Skip to content

Commit

Permalink
More outputs for exposed but infectious, etc... (#37)
Browse files Browse the repository at this point in the history
* more outputs

* fix problem with resetting hospitalizations and ventilations
atmyers authored Mar 28, 2024
1 parent a3d5c70 commit 931222d
Showing 3 changed files with 64 additions and 21 deletions.
4 changes: 2 additions & 2 deletions src/AgentContainer.H
Original file line number Diff line number Diff line change
@@ -64,7 +64,7 @@ struct IntIdx
workgroup, /*!< workgroup ID */
work_nborhood, /*!< work neighborhood ID */
withdrawn, /*!< quarantine status */
symptomatic, /*!< currently symptomatic? */
symptomatic, /*!< currently symptomatic? 0: no, but will be, 1: yes, 2: no, and will remain so until recovered */
nattribs /*!< number of integer-type attribute */
};
};
@@ -199,7 +199,7 @@ public:

void generateCellData (amrex::MultiFab& mf) const;

std::array<amrex::Long, 5> getTotals ();
std::array<amrex::Long, 9> getTotals ();

const DiseaseParm * getDiseaseParameters_h () const {
return h_parm;
70 changes: 53 additions & 17 deletions src/AgentContainer.cpp
Original file line number Diff line number Diff line change
@@ -805,10 +805,13 @@ void AgentContainer::updateStatus (MultiFab& disease_stats /*!< Community-wise d
auto timer_ptr = soa.GetRealData(RealIdx::treatment_timer).data();
auto prob_ptr = soa.GetRealData(RealIdx::prob).data();
auto withdrawn_ptr = soa.GetIntData(IntIdx::withdrawn).data();
auto symptomatic_ptr = soa.GetIntData(IntIdx::symptomatic).data();
auto incubation_period_ptr = soa.GetRealData(RealIdx::incubation_period).data();
auto infectious_period_ptr = soa.GetRealData(RealIdx::infectious_period).data();
auto symptomdev_period_ptr = soa.GetRealData(RealIdx::symptomdev_period).data();

auto* lparm = d_parm;

auto ds_arr = disease_stats[mfi].array();

struct DiseaseStats
@@ -838,8 +841,18 @@ void AgentContainer::updateStatus (MultiFab& disease_stats /*!< Community-wise d
}
else if (status_ptr[i] == Status::infected) {
counter_ptr[i] += 1;
if (counter_ptr[i] == amrex::Math::ceil(symptomdev_period_ptr[i]) && symptomatic_withdraw) {
withdrawn_ptr[i] = 1;
if (amrex::Random(engine) < lparm->p_asymp[0]) {
symptomatic_ptr[i] = 2;
} else {
symptomatic_ptr[i] = 0;
}
if (counter_ptr[i] == amrex::Math::ceil(symptomdev_period_ptr[i])) {
if (symptomatic_ptr[i] != 2) {
symptomatic_ptr[i] = 1;
}
if (symptomatic_ptr[i] && symptomatic_withdraw) {
withdrawn_ptr[i] = 1;
}
}
if (counter_ptr[i] < incubation_period_ptr[i]) {
// incubation phase
@@ -910,6 +923,9 @@ void AgentContainer::updateStatus (MultiFab& disease_stats /*!< Community-wise d
//pstruct_ptr[i].id() = -pstruct_ptr[i].id();
}
}
amrex::Gpu::Atomic::AddNoRet(
&ds_arr(home_i_ptr[i], home_j_ptr[i], 0,
DiseaseStats::hospitalization), -1.0_rt);
amrex::Gpu::Atomic::AddNoRet(
&ds_arr(home_i_ptr[i], home_j_ptr[i], 0,
DiseaseStats::ICU), -1.0_rt);
@@ -925,7 +941,12 @@ void AgentContainer::updateStatus (MultiFab& disease_stats /*!< Community-wise d
status_ptr[i] = Status::dead;
//pstruct_ptr[i].id() = -pstruct_ptr[i].id();
}

amrex::Gpu::Atomic::AddNoRet(
&ds_arr(home_i_ptr[i], home_j_ptr[i], 0,
DiseaseStats::hospitalization), -1.0_rt);
amrex::Gpu::Atomic::AddNoRet(
&ds_arr(home_i_ptr[i], home_j_ptr[i], 0,
DiseaseStats::ICU), -1.0_rt);
amrex::Gpu::Atomic::AddNoRet(
&ds_arr(home_i_ptr[i], home_j_ptr[i], 0,
DiseaseStats::ventilator), -1.0_rt);
@@ -937,6 +958,8 @@ void AgentContainer::updateStatus (MultiFab& disease_stats /*!< Community-wise d
else { // not hospitalized, recover once not infectious
if (counter_ptr[i] >= (incubation_period_ptr[i] + infectious_period_ptr[i])) {
status_ptr[i] = Status::immune;
counter_ptr[i] = 0.0;
symptomatic_ptr[i] = 0;
withdrawn_ptr[i] = 0;
}
}
@@ -1458,27 +1481,40 @@ void AgentContainer::generateCellData (MultiFab& mf /*!< MultiFab with at least
Returns a vector with 5 components corresponding to each value of #Status; each element is
the total number of agents at a step with the corresponding #Status (in that order).
*/
std::array<Long, 5> AgentContainer::getTotals () {
std::array<Long, 9> AgentContainer::getTotals () {
BL_PROFILE("getTotals");
amrex::ReduceOps<ReduceOpSum, ReduceOpSum, ReduceOpSum, ReduceOpSum, ReduceOpSum> reduce_ops;
auto r = amrex::ParticleReduce<ReduceData<int,int,int,int,int>> (
amrex::ReduceOps<ReduceOpSum, ReduceOpSum, ReduceOpSum, ReduceOpSum, ReduceOpSum, ReduceOpSum, ReduceOpSum, ReduceOpSum, ReduceOpSum> reduce_ops;
auto r = amrex::ParticleReduce<ReduceData<int,int,int,int,int,int,int,int,int>> (
*this, [=] AMREX_GPU_DEVICE (const SuperParticleType& p) noexcept
-> amrex::GpuTuple<int,int,int,int,int>
-> amrex::GpuTuple<int,int,int,int,int,int,int,int,int>
{
int s[5] = {0, 0, 0, 0, 0};
int s[9] = {0, 0, 0, 0, 0, 0, 0, 0, 0};
AMREX_ALWAYS_ASSERT(p.idata(IntIdx::status) >= 0);
AMREX_ALWAYS_ASSERT(p.idata(IntIdx::status) <= 4);
s[p.idata(IntIdx::status)] = 1;
return {s[0], s[1], s[2], s[3], s[4]};
if (p.idata(IntIdx::status) == 1) { // exposed
if (p.rdata(RealIdx::disease_counter) <= amrex::Math::ceil(p.rdata(RealIdx::incubation_period))) {
s[5] = 1; // exposed, but not infectious
} else { // infectious
if (p.idata(IntIdx::symptomatic) == 2) {
s[6] = 1; // asymptomatic and will remain so
}
else if (p.idata(IntIdx::symptomatic) == 0) {
s[7] = 1; // asymptomatic but will develop symptoms
}
else if (p.idata(IntIdx::symptomatic) == 1) {
s[8] = 1; // Infectious and symptomatic
} else {
amrex::Abort("how did I get here?");
}
}
}
return {s[0], s[1], s[2], s[3], s[4], s[5], s[6], s[7], s[8]};
}, reduce_ops);

std::array<Long, 5> counts = {amrex::get<0>(r), amrex::get<1>(r), amrex::get<2>(r), amrex::get<3>(r),
amrex::get<4>(r)};
ParallelDescriptor::ReduceLongSum(&counts[0], 5, ParallelDescriptor::IOProcessorNumber());
// amrex::Print() << "Never infected: " << counts[0] << "\n";
// amrex::Print() << "Infected: " << counts[1] << "\n";
// amrex::Print() << "Immune: " << counts[2] << "\n";
// amrex::Print() << "Previously infected: " << counts[3] << "\n";
// amrex::Print() << "Deaths: " << counts[4] << "\n";
std::array<Long, 9> counts = {amrex::get<0>(r), amrex::get<1>(r), amrex::get<2>(r), amrex::get<3>(r),
amrex::get<4>(r), amrex::get<5>(r), amrex::get<6>(r), amrex::get<7>(r),
amrex::get<8>(r)};
ParallelDescriptor::ReduceLongSum(&counts[0], 9, ParallelDescriptor::IOProcessorNumber());
return counts;
}
11 changes: 9 additions & 2 deletions src/main.cpp
Original file line number Diff line number Diff line change
@@ -125,7 +125,7 @@ void runAgent ()
amrex::FileOpenFailed(output_filename);
}

File << "Day Never Infected Immune Deaths Hospitalized Ventilated ICU\n";
File << "Day Never Infected Immune Deaths Hospitalized Ventilated ICU Exposed Asymptomatic Presymptomatic Symptomatic\n";

File.flush();

@@ -240,14 +240,21 @@ void runAgent ()
// total number of deaths computed on agents and on mesh should be the same...
AMREX_ALWAYS_ASSERT(mmc[3] == counts[4]);

// the total number of infected should equal the sum of
// exposed but not infectious
// infectious and asymptomatic
// infectious and pre-symptomatic
// infectious and symptomatic
AMREX_ALWAYS_ASSERT(counts[1] == counts[5] + counts[6] + counts[7] + counts[8]);

std::ofstream File;
File.open(output_filename.c_str(), std::ios::out|std::ios::app);

if (!File.good()) {
amrex::FileOpenFailed(output_filename);
}

File << i << " " << counts[0] << " " << counts[1] << " " << counts[2] << " " << counts[4] << " " << mmc[0] << " " << mmc[1] << " " << mmc[2] << "\n";
File << i << " " << counts[0] << " " << counts[1] << " " << counts[2] << " " << counts[4] << " " << mmc[0] << " " << mmc[1] << " " << mmc[2] << " " << counts[5] << " " << counts[6] << " " << counts[7] << " " << counts[8] << "\n";

File.flush();

0 comments on commit 931222d

Please sign in to comment.