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
00029
00030
00031
00032
00033
00034
00035
00036 #include "CellBasedPdeHandler.hpp"
00037 #include "MeshBasedCellPopulationWithGhostNodes.hpp"
00038 #include "NodeBasedCellPopulation.hpp"
00039 #include "PottsBasedCellPopulation.hpp"
00040 #include "VertexBasedCellPopulation.hpp"
00041 #include "CaBasedCellPopulation.hpp"
00042 #include "SimpleLinearEllipticSolver.hpp"
00043 #include "CellBasedPdeSolver.hpp"
00044 #include "Exception.hpp"
00045 #include "VtkMeshWriter.hpp"
00046
00047 template<unsigned DIM>
00048 CellBasedPdeHandler<DIM>::CellBasedPdeHandler(AbstractCellPopulation<DIM>* pCellPopulation,
00049 bool deleteMemberPointersInDestructor)
00050 : mpCellPopulation(pCellPopulation),
00051 mWriteAverageRadialPdeSolution(false),
00052 mWriteDailyAverageRadialPdeSolution(false),
00053 mAverageRadialSolutionVariableName(""),
00054 mSetBcsOnCoarseBoundary(true),
00055 mNumRadialIntervals(UNSIGNED_UNSET),
00056 mpCoarsePdeMesh(NULL),
00057 mDeleteMemberPointersInDestructor(deleteMemberPointersInDestructor)
00058 {
00059
00061 assert(mpCellPopulation->GetNumRealCells() != 0);
00062 }
00063
00064 template<unsigned DIM>
00065 CellBasedPdeHandler<DIM>::~CellBasedPdeHandler()
00066 {
00067
00068
00069
00070
00071
00072
00073 if (mDeleteMemberPointersInDestructor)
00074 {
00075 for (unsigned i=0; i<mPdeAndBcCollection.size(); i++)
00076 {
00077 delete mPdeAndBcCollection[i];
00078 }
00079 }
00080 if (mpCoarsePdeMesh)
00081 {
00082 delete mpCoarsePdeMesh;
00083 }
00084 }
00085
00086 template<unsigned DIM>
00087 const AbstractCellPopulation<DIM>* CellBasedPdeHandler<DIM>::GetCellPopulation() const
00088 {
00089 return mpCellPopulation;
00090 }
00091
00092 template<unsigned DIM>
00093 TetrahedralMesh<DIM,DIM>* CellBasedPdeHandler<DIM>::GetCoarsePdeMesh()
00094 {
00095 return mpCoarsePdeMesh;
00096 }
00097
00098 template<unsigned DIM>
00099 void CellBasedPdeHandler<DIM>::AddPdeAndBc(PdeAndBoundaryConditions<DIM>* pPdeAndBc)
00100 {
00101 if (!mPdeAndBcCollection.empty() && (pPdeAndBc->rGetDependentVariableName()==""))
00102 {
00103 EXCEPTION("When adding more than one PDE to CellBasedPdeHandler set the dependent variable name using SetDependentVariableName(name).");
00104 }
00105 for (unsigned i=0; i< mPdeAndBcCollection.size(); i++)
00106 {
00107 if (pPdeAndBc->rGetDependentVariableName() == mPdeAndBcCollection[i]->rGetDependentVariableName())
00108 {
00109 EXCEPTION("The name "+ pPdeAndBc->rGetDependentVariableName()+ " has already been used in the PDE collection");
00110 }
00111 }
00112 mPdeAndBcCollection.push_back(pPdeAndBc);
00113 }
00114
00115 template<unsigned DIM>
00116 Vec CellBasedPdeHandler<DIM>::GetPdeSolution(const std::string& rName)
00117 {
00118 if (rName != "")
00119 {
00120 for (unsigned i=0; i<mPdeAndBcCollection.size(); i++)
00121 {
00122 if (mPdeAndBcCollection[i]->rGetDependentVariableName() == rName)
00123 {
00124 return mPdeAndBcCollection[i]->GetSolution();
00125 }
00126 }
00127
00128 EXCEPTION("The PDE collection does not contain a PDE named " + rName);
00129 }
00130 else
00131 {
00132 assert(mPdeAndBcCollection.size()==1);
00133 return mPdeAndBcCollection[0]->GetSolution();
00134 }
00135 }
00136
00137 template<unsigned DIM>
00138 void CellBasedPdeHandler<DIM>::InitialiseCellPdeElementMap()
00139 {
00140 if (mpCoarsePdeMesh == NULL)
00141 {
00142 EXCEPTION("InitialiseCellPdeElementMap() should only be called if mpCoarsePdeMesh is set up.");
00143 }
00144
00145 mCellPdeElementMap.clear();
00146
00147
00148 for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
00149 cell_iter != mpCellPopulation->End();
00150 ++cell_iter)
00151 {
00152 const ChastePoint<DIM>& r_position_of_cell = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
00153 unsigned elem_index = mpCoarsePdeMesh->GetContainingElementIndex(r_position_of_cell);
00154 mCellPdeElementMap[*cell_iter] = elem_index;
00155 }
00156 }
00157
00158 template<unsigned DIM>
00159 void CellBasedPdeHandler<DIM>::UpdateCellPdeElementMap()
00160 {
00161
00162 for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
00163 cell_iter != mpCellPopulation->End();
00164 ++cell_iter)
00165 {
00166 const ChastePoint<DIM>& r_position_of_cell = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
00167 unsigned elem_index = mpCoarsePdeMesh->GetContainingElementIndexWithInitialGuess(r_position_of_cell, mCellPdeElementMap[*cell_iter]);
00168 mCellPdeElementMap[*cell_iter] = elem_index;
00169 }
00170 }
00171
00172 template<unsigned DIM>
00173 void CellBasedPdeHandler<DIM>::OpenResultsFiles(std::string outputDirectory)
00174 {
00175
00176 if ((dynamic_cast<CaBasedCellPopulation<DIM>*>(mpCellPopulation) != NULL) && mpCoarsePdeMesh==NULL)
00177 {
00178 assert(DIM ==2);
00179 ChasteCuboid<DIM> cuboid = mpCellPopulation->rGetMesh().CalculateBoundingBox();
00180
00181
00182 assert(cuboid.GetWidth(0) == cuboid.GetWidth(1));
00183
00184 UseCoarsePdeMesh(1, cuboid, false);
00185 }
00186
00187
00188 if (PdeSolveNeedsCoarseMesh() && mpCoarsePdeMesh==NULL)
00189 {
00190 EXCEPTION("Trying to solve a PDE on a cell population that doesn't have a mesh. Try calling UseCoarsePdeMesh().");
00191 }
00192
00193 if (mpCoarsePdeMesh != NULL)
00194 {
00195 InitialiseCellPdeElementMap();
00196
00197
00198 TrianglesMeshWriter<DIM,DIM> mesh_writer(outputDirectory+"/coarse_mesh_output", "coarse_mesh",false);
00199 mesh_writer.WriteFilesUsingMesh(*mpCoarsePdeMesh);
00200 }
00201
00202 if (PetscTools::AmMaster())
00203 {
00204 OutputFileHandler output_file_handler(outputDirectory+"/", false);
00205
00206 if (mpCoarsePdeMesh != NULL)
00207 {
00208 mpVizPdeSolutionResultsFile = output_file_handler.OpenOutputFile("results.vizcoarsepdesolution");
00209 }
00210 else
00211 {
00212 mpVizPdeSolutionResultsFile = output_file_handler.OpenOutputFile("results.vizpdesolution");
00213 }
00214
00215 if (mWriteAverageRadialPdeSolution)
00216 {
00217 mpAverageRadialPdeSolutionResultsFile = output_file_handler.OpenOutputFile("radial_dist.dat");
00218 }
00219 }
00220
00221 mDirPath = outputDirectory;
00222
00223
00224 }
00225
00226 template<unsigned DIM>
00227 void CellBasedPdeHandler<DIM>::CloseResultsFiles()
00228 {
00229
00230 if (PetscTools::AmMaster())
00231 {
00232 mpVizPdeSolutionResultsFile->close();
00233 if (mWriteAverageRadialPdeSolution)
00234 {
00235 WriteAverageRadialPdeSolution(SimulationTime::Instance()->GetTime());
00236 mpAverageRadialPdeSolutionResultsFile->close();
00237 }
00238 }
00239 }
00240
00241 template<unsigned DIM>
00242 void CellBasedPdeHandler<DIM>::UseCoarsePdeMesh(double stepSize, ChasteCuboid<DIM> meshCuboid, bool centreOnCellPopulation)
00243 {
00244 if (mPdeAndBcCollection.empty())
00245 {
00246 EXCEPTION("mPdeAndBcCollection should be populated prior to calling UseCoarsePdeMesh().");
00247 }
00248
00249 for (unsigned pde_index=0; pde_index<mPdeAndBcCollection.size(); pde_index++)
00250 {
00251 if (mPdeAndBcCollection[pde_index]->HasAveragedSourcePde() == false && !dynamic_cast<CaBasedCellPopulation<DIM>*>(mpCellPopulation))
00252 {
00253 EXCEPTION("UseCoarsePdeMesh() should only be called if averaged-source PDEs are specified.");
00254 }
00255 }
00256
00257
00258 mpCoarsePdeMesh = new TetrahedralMesh<DIM,DIM>;
00259 switch (DIM)
00260 {
00261 case 1:
00262 mpCoarsePdeMesh->ConstructRegularSlabMesh(stepSize, meshCuboid.GetWidth(0));
00263 break;
00264 case 2:
00265 mpCoarsePdeMesh->ConstructRegularSlabMesh(stepSize, meshCuboid.GetWidth(0), meshCuboid.GetWidth(1));
00266 break;
00267 case 3:
00268 mpCoarsePdeMesh->ConstructRegularSlabMesh(stepSize, meshCuboid.GetWidth(0), meshCuboid.GetWidth(1), meshCuboid.GetWidth(2));
00269 break;
00270 default:
00271 NEVER_REACHED;
00272 }
00273
00274 if (centreOnCellPopulation)
00275 {
00276
00277 c_vector<double,DIM> centre_of_coarse_mesh = zero_vector<double>(DIM);
00278 for (unsigned i=0; i<mpCoarsePdeMesh->GetNumNodes(); i++)
00279 {
00280 centre_of_coarse_mesh += mpCoarsePdeMesh->GetNode(i)->rGetLocation();
00281 }
00282 centre_of_coarse_mesh /= mpCoarsePdeMesh->GetNumNodes();
00283
00284
00285 c_vector<double,DIM> centre_of_cell_population = mpCellPopulation->GetCentroidOfCellPopulation();
00286 mpCoarsePdeMesh->Translate(centre_of_cell_population - centre_of_coarse_mesh);
00287 }
00288 else
00289 {
00290
00291 ChastePoint<DIM> upper = meshCuboid.rGetUpperCorner();
00292 ChastePoint<DIM> lower = meshCuboid.rGetLowerCorner();
00293 c_vector<double,DIM> centre_of_cuboid = 0.5*(upper.rGetLocation() + lower.rGetLocation());
00294
00295
00296 c_vector<double,DIM> centre_of_coarse_mesh = zero_vector<double>(DIM);
00297 for (unsigned i=0; i<mpCoarsePdeMesh->GetNumNodes(); i++)
00298 {
00299 centre_of_coarse_mesh += mpCoarsePdeMesh->GetNode(i)->rGetLocation();
00300 }
00301 centre_of_coarse_mesh /= mpCoarsePdeMesh->GetNumNodes();
00302
00303 mpCoarsePdeMesh->Translate(centre_of_cuboid - centre_of_coarse_mesh);
00304 }
00305 }
00306
00307 template<unsigned DIM>
00308 void CellBasedPdeHandler<DIM>::SolvePdeAndWriteResultsToFile(unsigned samplingTimestepMultiple)
00309 {
00310
00311 bool using_coarse_pde_mesh = (mpCoarsePdeMesh != NULL);
00312
00313
00314 assert(!mPdeAndBcCollection.empty());
00315 for (unsigned pde_index=0; pde_index<mPdeAndBcCollection.size(); pde_index++)
00316 {
00317 assert(mPdeAndBcCollection[pde_index]);
00318 assert(mPdeAndBcCollection[pde_index]->HasAveragedSourcePde() == using_coarse_pde_mesh || dynamic_cast<CaBasedCellPopulation<DIM>*>(mpCellPopulation));
00319 }
00320
00321
00322 mpCellPopulation->Update();
00323
00324
00325 TetrahedralMesh<DIM,DIM>* p_mesh;
00326 if (using_coarse_pde_mesh)
00327 {
00328 p_mesh = mpCoarsePdeMesh;
00329 }
00330 else
00331 {
00332
00333 p_mesh = &(static_cast<MeshBasedCellPopulation<DIM>*>(mpCellPopulation)->rGetMesh());
00334 }
00335
00336
00337 for (unsigned pde_index=0; pde_index<mPdeAndBcCollection.size(); pde_index++)
00338 {
00339
00340 PdeAndBoundaryConditions<DIM>* p_pde_and_bc = mPdeAndBcCollection[pde_index];
00341
00342
00343 std::auto_ptr<BoundaryConditionsContainer<DIM,DIM,1> > p_bcc = ConstructBoundaryConditionsContainer(p_pde_and_bc, p_mesh);
00344
00345
00346 PetscInt previous_solution_size = 0;
00347 if (p_pde_and_bc->GetSolution())
00348 {
00349 VecGetSize(p_pde_and_bc->GetSolution(), &previous_solution_size);
00350 }
00351
00352
00353 bool is_previous_solution_size_correct = (previous_solution_size == (int)p_mesh->GetNumNodes());
00354
00355
00356 Vec initial_guess;
00357 if (is_previous_solution_size_correct)
00358 {
00359
00360 VecDuplicate(p_pde_and_bc->GetSolution(), &initial_guess);
00361 VecCopy(p_pde_and_bc->GetSolution(), initial_guess);
00362 p_pde_and_bc->DestroySolution();
00363 }
00364 else
00365 {
00367 if (!using_coarse_pde_mesh && p_pde_and_bc->GetSolution())
00368 {
00369 assert(previous_solution_size != 0);
00370 p_pde_and_bc->DestroySolution();
00371 }
00372 }
00373
00374
00375 if (p_pde_and_bc->HasAveragedSourcePde())
00376 {
00377
00378
00379 this->UpdateCellPdeElementMap();
00380 p_pde_and_bc->SetUpSourceTermsForAveragedSourcePde(p_mesh, &mCellPdeElementMap);
00381
00382
00383 SimpleLinearEllipticSolver<DIM,DIM> solver(p_mesh, p_pde_and_bc->GetPde(), p_bcc.get());
00384
00385
00386 if (is_previous_solution_size_correct)
00387 {
00388 p_pde_and_bc->SetSolution(solver.Solve(initial_guess));
00389 PetscTools::Destroy(initial_guess);
00390 }
00391 else
00392 {
00393 p_pde_and_bc->SetSolution(solver.Solve());
00394 }
00395 }
00396 else
00397 {
00398 CellBasedPdeSolver<DIM> solver(p_mesh, p_pde_and_bc->GetPde(), p_bcc.get());
00399
00400
00401 if (is_previous_solution_size_correct)
00402 {
00403 p_pde_and_bc->SetSolution(solver.Solve(initial_guess));
00404 PetscTools::Destroy(initial_guess);
00405 }
00406 else
00407 {
00408 p_pde_and_bc->SetSolution(solver.Solve());
00409 }
00410 }
00411
00412
00413 ReplicatableVector solution_repl(p_pde_and_bc->GetSolution());
00414
00415
00416 for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
00417 cell_iter != mpCellPopulation->End();
00418 ++cell_iter)
00419 {
00420 unsigned node_index = mpCellPopulation->GetLocationIndexUsingCell(*cell_iter);
00421 double solution_at_node = 0.0;
00422
00423 if (using_coarse_pde_mesh)
00424 {
00425
00426
00427
00428 unsigned elem_index = mCellPdeElementMap[*cell_iter];
00429 Element<DIM,DIM>* p_element = mpCoarsePdeMesh->GetElement(elem_index);
00430
00431 const ChastePoint<DIM>& node_location = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
00432
00433 c_vector<double,DIM+1> weights = p_element->CalculateInterpolationWeights(node_location);
00434 for (unsigned i=0; i<DIM+1; i++)
00435 {
00436 double nodal_value = solution_repl[p_element->GetNodeGlobalIndex(i)];
00437 solution_at_node += nodal_value * weights(i);
00438 }
00439 }
00440 else
00441 {
00442 solution_at_node = solution_repl[node_index];
00443 }
00444 cell_iter->GetCellData()->SetItem(mPdeAndBcCollection[pde_index]->rGetDependentVariableName(), solution_at_node);
00445 }
00446 }
00447
00448
00449 SimulationTime* p_time = SimulationTime::Instance();
00450 if ((p_time->GetTimeStepsElapsed())%samplingTimestepMultiple == 0)
00451 {
00452 WritePdeSolution(p_time->GetTime());
00453 }
00454 #define COVERAGE_IGNORE
00456 if (!using_coarse_pde_mesh)
00457 {
00458 if (mWriteDailyAverageRadialPdeSolution)
00459 {
00461 p_time = SimulationTime::Instance();
00462 unsigned num_timesteps_per_day = (unsigned) (DBL_EPSILON + 24/SimulationTime::Instance()->GetTimeStep());
00463 if ((p_time->GetTimeStepsElapsed()) % num_timesteps_per_day == 0)
00464 {
00465 WriteAverageRadialPdeSolution(p_time->GetTime());
00466 }
00467 }
00468 }
00469 #undef COVERAGE_IGNORE
00470 }
00471
00472 template<unsigned DIM>
00473 std::auto_ptr<BoundaryConditionsContainer<DIM,DIM,1> > CellBasedPdeHandler<DIM>::ConstructBoundaryConditionsContainer(
00474 PdeAndBoundaryConditions<DIM>* pPdeAndBc,
00475 TetrahedralMesh<DIM,DIM>* pMesh)
00476 {
00477 std::auto_ptr<BoundaryConditionsContainer<DIM,DIM,1> > p_bcc(new BoundaryConditionsContainer<DIM,DIM,1>(false));
00478
00479 AbstractBoundaryCondition<DIM>* p_bc = pPdeAndBc->GetBoundaryCondition();
00480
00481 if (pPdeAndBc->IsNeumannBoundaryCondition())
00482 {
00483
00484 for (typename TetrahedralMesh<DIM,DIM>::BoundaryElementIterator elem_iter = pMesh->GetBoundaryElementIteratorBegin();
00485 elem_iter != pMesh->GetBoundaryElementIteratorEnd();
00486 ++elem_iter)
00487 {
00488 p_bcc->AddNeumannBoundaryCondition(*elem_iter, p_bc);
00489 }
00490 }
00491 else
00492 {
00493 bool using_coarse_pde_mesh = (mpCoarsePdeMesh != NULL);
00494
00495 if (using_coarse_pde_mesh && !mSetBcsOnCoarseBoundary)
00496 {
00497
00498 std::set<unsigned> coarse_element_indices_in_map;
00499 for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
00500 cell_iter != mpCellPopulation->End();
00501 ++cell_iter)
00502 {
00503 coarse_element_indices_in_map.insert(mCellPdeElementMap[*cell_iter]);
00504 }
00505
00506
00507 std::set<unsigned> coarse_mesh_boundary_node_indices;
00508 for (unsigned i=0; i<pMesh->GetNumElements(); i++)
00509 {
00510 if (coarse_element_indices_in_map.find(i) == coarse_element_indices_in_map.end())
00511 {
00512 Element<DIM,DIM>* p_element = pMesh->GetElement(i);
00513 for (unsigned j=0; j<DIM+1; j++)
00514 {
00515 unsigned node_index = p_element->GetNodeGlobalIndex(j);
00516 coarse_mesh_boundary_node_indices.insert(node_index);
00517 }
00518 }
00519 }
00520
00521
00522 for (std::set<unsigned>::iterator iter = coarse_mesh_boundary_node_indices.begin();
00523 iter != coarse_mesh_boundary_node_indices.end();
00524 ++iter)
00525 {
00526 p_bcc->AddDirichletBoundaryCondition(pMesh->GetNode(*iter), p_bc, 0, false);
00527 }
00528 }
00529 else
00530 {
00531 for (typename TetrahedralMesh<DIM,DIM>::BoundaryNodeIterator node_iter = pMesh->GetBoundaryNodeIteratorBegin();
00532 node_iter != pMesh->GetBoundaryNodeIteratorEnd();
00533 ++node_iter)
00534 {
00535 p_bcc->AddDirichletBoundaryCondition(*node_iter, p_bc);
00536 }
00537 }
00538 }
00539
00540 return p_bcc;
00541 }
00542
00543 template<unsigned DIM>
00544 unsigned CellBasedPdeHandler<DIM>::FindCoarseElementContainingCell(CellPtr pCell)
00545 {
00546
00547 unsigned old_element_index = mCellPdeElementMap[pCell];
00548
00549
00550 std::set<unsigned> test_elements;
00551 test_elements.insert(old_element_index);
00552
00553 Element<DIM,DIM>* p_element = mpCoarsePdeMesh->GetElement(old_element_index);
00554 for (unsigned local_index=0; local_index<DIM+1; local_index++)
00555 {
00556 std::set<unsigned> element_indices = p_element->GetNode(local_index)->rGetContainingElementIndices();
00557 for (std::set<unsigned>::iterator iter = element_indices.begin();
00558 iter != element_indices.end();
00559 ++iter)
00560 {
00561 test_elements.insert(*iter);
00562 }
00563 }
00564
00565
00566 const ChastePoint<DIM>& r_cell_position = mpCellPopulation->GetLocationOfCellCentre(pCell);
00567 unsigned new_element_index = mpCoarsePdeMesh->GetContainingElementIndex(r_cell_position, false, test_elements);
00568
00569
00570 mCellPdeElementMap[pCell] = new_element_index;
00571
00572 return new_element_index;
00573 }
00574
00575 template<unsigned DIM>
00576 void CellBasedPdeHandler<DIM>::WritePdeSolution(double time)
00577 {
00578 if (PetscTools::AmMaster())
00579 {
00580 (*mpVizPdeSolutionResultsFile) << time << "\t";
00581
00582 #ifdef CHASTE_VTK
00583
00584 VtkMeshWriter<DIM,DIM>* p_vtk_mesh_writer = NULL;
00585 if (DIM>1 && mpCoarsePdeMesh != NULL )
00586 {
00587 std::ostringstream time_string;
00588 time_string << SimulationTime::Instance()->GetTimeStepsElapsed()+1;
00589 std::string results_file = "pde_results_"+time_string.str();
00590
00591 p_vtk_mesh_writer = new VtkMeshWriter<DIM,DIM>(mDirPath, results_file, false);
00592 }
00593 #endif //CHASTE_VTK
00594 for (unsigned pde_index=0; pde_index<mPdeAndBcCollection.size(); pde_index++)
00595 {
00596 if (mpCoarsePdeMesh != NULL)
00597 {
00598 PdeAndBoundaryConditions<DIM>* p_pde_and_bc = mPdeAndBcCollection[pde_index];
00599 assert(p_pde_and_bc->rGetDependentVariableName() != "");
00600
00601 #ifdef CHASTE_VTK
00602 if (p_pde_and_bc->GetSolution())
00603 {
00604 if (DIM>1)
00605 {
00606 ReplicatableVector solution_repl(p_pde_and_bc->GetSolution());
00607 std::vector<double> pde_solution;
00608 for (unsigned i=0; i<mpCoarsePdeMesh->GetNumNodes(); i++)
00609 {
00610 pde_solution.push_back(solution_repl[i]);
00611 }
00612
00613 p_vtk_mesh_writer->AddPointData(p_pde_and_bc->rGetDependentVariableName(),pde_solution);
00614 }
00615 }
00616
00617 #endif //CHASTE_VTK
00618
00619 for (unsigned i=0; i<mpCoarsePdeMesh->GetNumNodes(); i++)
00620 {
00621 (*mpVizPdeSolutionResultsFile) << i << " ";
00622 c_vector<double,DIM> location = mpCoarsePdeMesh->GetNode(i)->rGetLocation();
00623 for (unsigned k=0; k<DIM; k++)
00624 {
00625 (*mpVizPdeSolutionResultsFile) << location[k] << " ";
00626 }
00627
00628 if (p_pde_and_bc->GetSolution())
00629 {
00630 ReplicatableVector solution_repl(p_pde_and_bc->GetSolution());
00631 (*mpVizPdeSolutionResultsFile) << solution_repl[i] << " ";
00632 }
00633 else
00634 {
00635
00636 NEVER_REACHED;
00637 }
00638 }
00639 }
00640 else
00641 {
00642 for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
00643 cell_iter != mpCellPopulation->End();
00644 ++cell_iter)
00645 {
00646 unsigned node_index = mpCellPopulation->GetLocationIndexUsingCell(*cell_iter);
00647 (*mpVizPdeSolutionResultsFile) << node_index << " ";
00648 const c_vector<double,DIM>& position = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
00649 for (unsigned i=0; i<DIM; i++)
00650 {
00651 (*mpVizPdeSolutionResultsFile) << position[i] << " ";
00652 }
00653 double solution = cell_iter->GetCellData()->GetItem(mPdeAndBcCollection[pde_index]->rGetDependentVariableName());
00654 (*mpVizPdeSolutionResultsFile) << solution << " ";
00655 }
00656 }
00657 }
00658 (*mpVizPdeSolutionResultsFile) << "\n";
00659 #ifdef CHASTE_VTK
00660 if (p_vtk_mesh_writer != NULL)
00661 {
00662 p_vtk_mesh_writer->WriteFilesUsingMesh(*mpCoarsePdeMesh);
00663 delete p_vtk_mesh_writer;
00664 }
00665 #endif //CHASTE_VTK
00666 }
00667 }
00668
00669 template<unsigned DIM>
00670 void CellBasedPdeHandler<DIM>::SetWriteAverageRadialPdeSolution(const std::string& rName, unsigned numRadialIntervals, bool writeDailyResults)
00671 {
00672 mAverageRadialSolutionVariableName = rName;
00673 mWriteAverageRadialPdeSolution = true;
00674 mNumRadialIntervals = numRadialIntervals;
00675 mWriteDailyAverageRadialPdeSolution = writeDailyResults;
00676 }
00677
00678 template<unsigned DIM>
00679 void CellBasedPdeHandler<DIM>::SetImposeBcsOnCoarseBoundary(bool setBcsOnCoarseBoundary)
00680 {
00681 mSetBcsOnCoarseBoundary = setBcsOnCoarseBoundary;
00682 }
00683
00684 template<unsigned DIM>
00685 void CellBasedPdeHandler<DIM>::WriteAverageRadialPdeSolution(double time)
00686 {
00687 (*mpAverageRadialPdeSolutionResultsFile) << time << " ";
00688
00689
00690 c_vector<double,DIM> centre = mpCellPopulation->GetCentroidOfCellPopulation();
00691
00692
00693 std::map<double, CellPtr> radius_cell_map;
00694 double max_distance_from_centre = 0.0;
00695 for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
00696 cell_iter != mpCellPopulation->End();
00697 ++cell_iter)
00698 {
00699 double distance = norm_2(mpCellPopulation->GetLocationOfCellCentre(*cell_iter) - centre);
00700 radius_cell_map[distance] = *cell_iter;
00701
00702 if (distance > max_distance_from_centre)
00703 {
00704 max_distance_from_centre = distance;
00705 }
00706 }
00707
00708
00709 std::vector<double> radius_intervals;
00710 for (unsigned i=0; i<mNumRadialIntervals; i++)
00711 {
00712 double upper_radius = max_distance_from_centre*((double) i+1)/((double) mNumRadialIntervals);
00713 radius_intervals.push_back(upper_radius);
00714 }
00715
00716
00717 double lower_radius = 0.0;
00718 for (unsigned i=0; i<mNumRadialIntervals; i++)
00719 {
00720 unsigned counter = 0;
00721 double average_solution = 0.0;
00722
00723 for (std::map<double, CellPtr>::iterator iter = radius_cell_map.begin(); iter != radius_cell_map.end(); ++iter)
00724 {
00725 if (iter->first > lower_radius && iter->first <= radius_intervals[i])
00726 {
00727 average_solution += (iter->second)->GetCellData()->GetItem(mAverageRadialSolutionVariableName);
00728 counter++;
00729 }
00730 }
00731 if (counter != 0)
00732 {
00733 average_solution /= (double) counter;
00734 }
00735
00736
00737 (*mpAverageRadialPdeSolutionResultsFile) << radius_intervals[i] << " " << average_solution << " ";
00738 lower_radius = radius_intervals[i];
00739 }
00740 (*mpAverageRadialPdeSolutionResultsFile) << "\n";
00741 }
00742
00743 template<unsigned DIM>
00744 double CellBasedPdeHandler<DIM>::GetPdeSolutionAtPoint(const c_vector<double,DIM>& rPoint, const std::string& rVariable)
00745 {
00746 double solution_at_point = 0.0;
00747
00748 unsigned pde_index = UINT_MAX;
00749
00750
00751 for (unsigned i=0; i<mPdeAndBcCollection.size(); i++)
00752 {
00753 if (mPdeAndBcCollection[i]->rGetDependentVariableName() == rVariable)
00754 {
00755 pde_index = i;
00756 break;
00757 }
00758 }
00759 if (pde_index == UINT_MAX)
00760 {
00761 EXCEPTION("Tried to get the solution of a variable name: " + rVariable + ". There is no PDE with that variable.");
00762 }
00763 PdeAndBoundaryConditions<DIM>* p_pde_and_bc = mPdeAndBcCollection[pde_index];
00764
00765 Element<DIM,DIM>* p_containing_element;
00766
00767 if (mpCoarsePdeMesh != NULL)
00768 {
00769
00770 unsigned elem_index = mpCoarsePdeMesh->GetContainingElementIndex(ChastePoint<DIM>(rPoint));
00771 p_containing_element = mpCoarsePdeMesh->GetElement(elem_index);
00772 }
00773 else
00774 {
00775
00776 TetrahedralMesh<DIM,DIM>* p_tetrahedral_mesh = &(static_cast<MeshBasedCellPopulation<DIM>*>(mpCellPopulation)->rGetMesh());
00777
00778 unsigned elem_index = p_tetrahedral_mesh->GetContainingElementIndex(ChastePoint<DIM>(rPoint));
00779 p_containing_element = p_tetrahedral_mesh->GetElement(elem_index);
00780 }
00781
00782
00783 if (p_pde_and_bc->GetSolution())
00784 {
00785 ReplicatableVector solution_repl(p_pde_and_bc->GetSolution());
00786 c_vector<double,DIM+1> weights = p_containing_element->CalculateInterpolationWeights(rPoint);
00787 for (unsigned i=0; i<DIM+1; i++)
00788 {
00789 double nodal_value = solution_repl[p_containing_element->GetNodeGlobalIndex(i)];
00790 solution_at_point += nodal_value * weights(i);
00791 }
00792 }
00793 else
00794 {
00795
00796 NEVER_REACHED;
00797 }
00798
00799 return solution_at_point;
00800 }
00801
00802 template<unsigned DIM>
00803 bool CellBasedPdeHandler<DIM>::GetWriteAverageRadialPdeSolution()
00804 {
00805 return mWriteAverageRadialPdeSolution;
00806 }
00807
00808 template<unsigned DIM>
00809 bool CellBasedPdeHandler<DIM>::GetWriteDailyAverageRadialPdeSolution()
00810 {
00811 return mWriteDailyAverageRadialPdeSolution;
00812 }
00813
00814 template<unsigned DIM>
00815 bool CellBasedPdeHandler<DIM>::GetImposeBcsOnCoarseBoundary()
00816 {
00817 return mSetBcsOnCoarseBoundary;
00818 }
00819
00820 template<unsigned DIM>
00821 unsigned CellBasedPdeHandler<DIM>::GetNumRadialIntervals()
00822 {
00823 return mNumRadialIntervals;
00824 }
00825
00826 template<unsigned DIM>
00827 void CellBasedPdeHandler<DIM>::OutputParameters(out_stream& rParamsFile)
00828 {
00829 std::string type = GetIdentifier();
00830
00831 *rParamsFile << "\t\t<" << type << ">\n";
00832 *rParamsFile << "\t\t<WriteAverageRadialPdeSolution>" << mWriteAverageRadialPdeSolution << "</WriteAverageRadialPdeSolution>\n";
00833 *rParamsFile << "\t\t<WriteDailyAverageRadialPdeSolution>" << mWriteDailyAverageRadialPdeSolution << "</WriteDailyAverageRadialPdeSolution>\n";
00834 *rParamsFile << "\t\t<SetBcsOnCoarseBoundary>" << mSetBcsOnCoarseBoundary << "</SetBcsOnCoarseBoundary>\n";
00835 *rParamsFile << "\t\t<NumRadialIntervals>" << mNumRadialIntervals << "</NumRadialIntervals>\n";
00836 *rParamsFile << "\t\t</" << type << ">\n";
00837 }
00838
00839 template<unsigned DIM>
00840 bool CellBasedPdeHandler<DIM>::PdeSolveNeedsCoarseMesh()
00841 {
00842 return ((dynamic_cast<NodeBasedCellPopulation<DIM>*>(mpCellPopulation) != NULL)
00843 || (dynamic_cast<PottsBasedCellPopulation<DIM>*>(mpCellPopulation) != NULL)
00844 || (dynamic_cast<CaBasedCellPopulation<DIM>*>(mpCellPopulation) != NULL)
00845 || (dynamic_cast<VertexBasedCellPopulation<DIM>*>(mpCellPopulation) != NULL));
00846 }
00847
00848
00849 #include "SerializationExportWrapperForCpp.hpp"
00850 EXPORT_TEMPLATE_CLASS_SAME_DIMS(CellBasedPdeHandler)
00851
00852
00853
00855
00856 template class CellBasedPdeHandler<1>;
00857 template class CellBasedPdeHandler<2>;
00858 template class CellBasedPdeHandler<3>;