00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 #include <cmath>
00029 #include <ctime>
00030 #include <iostream>
00031 #include <fstream>
00032 #include <set>
00033
00034 #include "TissueSimulation.hpp"
00035 #include "AbstractCellCentreBasedTissue.hpp"
00036 #include "CellBasedEventHandler.hpp"
00037 #include "LogFile.hpp"
00038
00039 template<unsigned DIM>
00040 TissueSimulation<DIM>::TissueSimulation(AbstractTissue<DIM>& rTissue,
00041 std::vector<AbstractForce<DIM>*> forceCollection,
00042 bool deleteTissueAndForceCollection,
00043 bool initialiseCells)
00044 : mEndTime(0.0),
00045 mrTissue(rTissue),
00046 mDeleteTissue(deleteTissueAndForceCollection),
00047 mAllocatedMemoryForForceCollection(deleteTissueAndForceCollection),
00048 mInitialiseCells(initialiseCells),
00049 mNoBirth(false),
00050 mUpdateTissue(true),
00051 mOutputDirectory(""),
00052 mSimulationOutputDirectory(mOutputDirectory),
00053 mNumBirths(0),
00054 mNumDeaths(0),
00055 mSamplingTimestepMultiple(1),
00056 mForceCollection(forceCollection)
00057 {
00058 mpConfig = TissueConfig::Instance();
00059
00060
00061 mpRandomGenerator = RandomNumberGenerator::Instance();
00062
00063
00064 if (dynamic_cast<AbstractCellCentreBasedTissue<DIM>*>(&rTissue))
00065 {
00066 mDt = 1.0/120.0;
00067 }
00068 else
00069 {
00070 mDt = 0.002;
00071 }
00072
00073 if (mInitialiseCells)
00074 {
00075 mrTissue.InitialiseCells();
00076 }
00077 }
00078
00079
00080 template<unsigned DIM>
00081 TissueSimulation<DIM>::~TissueSimulation()
00082 {
00083 if (mAllocatedMemoryForForceCollection)
00084 {
00085 for (typename std::vector<AbstractForce<DIM>*>::iterator force_iter = mForceCollection.begin();
00086 force_iter != mForceCollection.end();
00087 ++force_iter)
00088 {
00089 delete *force_iter;
00090 }
00091 }
00092
00093 if (mDeleteTissue)
00094 {
00095 for (typename std::vector<AbstractCellKiller<DIM>*>::iterator it=mCellKillers.begin();
00096 it != mCellKillers.end();
00097 ++it)
00098 {
00099 delete *it;
00100 }
00101 delete &mrTissue;
00102 }
00103 }
00104
00105
00106 template<unsigned DIM>
00107 unsigned TissueSimulation<DIM>::DoCellBirth()
00108 {
00109 if (mNoBirth)
00110 {
00111 return 0;
00112 }
00113
00114 unsigned num_births_this_step = 0;
00115
00116
00117 for (typename AbstractTissue<DIM>::Iterator cell_iter = mrTissue.Begin();
00118 cell_iter != mrTissue.End();
00119 ++cell_iter)
00120 {
00121
00122 if (cell_iter->GetAge() > 0.0)
00123 {
00124 if (cell_iter->ReadyToDivide())
00125 {
00126
00127 TissueCell new_cell = cell_iter->Divide();
00128
00129
00130 c_vector<double, DIM> new_location = CalculateCellDivisionVector(*cell_iter);
00131
00132
00133 mrTissue.AddCell(new_cell, new_location, &(*cell_iter));
00134
00135
00136 num_births_this_step++;
00137 }
00138 }
00139 }
00140 return num_births_this_step;
00141 }
00142
00143
00144 template<unsigned DIM>
00145 unsigned TissueSimulation<DIM>::DoCellRemoval()
00146 {
00147 unsigned num_deaths_this_step = 0;
00148
00149
00150
00151 for (typename std::vector<AbstractCellKiller<DIM>*>::iterator killer_iter = mCellKillers.begin();
00152 killer_iter != mCellKillers.end();
00153 ++killer_iter)
00154 {
00155 (*killer_iter)->TestAndLabelCellsForApoptosisOrDeath();
00156 }
00157
00158 num_deaths_this_step += mrTissue.RemoveDeadCells();
00159
00160 return num_deaths_this_step;
00161 }
00162
00163
00164 template<unsigned DIM>
00165 const std::vector<AbstractForce<DIM>*> TissueSimulation<DIM>::rGetForceCollection() const
00166 {
00167 return mForceCollection;
00168 }
00169
00170
00171 template<unsigned DIM>
00172 c_vector<double, DIM> TissueSimulation<DIM>::CalculateCellDivisionVector(TissueCell& rParentCell)
00173 {
00179 if (dynamic_cast<AbstractCellCentreBasedTissue<DIM>*>(&mrTissue))
00180 {
00181
00182 c_vector<double, DIM> parent_coords = mrTissue.GetLocationOfCellCentre(rParentCell);
00183 c_vector<double, DIM> daughter_coords;
00184
00185
00186 double separation = TissueConfig::Instance()->GetDivisionSeparation();
00187
00188
00189 c_vector<double, DIM> random_vector;
00190
00191
00192
00193
00194
00195
00196 switch (DIM)
00197 {
00198 case 1:
00199 {
00200 double random_direction = -1.0 + 2.0*(RandomNumberGenerator::Instance()->ranf() < 0.5);
00201
00202 random_vector(0) = 0.5*separation*random_direction;
00203 break;
00204 }
00205 case 2:
00206 {
00207 double random_angle = 2.0*M_PI*RandomNumberGenerator::Instance()->ranf();
00208
00209 random_vector(0) = 0.5*separation*cos(random_angle);
00210 random_vector(1) = 0.5*separation*sin(random_angle);
00211 break;
00212 }
00213 case 3:
00214 {
00215 double random_zenith_angle = M_PI*RandomNumberGenerator::Instance()->ranf();
00216 double random_azimuth_angle = 2*M_PI*RandomNumberGenerator::Instance()->ranf();
00217
00218 random_vector(0) = 0.5*separation*cos(random_azimuth_angle)*sin(random_zenith_angle);
00219 random_vector(1) = 0.5*separation*sin(random_azimuth_angle)*sin(random_zenith_angle);
00220 random_vector(2) = 0.5*separation*cos(random_zenith_angle);
00221 break;
00222 }
00223 default:
00224
00225 NEVER_REACHED;
00226 }
00227
00228 parent_coords = parent_coords - random_vector;
00229 daughter_coords = parent_coords + random_vector;
00230
00231
00232 ChastePoint<DIM> parent_coords_point(parent_coords);
00233 unsigned node_index = mrTissue.GetLocationIndexUsingCell(rParentCell);
00234 mrTissue.SetNode(node_index, parent_coords_point);
00235
00236 return daughter_coords;
00237 }
00238 else
00239 {
00240 return zero_vector<double>(DIM);
00241 }
00242 }
00243
00244
00245 template<unsigned DIM>
00246 void TissueSimulation<DIM>::UpdateNodePositions(const std::vector< c_vector<double, DIM> >& rNodeForces)
00247 {
00248
00249
00250
00251
00252 std::vector<c_vector<double, DIM> > old_node_locations;
00253 old_node_locations.reserve(mrTissue.GetNumNodes());
00254 for (unsigned node_index=0; node_index<mrTissue.GetNumNodes(); node_index++)
00255 {
00256 old_node_locations[node_index] = mrTissue.GetNode(node_index)->rGetLocation();
00257 }
00258
00259
00260 mrTissue.UpdateNodeLocations(rNodeForces, mDt);
00261
00262
00263 ApplyTissueBoundaryConditions(old_node_locations);
00264
00265
00266 if ( TissueConfig::Instance()->GetOutputNodeVelocities() )
00267 {
00268 if (SimulationTime::Instance()->GetTimeStepsElapsed()%mSamplingTimestepMultiple==0)
00269 {
00270 *mpNodeVelocitiesFile << SimulationTime::Instance()->GetTime() << "\t";
00271 for (unsigned node_index=0; node_index<mrTissue.GetNumNodes(); node_index++)
00272 {
00273
00274 assert(!mrTissue.GetNode(node_index)->IsDeleted());
00275
00276
00277 bool is_real_node = true;
00278 if (dynamic_cast<AbstractCellCentreBasedTissue<DIM>*>(&mrTissue))
00279 {
00280 if (static_cast<AbstractCellCentreBasedTissue<DIM>*>(&mrTissue)->IsGhostNode(node_index))
00281 {
00282
00283 is_real_node = false;
00284 }
00285 else
00286 {
00287
00288 assert(!mrTissue.rGetCellUsingLocationIndex(node_index).IsDead());
00289 }
00290 }
00291
00292
00293 if (is_real_node)
00294 {
00295 const c_vector<double,DIM>& position = mrTissue.GetNode(node_index)->rGetLocation();
00296 c_vector<double, 2> velocity = this->mDt * rNodeForces[node_index] / this->mrTissue.GetDampingConstant(node_index);
00297
00298 *mpNodeVelocitiesFile << node_index << " ";
00299 for (unsigned i=0; i<DIM; i++)
00300 {
00301 *mpNodeVelocitiesFile << position[i] << " ";
00302 }
00303 for (unsigned i=0; i<DIM; i++)
00304 {
00305 *mpNodeVelocitiesFile << velocity[i] << " ";
00306 }
00307 }
00308 }
00309 *mpNodeVelocitiesFile << "\n";
00310 }
00311 }
00312 }
00313
00314
00315 template<unsigned DIM>
00316 void TissueSimulation<DIM>::SetDt(double dt)
00317 {
00318 assert(dt > 0);
00319 mDt = dt;
00320 }
00321
00322
00323 template<unsigned DIM>
00324 double TissueSimulation<DIM>::GetDt()
00325 {
00326 return mDt;
00327 }
00328
00329
00330 template<unsigned DIM>
00331 unsigned TissueSimulation<DIM>::GetNumBirths()
00332 {
00333 return mNumBirths;
00334 }
00335
00336
00337 template<unsigned DIM>
00338 unsigned TissueSimulation<DIM>::GetNumDeaths()
00339 {
00340 return mNumDeaths;
00341 }
00342
00343
00344 template<unsigned DIM>
00345 void TissueSimulation<DIM>::SetEndTime(double endTime)
00346 {
00347 assert(endTime > 0);
00348 mEndTime = endTime;
00349 }
00350
00351
00352 template<unsigned DIM>
00353 void TissueSimulation<DIM>::SetOutputDirectory(std::string outputDirectory)
00354 {
00355 mOutputDirectory = outputDirectory;
00356 mSimulationOutputDirectory = mOutputDirectory;
00357 }
00358
00359
00360 template<unsigned DIM>
00361 std::string TissueSimulation<DIM>::GetOutputDirectory()
00362 {
00363 return mOutputDirectory;
00364 }
00365
00366 template<unsigned DIM>
00367 void TissueSimulation<DIM>::SetSamplingTimestepMultiple(unsigned samplingTimestepMultiple)
00368 {
00369 assert(samplingTimestepMultiple > 0);
00370 mSamplingTimestepMultiple = samplingTimestepMultiple;
00371 }
00372
00373
00374 template<unsigned DIM>
00375 AbstractTissue<DIM>& TissueSimulation<DIM>::rGetTissue()
00376 {
00377 return mrTissue;
00378 }
00379
00380
00381 template<unsigned DIM>
00382 const AbstractTissue<DIM>& TissueSimulation<DIM>::rGetTissue() const
00383 {
00384 return mrTissue;
00385 }
00386
00387
00388 template<unsigned DIM>
00389 void TissueSimulation<DIM>::SetUpdateTissueRule(bool updateTissue)
00390 {
00391 mUpdateTissue = updateTissue;
00392 }
00393
00394
00395 template<unsigned DIM>
00396 void TissueSimulation<DIM>::SetNoBirth(bool noBirth)
00397 {
00398 mNoBirth = noBirth;
00399 }
00400
00401
00402 template<unsigned DIM>
00403 void TissueSimulation<DIM>::AddCellKiller(AbstractCellKiller<DIM>* pCellKiller)
00404 {
00405 mCellKillers.push_back(pCellKiller);
00406 }
00407
00408
00409 template<unsigned DIM>
00410 std::vector<double> TissueSimulation<DIM>::GetNodeLocation(const unsigned& rNodeIndex)
00411 {
00412 std::vector<double> location;
00413 for (unsigned i=0; i<DIM; i++)
00414 {
00415 location.push_back( mrTissue.GetNode(rNodeIndex)->rGetLocation()[i] );
00416 }
00417 return location;
00418 }
00419
00420
00421 template<unsigned DIM>
00422 void TissueSimulation<DIM>::Solve()
00423 {
00424 CellBasedEventHandler::BeginEvent(CellBasedEventHandler::EVERYTHING);
00425 CellBasedEventHandler::BeginEvent(CellBasedEventHandler::SETUP);
00426
00427
00428 SimulationTime* p_simulation_time = SimulationTime::Instance();
00429 double current_time = p_simulation_time->GetTime();
00430
00431 unsigned num_time_steps = (unsigned) ((mEndTime-current_time)/mDt+0.5);
00432
00433 if (current_time > 0)
00434 {
00435 p_simulation_time->ResetEndTimeAndNumberOfTimeSteps(mEndTime, num_time_steps);
00436 }
00437 else
00438 {
00439 p_simulation_time->SetEndTimeAndNumberOfTimeSteps(mEndTime, num_time_steps);
00440 }
00441
00442 if (mOutputDirectory=="")
00443 {
00444 EXCEPTION("OutputDirectory not set");
00445 }
00446
00447 double time_now = p_simulation_time->GetTime();
00448 std::ostringstream time_string;
00449 time_string << time_now;
00450
00451 std::string results_directory = mOutputDirectory +"/results_from_time_" + time_string.str();
00452 mSimulationOutputDirectory = results_directory;
00453
00455
00457
00458
00459 OutputFileHandler output_file_handler(results_directory+"/", true);
00460
00461 mrTissue.CreateOutputFiles(results_directory+"/", false);
00462
00463 mpSetupFile = output_file_handler.OpenOutputFile("results.vizsetup");
00464
00465 if (TissueConfig::Instance()->GetOutputNodeVelocities())
00466 {
00467 OutputFileHandler output_file_handler2(results_directory+"/", false);
00468 mpNodeVelocitiesFile = output_file_handler2.OpenOutputFile("nodevelocities.dat");
00469 }
00470
00471 SetupSolve();
00472
00473
00474
00475 LOG(1, "Setting up cells...");
00476 for (typename AbstractTissue<DIM>::Iterator cell_iter = mrTissue.Begin();
00477 cell_iter != mrTissue.End();
00478 ++cell_iter)
00479 {
00480
00481
00482 cell_iter->ReadyToDivide();
00483 }
00484 LOG(1, "\tdone\n");
00485
00486
00487 WriteVisualizerSetupFile();
00488
00489 *mpSetupFile << std::flush;
00490
00491 mrTissue.WriteResultsToFiles();
00492
00493 CellBasedEventHandler::EndEvent(CellBasedEventHandler::SETUP);
00494
00495
00496 std::vector<c_vector<double, DIM> > forces(mrTissue.GetNumNodes(),zero_vector<double>(DIM));
00497
00499
00501
00502 while ((p_simulation_time->GetTimeStepsElapsed() < num_time_steps) && !(StoppingEventHasOccurred()) )
00503 {
00504 LOG(1, "--TIME = " << p_simulation_time->GetTime() << "\n");
00505
00512 UpdateTissue();
00513
00515
00517 CellBasedEventHandler::BeginEvent(CellBasedEventHandler::FORCE);
00518
00519
00520 for (unsigned i=0; i<forces.size(); i++)
00521 {
00522 forces[i].clear();
00523 }
00524
00525
00526
00527 if (mrTissue.GetNumNodes()!=forces.size())
00528 {
00529 forces.resize(mrTissue.GetNumNodes(), zero_vector<double>(DIM));
00530 }
00531
00532
00533 for (typename std::vector<AbstractForce<DIM>*>::iterator iter = mForceCollection.begin();
00534 iter !=mForceCollection.end();
00535 iter++)
00536 {
00537 (*iter)->AddForceContribution(forces, mrTissue);
00538 }
00539 CellBasedEventHandler::EndEvent(CellBasedEventHandler::FORCE);
00540
00542
00544 CellBasedEventHandler::BeginEvent(CellBasedEventHandler::POSITION);
00545 UpdateNodePositions(forces);
00546 CellBasedEventHandler::EndEvent(CellBasedEventHandler::POSITION);
00547
00549
00550
00552 PostSolve();
00553
00554
00555 p_simulation_time->IncrementTimeOneStep();
00556
00558
00560 CellBasedEventHandler::BeginEvent(CellBasedEventHandler::OUTPUT);
00561
00562
00563 if (p_simulation_time->GetTimeStepsElapsed()%mSamplingTimestepMultiple==0)
00564 {
00565 mrTissue.WriteResultsToFiles();
00566 }
00567
00568 CellBasedEventHandler::EndEvent(CellBasedEventHandler::OUTPUT);
00569 }
00570
00571 LOG(1, "--END TIME = " << SimulationTime::Instance()->GetTime() << "\n");
00572
00573
00574 UpdateTissue();
00575
00576 AfterSolve();
00577
00578 CellBasedEventHandler::BeginEvent(CellBasedEventHandler::OUTPUT);
00579
00580 mrTissue.CloseOutputFiles();
00581
00582 if (TissueConfig::Instance()->GetOutputNodeVelocities())
00583 {
00584 mpNodeVelocitiesFile->close();
00585 }
00586
00587 *mpSetupFile << "Complete\n";
00588 mpSetupFile->close();
00589
00590 CellBasedEventHandler::EndEvent(CellBasedEventHandler::OUTPUT);
00591
00592 CellBasedEventHandler::EndEvent(CellBasedEventHandler::EVERYTHING);
00593 }
00594
00595
00596 template<unsigned DIM>
00597 bool TissueSimulation<DIM>::StoppingEventHasOccurred()
00598 {
00599 return false;
00600 }
00601
00602
00603 template<unsigned DIM>
00604 void TissueSimulation<DIM>::UpdateTissue()
00605 {
00607
00609 CellBasedEventHandler::BeginEvent(CellBasedEventHandler::DEATH);
00610 unsigned deaths_this_step = DoCellRemoval();
00611 mNumDeaths += deaths_this_step;
00612 LOG(1, "\tNum deaths = " << mNumDeaths << "\n");
00613 CellBasedEventHandler::EndEvent(CellBasedEventHandler::DEATH);
00614
00616
00618 CellBasedEventHandler::BeginEvent(CellBasedEventHandler::BIRTH);
00619 unsigned births_this_step = DoCellBirth();
00620 mNumBirths += births_this_step;
00621 LOG(1, "\tNum births = " << mNumBirths << "\n");
00622 CellBasedEventHandler::EndEvent(CellBasedEventHandler::BIRTH);
00623
00624
00625 bool births_or_death_occurred = ((births_this_step>0) || (deaths_this_step>0));
00626
00628
00630 CellBasedEventHandler::BeginEvent(CellBasedEventHandler::UPDATETISSUE);
00631 if (mUpdateTissue)
00632 {
00633 LOG(1, "\tUpdating tissue...");
00634 mrTissue.Update(births_or_death_occurred);
00635 LOG(1, "\tdone.\n");
00636 }
00637 else if (births_or_death_occurred)
00638 {
00639 EXCEPTION("Tissue has had births or deaths but mUpdateTissue is set to false, please set it to true.");
00640 }
00641 CellBasedEventHandler::EndEvent(CellBasedEventHandler::UPDATETISSUE);
00642 }
00643
00645
00647
00648
00649 template class TissueSimulation<1>;
00650 template class TissueSimulation<2>;
00651 template class TissueSimulation<3>;