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 #include "CardiacElectroMechanicsProblem.hpp"
00030
00031 #include "OutputFileHandler.hpp"
00032 #include "ReplicatableVector.hpp"
00033 #include "HeartConfig.hpp"
00034 #include "LogFile.hpp"
00035 #include "ChastePoint.hpp"
00036 #include "Element.hpp"
00037 #include "BoundaryConditionsContainer.hpp"
00038 #include "AbstractDynamicLinearPdeSolver.hpp"
00039 #include "TimeStepper.hpp"
00040 #include "QuadraturePointsGroup.hpp"
00041 #include "TrianglesMeshWriter.hpp"
00042 #include "Hdf5ToMeshalyzerConverter.hpp"
00043 #include "Hdf5ToCmguiConverter.hpp"
00044 #include "MeshalyzerMeshWriter.hpp"
00045 #include "PetscTools.hpp"
00046 #include "NashHunterPoleZeroLaw.hpp"
00047 #include "ImplicitCardiacMechanicsSolver.hpp"
00048 #include "ExplicitCardiacMechanicsSolver.hpp"
00049 #include "MooneyRivlinMaterialLaw.hpp"
00050 #include "CmguiDeformedSolutionsWriter.hpp"
00051 #include "VoltageInterpolaterOntoMechanicsMesh.hpp"
00052
00053
00054 template<unsigned DIM>
00055 void CardiacElectroMechanicsProblem<DIM>::DetermineWatchedNodes()
00056 {
00057 assert(mIsWatchedLocation);
00058
00059
00060 double min_dist = DBL_MAX;
00061 unsigned node_index = UNSIGNED_UNSET;
00062 for(unsigned i=0; i<mpElectricsMesh->GetNumNodes(); i++)
00063 {
00064 double dist = norm_2(mWatchedLocation - mpElectricsMesh->GetNode(i)->rGetLocation());
00065 if(dist < min_dist)
00066 {
00067 min_dist = dist;
00068 node_index = i;
00069 }
00070 }
00071
00072
00073 assert(node_index != UNSIGNED_UNSET);
00074 c_vector<double,DIM> pos = mpElectricsMesh->GetNode(node_index)->rGetLocation();
00075
00076 if(min_dist > 1e-8)
00077 {
00078 #define COVERAGE_IGNORE
00079 std::cout << "ERROR: Could not find an electrics node very close to requested watched location - "
00080 << "min distance was " << min_dist << " for node " << node_index
00081 << " at location " << pos << std::flush;;
00082
00084
00085 NEVER_REACHED;
00086 #undef COVERAGE_IGNORE
00087 }
00088 else
00089 {
00090 LOG(2,"Chosen electrics node "<<node_index<<" at location " << pos << " to be watched");
00091 mWatchedElectricsNodeIndex = node_index;
00092 }
00093
00094
00095 min_dist = DBL_MAX;
00096 node_index = UNSIGNED_UNSET;
00097 c_vector<double,DIM> pos_at_min;
00098
00099 for(unsigned i=0; i<mpMechanicsMesh->GetNumNodes(); i++)
00100 {
00101 c_vector<double,DIM> position = mpMechanicsMesh->GetNode(i)->rGetLocation();
00102
00103 double dist = norm_2(position-mWatchedLocation);
00104
00105 if(dist < min_dist)
00106 {
00107 min_dist = dist;
00108 node_index = i;
00109 pos_at_min = position;
00110 }
00111 }
00112
00113
00114 assert(node_index != UNSIGNED_UNSET);
00115
00116 if(min_dist > 1e-8)
00117 {
00118 #define COVERAGE_IGNORE
00119 std::cout << "ERROR: Could not find a mechanics node very close to requested watched location - "
00120 << "min distance was " << min_dist << " for node " << node_index
00121 << " at location " << pos_at_min;
00122
00124
00125 assert(0);
00126 #undef COVERAGE_IGNORE
00127 }
00128 else
00129 {
00130 LOG(2,"Chosen mechanics node "<<node_index<<" at location " << pos << " to be watched");
00131 mWatchedMechanicsNodeIndex = node_index;
00132 }
00133
00134 OutputFileHandler handler(mOutputDirectory,false);
00135 mpWatchedLocationFile = handler.OpenOutputFile("watched.txt");
00136 }
00137
00138 template<unsigned DIM>
00139 void CardiacElectroMechanicsProblem<DIM>::WriteWatchedLocationData(double time, Vec voltage)
00140 {
00141 assert(mIsWatchedLocation);
00142
00143 std::vector<c_vector<double,DIM> >& deformed_position = mpMechanicsSolver->rGetDeformedPosition();
00144
00146 ReplicatableVector voltage_replicated(voltage);
00147 double V = voltage_replicated[mWatchedElectricsNodeIndex];
00148
00151
00152
00153
00154 *mpWatchedLocationFile << time << " ";
00155 for(unsigned i=0; i<DIM; i++)
00156 {
00157 *mpWatchedLocationFile << deformed_position[mWatchedMechanicsNodeIndex](i) << " ";
00158 }
00159 *mpWatchedLocationFile << V << "\n";
00160 mpWatchedLocationFile->flush();
00161 }
00162
00163 template<unsigned DIM>
00164 c_matrix<double,DIM,DIM>& CardiacElectroMechanicsProblem<DIM>::rGetModifiedConductivityTensor(unsigned elementIndex, const c_matrix<double,DIM,DIM>& rOriginalConductivity)
00165 {
00166 if(mLastModifiedConductivity.first==elementIndex)
00167 {
00168
00169 return mLastModifiedConductivity.second;
00170 }
00171 else
00172 {
00173
00174
00175
00176 unsigned containing_mechanics_elem = mpMeshPair->rGetCoarseElementsForFineElementCentroids()[elementIndex];
00177 c_matrix<double,DIM,DIM>& r_deformation_gradient = mDeformationGradientsForEachMechanicsElement[containing_mechanics_elem];
00178
00179
00180 c_matrix<double,DIM,DIM> inv_F = Inverse(r_deformation_gradient);
00181 mLastModifiedConductivity.second = prod(inv_F, rOriginalConductivity);
00182 mLastModifiedConductivity.second = prod(mLastModifiedConductivity.second, trans(inv_F));
00183
00184
00185 mLastModifiedConductivity.first = elementIndex;
00186 return mLastModifiedConductivity.second;
00187 }
00188 }
00189
00190
00191
00192
00194
00195
00196
00197
00198
00199
00200
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219
00220
00221
00222
00223
00224
00225 template<unsigned DIM>
00226 CardiacElectroMechanicsProblem<DIM>::CardiacElectroMechanicsProblem(
00227 ContractionModel contractionModel,
00228 TetrahedralMesh<DIM,DIM>* pElectricsMesh,
00229 QuadraticMesh<DIM>* pMechanicsMesh,
00230 std::vector<unsigned> fixedMechanicsNodes,
00231 AbstractCardiacCellFactory<DIM>* pCellFactory,
00232 double endTime,
00233 double electricsPdeTimeStep,
00234 double mechanicsSolveTimestep,
00235 double contractionModelOdeTimeStep,
00236 std::string outputDirectory = "") :
00237 mpMeshPair(NULL)
00238 {
00239
00240 MechanicsEventHandler::Reset();
00241 MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ALL);
00242
00243
00244
00245
00246 HeartEventHandler::Disable();
00247
00248 mContractionModel = contractionModel;
00249
00250
00251 assert(endTime > 0);
00252 mEndTime = endTime;
00253
00254
00255 HeartConfig::Instance()->SetSimulationDuration(mEndTime);
00256
00257 assert(electricsPdeTimeStep>0);
00258 mElectricsTimeStep = electricsPdeTimeStep;
00259
00260 mNumElecTimestepsPerMechTimestep = (unsigned) floor((mechanicsSolveTimestep/electricsPdeTimeStep)+0.5);
00261 if(fabs(mNumElecTimestepsPerMechTimestep*electricsPdeTimeStep - mechanicsSolveTimestep) > 1e-6)
00262 {
00263 EXCEPTION("Electrics PDE timestep does not divide mechanics solve timestep");
00264 }
00265 mMechanicsTimeStep = mechanicsSolveTimestep;
00266
00267
00268
00269 assert(mechanicsSolveTimestep <= 10);
00270
00271 assert(contractionModelOdeTimeStep <= mMechanicsTimeStep+1e-14);
00272 mContractionModelOdeTimeStep = contractionModelOdeTimeStep;
00273
00274
00275
00276
00277 assert(pCellFactory != NULL);
00278 mpMonodomainProblem = new MonodomainProblem<DIM>(pCellFactory);
00279
00280
00281 mWriteOutput = (outputDirectory!="");
00282 if(mWriteOutput)
00283 {
00284 mOutputDirectory = outputDirectory;
00285
00286 OutputFileHandler handler(mOutputDirectory);
00287 mDeformationOutputDirectory = mOutputDirectory + "/deformation";
00288 HeartConfig::Instance()->SetOutputDirectory(mOutputDirectory + "/electrics");
00289 HeartConfig::Instance()->SetOutputFilenamePrefix("voltage");
00290 }
00291 else
00292 {
00293 mDeformationOutputDirectory = "";
00294 }
00295 mNoElectricsOutput = false;
00296
00297
00298 mpElectricsMesh = pElectricsMesh;
00299 mpMechanicsMesh = pMechanicsMesh;
00300 mFixedNodes = fixedMechanicsNodes;
00301
00302 mpCardiacMechSolver = NULL;
00303 mpMechanicsSolver = NULL;
00304
00305 mpMaterialLaw = NULL;
00306 mAllocatedMaterialLawMemory = false;
00307
00308
00309
00310
00311 std::string log_dir = mOutputDirectory;
00312 LogFile::Instance()->Set(2, mOutputDirectory);
00313 LogFile::Instance()->WriteHeader("Electromechanics");
00314 LOG(2, DIM << "d Implicit CardiacElectroMechanics Simulation:");
00315 LOG(2, "End time = " << mEndTime << ", electrics time step = " << mElectricsTimeStep << ", mechanics timestep = " << mMechanicsTimeStep << "\n");
00316 LOG(2, "Contraction model ode timestep " << mContractionModelOdeTimeStep);
00317 LOG(2, "Output is written to " << mOutputDirectory << "/[deformation/electrics]");
00318
00319 if(mpElectricsMesh != NULL)
00320 {
00321 LOG(2, "Electrics mesh has " << mpElectricsMesh->GetNumNodes() << " nodes");
00322 }
00323 if(mpMechanicsMesh != NULL)
00324 {
00325 LOG(2, "Mechanics mesh has " << mpMechanicsMesh->GetNumNodes() << " nodes");
00326 }
00327
00328 mIsWatchedLocation = false;
00329 mWatchedElectricsNodeIndex = UNSIGNED_UNSET;
00330 mWatchedMechanicsNodeIndex = UNSIGNED_UNSET;
00331
00332 mFibreSheetDirectionsFile = "";
00333
00334 mConductivityAffectedByDeformationMef = false;
00335 mCellModelsAffectedByDeformationMef = false;
00336
00337 mLastModifiedConductivity.first = UNSIGNED_UNSET;
00338
00339
00340 }
00341
00342 template<unsigned DIM>
00343 CardiacElectroMechanicsProblem<DIM>::~CardiacElectroMechanicsProblem()
00344 {
00350 if(mIsWatchedLocation && mpMechanicsMesh)
00351 {
00352 mpWatchedLocationFile->close();
00353 }
00354
00355 delete mpMonodomainProblem;
00356
00357 delete mpCardiacMechSolver;
00358
00359 delete mpMeshPair;
00360
00361 if(mAllocatedMaterialLawMemory && mpMaterialLaw)
00362 {
00363 delete mpMaterialLaw;
00364 }
00365
00366 LogFile::Close();
00367 }
00368
00369 template<unsigned DIM>
00370 void CardiacElectroMechanicsProblem<DIM>::Initialise()
00371 {
00372 LOG(2, "Initialising..");
00373
00374 assert(mpElectricsMesh!=NULL);
00375 assert(mpMechanicsMesh!=NULL);
00376 assert(mpCardiacMechSolver==NULL);
00377
00378 if(mIsWatchedLocation)
00379 {
00380 DetermineWatchedNodes();
00381 }
00382
00383
00384
00385 mpMonodomainProblem->SetMesh(mpElectricsMesh);
00386 HeartConfig::Instance()->SetIntracellularConductivities(Create_c_vector(1.75,1.75,1.75));
00387 mpMonodomainProblem->Initialise();
00388
00389
00390 if(mpMaterialLaw == NULL)
00391 {
00392 mpMaterialLaw = new NashHunterPoleZeroLaw<DIM>();
00393 mAllocatedMaterialLawMemory = true;
00394 }
00395
00396
00397
00398
00399 switch(mContractionModel)
00400 {
00401 case NASH2004:
00402
00403 mpCardiacMechSolver = new ExplicitCardiacMechanicsSolver<NonlinearElasticitySolver<DIM>,DIM>(
00404 mContractionModel,mpMechanicsMesh,mDeformationOutputDirectory,mFixedNodes,mpMaterialLaw);
00405 break;
00406 case KERCHOFFS2003:
00407
00408 mpCardiacMechSolver = new ImplicitCardiacMechanicsSolver<NonlinearElasticitySolver<DIM>,DIM>(
00409 mContractionModel,mpMechanicsMesh,mDeformationOutputDirectory,mFixedNodes,mpMaterialLaw);
00410 break;
00411 case NHS:
00412
00413 mpCardiacMechSolver = new ImplicitCardiacMechanicsSolver<NonlinearElasticitySolver<DIM>,DIM>(
00414 mContractionModel,mpMechanicsMesh,mDeformationOutputDirectory,mFixedNodes,mpMaterialLaw);
00415 break;
00416 default:
00417 EXCEPTION("Invalid contraction model, options are: KERCHOFFS2003 or NHS");
00418 }
00419
00420 mpMechanicsSolver = dynamic_cast<AbstractNonlinearElasticitySolver<DIM>*>(mpCardiacMechSolver);
00421 assert(mpMechanicsSolver);
00422
00423 if(mFibreSheetDirectionsFile!="")
00424 {
00425 mpCardiacMechSolver->SetVariableFibreSheetDirections(mFibreSheetDirectionsFile, mFibreSheetDirectionsDefinedPerQuadraturePoint);
00426 }
00427
00428
00429
00430 mpMeshPair = new FineCoarseMeshPair<DIM>(*mpElectricsMesh, *mpMechanicsMesh);
00431 mpMeshPair->SetUpBoxesOnFineMesh();
00432 mpMeshPair->ComputeFineElementsAndWeightsForCoarseQuadPoints(*(mpCardiacMechSolver->GetQuadratureRule()), false);
00433 mpMeshPair->DeleteFineBoxCollection();
00434
00435
00436
00437
00438 if(mConductivityAffectedByDeformationMef || mCellModelsAffectedByDeformationMef)
00439 {
00440 mpMeshPair->SetUpBoxesOnCoarseMesh();
00441 }
00442
00443 if(mCellModelsAffectedByDeformationMef)
00444 {
00445
00446
00447 mpMeshPair->ComputeCoarseElementsForFineNodes(false);
00448
00449
00450 mStretchesForEachMechanicsElement.resize(mpMechanicsMesh->GetNumElements(),1.0);
00451 }
00452
00453 if(mConductivityAffectedByDeformationMef)
00454 {
00455
00456
00457 mpMeshPair->ComputeCoarseElementsForFineElementCentroids(false);
00458
00459
00460 mDeformationGradientsForEachMechanicsElement.resize(mpMechanicsMesh->GetNumElements(),identity_matrix<double>(DIM));
00461
00462
00463
00464 mpMonodomainProblem->GetMonodomainTissue()->SetConductivityModifier(this);
00465 }
00466
00467 if(mWriteOutput)
00468 {
00469 TrianglesMeshWriter<DIM,DIM> mesh_writer(mOutputDirectory,"electrics_mesh",false);
00470 mesh_writer.WriteFilesUsingMesh(*mpElectricsMesh);
00471 }
00472 }
00473
00474 template<unsigned DIM>
00475 void CardiacElectroMechanicsProblem<DIM>::Solve()
00476 {
00477
00478
00479 if(mpCardiacMechSolver==NULL)
00480 {
00481 Initialise();
00482 }
00483
00484 boost::shared_ptr<BoundaryConditionsContainer<DIM,DIM,1> > p_bcc(new BoundaryConditionsContainer<DIM,DIM,1>);
00485 p_bcc->DefineZeroNeumannOnMeshBoundary(mpElectricsMesh, 0);
00486 mpMonodomainProblem->SetBoundaryConditionsContainer(p_bcc);
00487
00488
00489
00490 AbstractDynamicLinearPdeSolver<DIM,DIM,1>* p_electrics_solver
00491 = mpMonodomainProblem->CreateSolver();
00492
00493
00494 Vec voltage=NULL;
00495 Vec initial_voltage = mpMonodomainProblem->CreateInitialCondition();
00496
00497 unsigned num_quad_points = mpCardiacMechSolver->GetTotalNumQuadPoints();
00498 std::vector<double> interpolated_calcium_concs(num_quad_points, 0.0);
00499 std::vector<double> interpolated_voltages(num_quad_points, 0.0);
00500
00501
00502 unsigned counter = 0;
00503
00504 TimeStepper stepper(0.0, mEndTime, mMechanicsTimeStep);
00505
00506 CmguiDeformedSolutionsWriter<DIM>* p_cmgui_writer = NULL;
00507
00508 unsigned mech_writer_counter = 0;
00509 if (mWriteOutput)
00510 {
00511 mpMechanicsSolver->SetWriteOutput();
00512 mpMechanicsSolver->WriteCurrentDeformation("solution",mech_writer_counter);
00513
00514 p_cmgui_writer = new CmguiDeformedSolutionsWriter<DIM>(mOutputDirectory+"/deformation/cmgui",
00515 "solution",
00516 *(this->mpMechanicsMesh),
00517 WRITE_QUADRATIC_MESH);
00518 std::vector<std::string> fields;
00519 fields.push_back("V");
00520 p_cmgui_writer->SetAdditionalFieldNames(fields);
00521 p_cmgui_writer->WriteInitialMesh();
00522
00523
00524 if(!mNoElectricsOutput)
00525 {
00526
00527
00528
00529 HeartConfig::Instance()->SetPrintingTimeStep(mMechanicsTimeStep);
00530 mpMonodomainProblem->InitialiseWriter();
00531 mpMonodomainProblem->WriteOneStep(stepper.GetTime(), initial_voltage);
00532 }
00533
00534 if(mIsWatchedLocation)
00535 {
00536 WriteWatchedLocationData(stepper.GetTime(), initial_voltage);
00537 }
00538 }
00539
00540 PrepareForSolve();
00541
00542 while (!stepper.IsTimeAtEnd())
00543 {
00544 LOG(2, "\nCurrent time = " << stepper.GetTime());
00545 #ifdef MECH_VERBOSE // defined in AbstractNonlinearElasticitySolver
00546
00547 std::cout << "\n\n ** Current time = " << stepper.GetTime() << "\n";
00548 #endif
00549
00556 if(mCellModelsAffectedByDeformationMef)
00557 {
00558
00559
00560 mpCardiacMechSolver->ComputeDeformationGradientAndStretchInEachElement(mDeformationGradientsForEachMechanicsElement, mStretchesForEachMechanicsElement);
00561
00562
00563 for(unsigned global_index = mpElectricsMesh->GetDistributedVectorFactory()->GetLow();
00564 global_index < mpElectricsMesh->GetDistributedVectorFactory()->GetHigh();
00565 global_index++)
00566 {
00567 unsigned containing_elem = mpMeshPair->rGetCoarseElementsForFineNodes()[global_index];
00568 double stretch = mStretchesForEachMechanicsElement[containing_elem];
00569 mpMonodomainProblem->GetTissue()->GetCardiacCell(global_index)->SetStretch(stretch);
00570 }
00571
00572
00573
00574
00575 }
00576
00577 p_electrics_solver->SetTimeStep(mElectricsTimeStep);
00578
00584 LOG(2, " Solving electrics");
00585 MechanicsEventHandler::BeginEvent(MechanicsEventHandler::NON_MECH);
00586 for(unsigned i=0; i<mNumElecTimestepsPerMechTimestep; i++)
00587 {
00588 double current_time = stepper.GetTime() + i*mElectricsTimeStep;
00589 double next_time = stepper.GetTime() + (i+1)*mElectricsTimeStep;
00590
00591
00592 p_electrics_solver->SetTimes(current_time, next_time);
00593 p_electrics_solver->SetInitialCondition( initial_voltage );
00594
00595 voltage = p_electrics_solver->Solve();
00596
00597 PetscReal min_voltage, max_voltage;
00598 VecMax(voltage,PETSC_NULL,&max_voltage);
00599 VecMin(voltage,PETSC_NULL,&min_voltage);
00600 if(i==0)
00601 {
00602 LOG(2, " minimum and maximum voltage is " << min_voltage <<", "<<max_voltage);
00603 }
00604 else if(i==1)
00605 {
00606 LOG(2, " ..");
00607 }
00608
00609 VecDestroy(initial_voltage);
00610 initial_voltage = voltage;
00611 }
00612
00613 if(mConductivityAffectedByDeformationMef)
00614 {
00615 p_electrics_solver->SetMatrixIsNotAssembled();
00616 }
00617
00618
00624
00625
00626
00627 LOG(2, " Interpolating Ca_I and voltage");
00628
00629 ReplicatableVector voltage_repl(voltage);
00630
00631
00632
00633 ReplicatableVector calcium_repl(mpElectricsMesh->GetNumNodes());
00634 PetscInt lo, hi;
00635 VecGetOwnershipRange(voltage, &lo, &hi);
00636 for(int index=lo; index<hi; index++)
00637 {
00638 calcium_repl[index] = mpMonodomainProblem->GetTissue()->GetCardiacCell(index)->GetIntracellularCalciumConcentration();
00639 }
00640 PetscTools::Barrier();
00641 calcium_repl.Replicate(lo,hi);
00642
00643
00644 for(unsigned i=0; i<mpMeshPair->rGetElementsAndWeights().size(); i++)
00645 {
00646 double interpolated_CaI = 0;
00647 double interpolated_voltage = 0;
00648
00649 Element<DIM,DIM>& element = *(mpElectricsMesh->GetElement(mpMeshPair->rGetElementsAndWeights()[i].ElementNum));
00650 for(unsigned node_index = 0; node_index<element.GetNumNodes(); node_index++)
00651 {
00652 unsigned global_node_index = element.GetNodeGlobalIndex(node_index);
00653 double CaI_at_node = calcium_repl[global_node_index];
00654 interpolated_CaI += CaI_at_node*mpMeshPair->rGetElementsAndWeights()[i].Weights(node_index);
00655 interpolated_voltage += voltage_repl[global_node_index]*mpMeshPair->rGetElementsAndWeights()[i].Weights(node_index);
00656 }
00657
00658 interpolated_calcium_concs[i] = interpolated_CaI;
00659 interpolated_voltages[i] = interpolated_voltage;
00660 }
00661
00662 LOG(2, " Setting Ca_I. max value = " << Max(interpolated_calcium_concs));
00663
00664
00665
00666
00667
00668 mpCardiacMechSolver->SetCalciumAndVoltage(interpolated_calcium_concs, interpolated_voltages);
00669 MechanicsEventHandler::EndEvent(MechanicsEventHandler::NON_MECH);
00670
00671
00672
00678 LOG(2, " Solving mechanics ");
00679 mpMechanicsSolver->SetWriteOutput(false);
00680
00681
00682
00683 mpMechanicsSolver->SetCurrentTime(stepper.GetTime());
00684
00686
00687
00688 MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ALL_MECH);
00689 mpCardiacMechSolver->Solve(stepper.GetTime(), stepper.GetNextTime(), mContractionModelOdeTimeStep);
00690 MechanicsEventHandler::EndEvent(MechanicsEventHandler::ALL_MECH);
00691
00692 LOG(2, " Number of newton iterations = " << mpMechanicsSolver->GetNumNewtonIterations());
00693
00694
00695 stepper.AdvanceOneTimeStep();
00696 counter++;
00697
00698
00704 MechanicsEventHandler::BeginEvent(MechanicsEventHandler::OUTPUT);
00705 if(mWriteOutput && (counter%WRITE_EVERY_NTH_TIME==0))
00706 {
00707 LOG(2, " Writing output");
00708
00709 mech_writer_counter++;
00710 mpMechanicsSolver->SetWriteOutput();
00711 mpMechanicsSolver->WriteCurrentDeformation("solution",mech_writer_counter);
00712
00713 p_cmgui_writer->WriteDeformationPositions(rGetDeformedPosition(), counter);
00714
00715 if(!mNoElectricsOutput)
00716 {
00717 mpMonodomainProblem->mpWriter->AdvanceAlongUnlimitedDimension();
00718 mpMonodomainProblem->WriteOneStep(stepper.GetTime(), voltage);
00719 }
00720
00721 if(mIsWatchedLocation)
00722 {
00723 WriteWatchedLocationData(stepper.GetTime(), voltage);
00724 }
00725 OnEndOfTimeStep(counter);
00726 }
00727 MechanicsEventHandler::EndEvent(MechanicsEventHandler::OUTPUT);
00728
00729
00730 LogFile::Instance()->WriteElapsedTime(" ");
00731 }
00732
00733
00734
00735 if ((mWriteOutput) && (!mNoElectricsOutput))
00736 {
00737 HeartConfig::Instance()->Reset();
00738 mpMonodomainProblem->mpWriter->Close();
00739 delete mpMonodomainProblem->mpWriter;
00740
00741 std::string input_dir = mOutputDirectory+"/electrics";
00742
00743
00744 std::string config_directory = HeartConfig::Instance()->GetOutputDirectory();
00745 HeartConfig::Instance()->SetOutputDirectory(input_dir);
00746
00747
00748
00749
00750 Hdf5ToCmguiConverter<DIM,DIM> cmgui_converter(input_dir,"voltage",mpElectricsMesh);
00751
00752
00753
00754
00755
00756
00757
00758
00759
00760
00761
00762
00763
00764 VoltageInterpolaterOntoMechanicsMesh<DIM> cnverter(*mpElectricsMesh,*mpMechanicsMesh,input_dir,"voltage");
00765
00766
00767 HeartConfig::Instance()->SetOutputDirectory(config_directory);
00768 }
00769
00770 if(p_cmgui_writer)
00771 {
00772 if(mNoElectricsOutput)
00773 {
00774 p_cmgui_writer->WriteCmguiScript();
00775 }
00776 else
00777 {
00778 p_cmgui_writer->WriteCmguiScript("../../electrics/cmgui_output/voltage_mechanics_mesh");
00779 }
00780 delete p_cmgui_writer;
00781 }
00782 VecDestroy(voltage);
00783 delete p_electrics_solver;
00784
00785 MechanicsEventHandler::EndEvent(MechanicsEventHandler::ALL);
00786 }
00787
00788
00789
00790 template<unsigned DIM>
00791 double CardiacElectroMechanicsProblem<DIM>::Max(std::vector<double>& vec)
00792 {
00793 double max = -1e200;
00794 for(unsigned i=0; i<vec.size(); i++)
00795 {
00796 if(vec[i]>max) max=vec[i];
00797 }
00798 return max;
00799 }
00800
00801 template<unsigned DIM>
00802 void CardiacElectroMechanicsProblem<DIM>::SetNoElectricsOutput()
00803 {
00804 mNoElectricsOutput = true;
00805 }
00806
00807 template<unsigned DIM>
00808 void CardiacElectroMechanicsProblem<DIM>::SetWatchedPosition(c_vector<double,DIM> watchedLocation)
00809 {
00810 mIsWatchedLocation = true;
00811 mWatchedLocation = watchedLocation;
00812 }
00813
00814 template<unsigned DIM>
00815 void CardiacElectroMechanicsProblem<DIM>::SetVariableFibreSheetDirectionsFile(std::string fibreSheetDirectionsFile, bool definedPerQuadraturePoint)
00816 {
00817 mFibreSheetDirectionsFile = fibreSheetDirectionsFile;
00818 mFibreSheetDirectionsDefinedPerQuadraturePoint = definedPerQuadraturePoint;
00819 }
00820
00821 template<unsigned DIM>
00822 std::vector<c_vector<double,DIM> >& CardiacElectroMechanicsProblem<DIM>::rGetDeformedPosition()
00823 {
00824 return mpMechanicsSolver->rGetDeformedPosition();
00825 }
00826
00827
00828
00829
00831
00833
00834
00835 template class CardiacElectroMechanicsProblem<2>;
00836 template class CardiacElectroMechanicsProblem<3>;