Skip to content

Commit

Permalink
refactor how intracellular updates are handled
Browse files Browse the repository at this point in the history
- first do all pre updates, then all updates, then all post updates
- this ensures that any changes manifested in updates/post-updates are not immediately
  visible to the pre-updates of other cells
- with the more verbose logic controlling intracellular updates, make a separate function to handle it
    - in that function, encode directly dependence on `diffusion_dt` rather than the argument `diffusion_dt_`
    - `initialzed` used to determine when to call `intracellular` rather than passed in
    - if after first looping over all cells, none were found to require intracellular updates, then return early
  • Loading branch information
drbergman committed Feb 18, 2025
1 parent cbf625e commit e7fdcaf
Show file tree
Hide file tree
Showing 2 changed files with 47 additions and 19 deletions.
62 changes: 44 additions & 18 deletions core/PhysiCell_cell_container.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -140,24 +140,8 @@ void Cell_Container::update_all_cells(double t, double phenotype_dt_ , double me
static double mechanics_dt_tolerance = 0.001 * mechanics_dt_;

// intracellular update. called for every diffusion_dt, but actually depends on the intracellular_dt of each cell (as it can be noisy)

#pragma omp parallel for
for( int i=0; i < (*all_cells).size(); i++ )
{
if( (*all_cells)[i]->is_out_of_domain == false && initialzed ) {

if( (*all_cells)[i]->phenotype.intracellular != NULL && (*all_cells)[i]->phenotype.intracellular->need_update())
{
if ((*all_cells)[i]->functions.pre_update_intracellular != NULL)
(*all_cells)[i]->functions.pre_update_intracellular( (*all_cells)[i], (*all_cells)[i]->phenotype , diffusion_dt_ );

(*all_cells)[i]->phenotype.intracellular->update( (*all_cells)[i], (*all_cells)[i]->phenotype , diffusion_dt_ );

if ((*all_cells)[i]->functions.post_update_intracellular != NULL)
(*all_cells)[i]->functions.post_update_intracellular( (*all_cells)[i], (*all_cells)[i]->phenotype , diffusion_dt_ );
}
}
}
if (initialzed)
{ update_all_cells_intracellular(); }

if( time_since_last_cycle > phenotype_dt_ - 0.5 * diffusion_dt_ || !initialzed )
{
Expand Down Expand Up @@ -318,6 +302,48 @@ void Cell_Container::update_all_cells(double t, double phenotype_dt_ , double me
return;
}

void Cell_Container::update_all_cells_intracellular( void )
{
bool anything_to_update = false;
std::vector<bool> ready_to_update_intracellular = std::vector<bool>( (*all_cells).size(), false );
#pragma omp parallel for
for( int i=0; i < (*all_cells).size(); i++ )
{
if( (*all_cells)[i]->is_out_of_domain == false ) {
if( (*all_cells)[i]->phenotype.intracellular != NULL && (*all_cells)[i]->phenotype.intracellular->need_update())
{
ready_to_update_intracellular[i] = true;
anything_to_update = true;
if ((*all_cells)[i]->functions.pre_update_intracellular != NULL)
{
(*all_cells)[i]->functions.pre_update_intracellular((*all_cells)[i], (*all_cells)[i]->phenotype, diffusion_dt);
}
}
}
}

if (!anything_to_update)
{ return; }

#pragma omp parallel for
for ( int i=0; i < (*all_cells).size(); i++ )
{
if (ready_to_update_intracellular[i])
{
(*all_cells)[i]->phenotype.intracellular->update((*all_cells)[i], (*all_cells)[i]->phenotype, diffusion_dt);
}
}

#pragma omp parallel for
for ( int i=0; i < (*all_cells).size(); i++ )
{
if (ready_to_update_intracellular[i] && (*all_cells)[i]->functions.post_update_intracellular != NULL)
{
(*all_cells)[i]->functions.post_update_intracellular((*all_cells)[i], (*all_cells)[i]->phenotype, diffusion_dt);
}
}
}

void Cell_Container::register_agent( Cell* agent )
{
agent_grid[agent->get_current_mechanics_voxel_index()].push_back(agent);
Expand Down
4 changes: 3 additions & 1 deletion core/PhysiCell_cell_container.h
Original file line number Diff line number Diff line change
Expand Up @@ -104,7 +104,9 @@ class Cell_Container : public BioFVM::Agent_Container
void update_all_cells(double t);
void update_all_cells(double t, double dt);
void update_all_cells(double t, double phenotype_dt, double mechanics_dt);
void update_all_cells(double t, double phenotype_dt, double mechanics_dt, double diffusion_dt );
void update_all_cells(double t, double phenotype_dt, double mechanics_dt, double diffusion_dt );

void update_all_cells_intracellular( void );

void register_agent( Cell* agent );
void add_agent_to_outer_voxel(Cell* agent);
Expand Down

0 comments on commit e7fdcaf

Please sign in to comment.