CardiacElectroMechanicsProblem.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2010
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
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 
00047 #include "ImplicitCardiacMechanicsSolver.hpp"
00048 #include "ExplicitCardiacMechanicsSolver.hpp"
00049 // if including Cinv in monobidomain equations
00050 //#include "NodewiseData.hpp"
00051 
00052 #include "MooneyRivlinMaterialLaw.hpp"
00053 #include "CmguiDeformedSolutionsWriter.hpp"
00054 #include "VoltageInterpolaterOntoMechanicsMesh.hpp"
00055 
00056 template<unsigned DIM>
00057 void CardiacElectroMechanicsProblem<DIM>::DetermineWatchedNodes()
00058 {
00059     assert(mIsWatchedLocation);
00060 
00061     // find the nearest electrics mesh node
00062     double min_dist = DBL_MAX;
00063     unsigned node_index = UNSIGNED_UNSET;
00064     for(unsigned i=0; i<mpElectricsMesh->GetNumNodes(); i++)
00065     {
00066         double dist = norm_2(mWatchedLocation - mpElectricsMesh->GetNode(i)->rGetLocation());
00067         if(dist < min_dist)
00068         {
00069             min_dist = dist;
00070             node_index = i;
00071         }
00072     }
00073 
00074     // set up watched node, if close enough
00075     assert(node_index != UNSIGNED_UNSET); // should def have found something
00076     c_vector<double,DIM> pos = mpElectricsMesh->GetNode(node_index)->rGetLocation();
00077 
00078     if(min_dist > 1e-8)
00079     {
00080         #define COVERAGE_IGNORE
00081         std::cout << "ERROR: Could not find an electrics node very close to requested watched location - "
00082                   << "min distance was " << min_dist << " for node " << node_index
00083                   << " at location " << pos << std::flush;;
00084 
00086         //EXCEPTION("Could not find an electrics node very close to requested watched location");
00087         NEVER_REACHED;
00088         #undef COVERAGE_IGNORE
00089     }
00090     else
00091     {
00092         LOG(2,"Chosen electrics node "<<node_index<<" at location " << pos << " to be watched");
00093         mWatchedElectricsNodeIndex = node_index;
00094     }
00095 
00096     // find nearest mechanics mesh
00097     min_dist = DBL_MAX;
00098     node_index = UNSIGNED_UNSET;
00099     c_vector<double,DIM> pos_at_min;
00100 
00101     for(unsigned i=0; i<mpMechanicsMesh->GetNumNodes(); i++)
00102     {
00103         c_vector<double,DIM> position = mpMechanicsMesh->GetNode(i)->rGetLocation();
00104 
00105         double dist = norm_2(position-mWatchedLocation);
00106 
00107         if(dist < min_dist)
00108         {
00109             min_dist = dist;
00110             node_index = i;
00111             pos_at_min = position;
00112         }
00113     }
00114 
00115     // set up watched node, if close enough
00116     assert(node_index != UNSIGNED_UNSET); // should def have found something
00117 
00118     if(min_dist > 1e-8)
00119     {
00120         #define COVERAGE_IGNORE
00121         std::cout << "ERROR: Could not find a mechanics node very close to requested watched location - "
00122                   << "min distance was " << min_dist << " for node " << node_index
00123                   << " at location " << pos_at_min;
00124 
00126         //EXCEPTION("Could not find a mechanics node very close to requested watched location");
00127         assert(0);
00128         #undef COVERAGE_IGNORE
00129     }
00130     else
00131     {
00132         LOG(2,"Chosen mechanics node "<<node_index<<" at location " << pos << " to be watched");
00133         mWatchedMechanicsNodeIndex = node_index;
00134     }
00135 
00136     OutputFileHandler handler(mOutputDirectory,false);
00137     mpWatchedLocationFile = handler.OpenOutputFile("watched.txt");
00138 }
00139 
00140 template<unsigned DIM>
00141 void CardiacElectroMechanicsProblem<DIM>::WriteWatchedLocationData(double time, Vec voltage)
00142 {
00143     assert(mIsWatchedLocation);
00144 
00145     std::vector<c_vector<double,DIM> >& deformed_position = mpCardiacMechSolver->rGetDeformedPosition();
00146 
00148     ReplicatableVector voltage_replicated(voltage);
00149     double V = voltage_replicated[mWatchedElectricsNodeIndex];
00150 
00153     //     // Metadata is currently being added to CellML models and then this will be avoided by asking for Calcium.
00154     //    double Ca = mpMonodomainProblem->GetMonodomainTissue()->GetCardiacCell(mWatchedElectricsNodeIndex)->GetIntracellularCalciumConcentration();
00155 
00156     *mpWatchedLocationFile << time << " ";
00157     for(unsigned i=0; i<DIM; i++)
00158     {
00159         *mpWatchedLocationFile << deformed_position[mWatchedMechanicsNodeIndex](i) << " ";
00160     }
00161     *mpWatchedLocationFile << V << "\n";
00162     mpWatchedLocationFile->flush();
00163 }
00164 
00165 //
00166 //
00168 //template<unsigned DIM>
00169 //void CardiacElectroMechanicsProblem<DIM>::SetImpactRegion(std::vector<BoundaryElement<DIM-1,DIM>*>& rImpactRegion)
00170 //{
00171 //    assert(mpImpactRegion == NULL);
00172 //    mpImpactRegion = &rImpactRegion; 
00173 //}
00174 //
00176 //template<unsigned DIM>
00177 //void CardiacElectroMechanicsProblem<DIM>::ApplyImpactTractions(double time)
00178 //{    
00179 //    if(mpImpactRegion==NULL)
00180 //    {
00181 //        return;
00182 //    }
00183 //
00184 //    double start_time = 10;
00185 //    double end_time = 20;
00186 //    double magnitude = 1.5;
00187 //    unsigned direction = 1;
00188 //    
00189 //    c_vector<double,DIM> traction = zero_vector<double>(DIM);
00190 //    if( (time>=start_time) && (time<=end_time) )
00191 //    {
00192 //        traction(direction) = magnitude;
00193 //    }
00194 //    mImpactTractions.clear();
00195 //    mImpactTractions.resize(mpImpactRegion->size(), traction);
00196 //    mpCardiacMechSolver->SetSurfaceTractionBoundaryConditions(*mpImpactRegion, mImpactTractions);
00197 //}
00198 
00199 template<unsigned DIM>
00200 CardiacElectroMechanicsProblem<DIM>::CardiacElectroMechanicsProblem(
00201             ContractionModel contractionModel,
00202             TetrahedralMesh<DIM,DIM>* pElectricsMesh,
00203             QuadraticMesh<DIM>* pMechanicsMesh,
00204             std::vector<unsigned> fixedMechanicsNodes,
00205             AbstractCardiacCellFactory<DIM>* pCellFactory,
00206             double endTime,
00207             double electricsPdeTimeStep,
00208             unsigned numElecTimeStepsPerMechTimestep,
00209             double contractionModelOdeTimeStep,
00210             std::string outputDirectory = "") :
00211         mpMeshPair(NULL)
00212 {
00213     // Start-up mechanics event handler..
00214     MechanicsEventHandler::Reset();
00215     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ALL);
00216     // disable the electric event handler, because we use a problem class but
00217     // don't call Solve, so we would have to worry about starting and ending any
00218     // events in AbstractCardiacProblem::Solve() (esp. calling EndEvent(EVERYTHING))
00219     // if we didn't disable it.
00220     HeartEventHandler::Disable();
00221 
00222     mContractionModel = contractionModel;
00223 
00224     // create the monodomain problem. Note the we use this to set up the cells,
00225     // get an initial condition (voltage) vector, and get an solver. We won't
00226     // ever call solve on the MonodomainProblem
00227     assert(pCellFactory != NULL);
00228     mpMonodomainProblem = new MonodomainProblem<DIM>(pCellFactory);
00229 
00230     // save time infomation
00231     assert(endTime > 0);
00232     mEndTime = endTime;
00233 
00234     assert(electricsPdeTimeStep>0);
00235     mElectricsTimeStep = electricsPdeTimeStep;
00236 
00237     assert(numElecTimeStepsPerMechTimestep>0);
00238 
00239     mNumElecTimestepsPerMechTimestep = numElecTimeStepsPerMechTimestep;
00240 
00241     mMechanicsTimeStep = mElectricsTimeStep*mNumElecTimestepsPerMechTimestep;
00242     assert(contractionModelOdeTimeStep <= mMechanicsTimeStep+1e-14);
00243     mContractionModelOdeTimeStep = contractionModelOdeTimeStep;
00244 
00245     // check whether output is required
00246     mWriteOutput = (outputDirectory!="");
00247     if(mWriteOutput)
00248     {
00249         mOutputDirectory = outputDirectory;
00250         // create the directory
00251         OutputFileHandler handler(mOutputDirectory);
00252         mDeformationOutputDirectory = mOutputDirectory + "/deformation";
00253         HeartConfig::Instance()->SetOutputDirectory(mOutputDirectory + "/electrics");
00254         HeartConfig::Instance()->SetOutputFilenamePrefix("voltage");
00255     }
00256     else
00257     {
00258         mDeformationOutputDirectory = "";
00259     }
00260     mNoElectricsOutput = false;
00261 
00262     // initialise all the pointers
00263     mpElectricsMesh = pElectricsMesh; // note these are allowed to be null, in case a child constructor wants to create them
00264     mpMechanicsMesh = pMechanicsMesh;
00265     mFixedNodes = fixedMechanicsNodes;
00266 
00267     mpCardiacMechSolver = NULL;
00268 
00269     // Create the Logfile (note we have to do this after the output dir has been
00270     // created, else the log file might get cleaned away
00271     std::string log_dir = mOutputDirectory; // just the TESTOUTPUT dir if mOutputDir="";
00272     LogFile::Instance()->Set(2, mOutputDirectory);
00273     LogFile::Instance()->WriteHeader("Electromechanics");
00274     LOG(2, DIM << "d Implicit CardiacElectroMechanics Simulation:");
00275     LOG(2, "End time = " << mEndTime << ", electrics time step = " << mElectricsTimeStep << ", mechanics timestep = " << mMechanicsTimeStep << "\n");
00276     LOG(2, "Contraction model ode timestep " << mContractionModelOdeTimeStep);
00277     LOG(2, "Output is written to " << mOutputDirectory << "/[deformation/electrics]");
00278 #define COVERAGE_IGNORE
00280     if(mpElectricsMesh != NULL)
00281     {
00282         LOG(2, "Electrics mesh has " << mpElectricsMesh->GetNumNodes() << " nodes");
00283     }
00284     if(mpMechanicsMesh != NULL)
00285     {
00286         LOG(2, "Mechanics mesh has " << mpMechanicsMesh->GetNumNodes() << " nodes");
00287     }
00288 #undef COVERAGE_IGNORE
00289     mIsWatchedLocation = false;
00290     mWatchedElectricsNodeIndex = UNSIGNED_UNSET;
00291     mWatchedMechanicsNodeIndex = UNSIGNED_UNSET;
00292 
00293     mFibreSheetDirectionsFile = "";
00294     mNoMechanoElectricFeedback = true;
00295 
00296 //    mpImpactRegion=NULL;
00297 }
00298 
00299 template<unsigned DIM>
00300 CardiacElectroMechanicsProblem<DIM>::~CardiacElectroMechanicsProblem()
00301 {   
00307     if(mIsWatchedLocation && mpMechanicsMesh)
00308     {
00309         mpWatchedLocationFile->close();
00310     }
00311 
00312     delete mpMonodomainProblem;
00313 
00314     delete mpCardiacMechSolver;
00315 
00316     delete mpMeshPair;
00317 
00318     LogFile::Close();
00319 }
00320 
00321 template<unsigned DIM>
00322 void CardiacElectroMechanicsProblem<DIM>::Initialise()
00323 {
00324     LOG(2, "Initialising..");
00325 
00326     assert(mpElectricsMesh!=NULL);
00327     assert(mpMechanicsMesh!=NULL);
00328     assert(mpCardiacMechSolver==NULL);
00329 
00330     if(mIsWatchedLocation)
00331     {
00332         DetermineWatchedNodes();
00333     }
00334 
00335 
00336     // initialise monodomain problem
00337     mpMonodomainProblem->SetMesh(mpElectricsMesh);
00338     HeartConfig::Instance()->SetIntracellularConductivities(Create_c_vector(1.75,1.75,1.75));
00339     mpMonodomainProblem->Initialise();
00340 
00341     // Construct mechanics solver
00342     // Here we pick the best solver for each particular contraction model. Commented out versions are
00343     // for experimentation.
00344     switch(mContractionModel)
00345     {
00346         case NASH2004:
00347             // stretch and stretch-rate independent, so should use explicit
00348             mpCardiacMechSolver = new ExplicitCardiacMechanicsSolver<DIM>(mContractionModel,mpMechanicsMesh,mDeformationOutputDirectory,mFixedNodes);
00349             break;
00350         case KERCHOFFS2003:
00351             // stretch independent, so should use implicit (explicit may be unstable)
00352             mpCardiacMechSolver = new ImplicitCardiacMechanicsSolver<DIM>(mContractionModel,mpMechanicsMesh,mDeformationOutputDirectory,mFixedNodes);
00353             break;
00354         case NHS:
00355             // stretch and stretch-rate independent, so should definitely use implicit
00356             mpCardiacMechSolver = new ImplicitCardiacMechanicsSolver<DIM>(mContractionModel,mpMechanicsMesh,mDeformationOutputDirectory,mFixedNodes);
00357             break;
00358         default:
00359             EXCEPTION("Invalid contraction model, options are: KERCHOFFS2003 or NHS");
00360     }
00361 
00362     if(mFibreSheetDirectionsFile!="")
00363     {
00364        mpCardiacMechSolver->SetVariableFibreSheetDirections(mFibreSheetDirectionsFile, mFibreSheetDirectionsDefinedPerQuadraturePoint);
00365     }
00366 
00367     // set up mesh pair and determine the fine mesh elements and corresponding weights for each
00368     // quadrature point in the coarse mesh
00369     mpMeshPair = new FineCoarseMeshPair<DIM>(*mpElectricsMesh, *mpMechanicsMesh);
00370     mpMeshPair->SetUpBoxesOnFineMesh();
00371     mpMeshPair->ComputeFineElementsAndWeightsForCoarseQuadPoints(*(mpCardiacMechSolver->GetQuadratureRule()), false);
00372     mpMeshPair->DeleteFineBoxCollection();
00373     
00374     if(!mNoMechanoElectricFeedback)
00375     {
00376         mpMeshPair->SetUpBoxesOnCoarseMesh();
00377         
00378         // compute the coarse elements which contain each fine node -- for transferring stretch from 
00379         // mechanics solve electrics cell models
00380         mpMeshPair->ComputeCoarseElementsForFineNodes(false);
00381 
00382         // compute the coarse elements which contain each fine element centroid -- for transferring F from
00383         // mechanics solve to electrics mesh elements 
00384         mpMeshPair->ComputeCoarseElementsForFineElementCentroids(false);
00385 
00386         // initialise the stretches saved for each mechanics element
00387         mStretchesForEachMechanicsElement.resize(mpMechanicsMesh->GetNumElements(),1.0);
00388         
00389         // initialise the store of the F in each mechanics element (one constant value of F) in each 
00390         mDeformationGradientsForEachMechanicsElement.resize(mpMechanicsMesh->GetNumElements(),identity_matrix<double>(DIM));
00391     }
00392         
00393     if(mWriteOutput)
00394     {
00395         TrianglesMeshWriter<DIM,DIM> mesh_writer(mOutputDirectory,"electrics_mesh",false);
00396         mesh_writer.WriteFilesUsingMesh(*mpElectricsMesh);
00397     }
00398 }
00399 
00400 template<unsigned DIM>
00401 void CardiacElectroMechanicsProblem<DIM>::Solve()
00402 {
00403 
00404     // initialise the meshes and mechanics solver
00405     if(mpCardiacMechSolver==NULL)
00406     {
00407         Initialise();
00408     }
00409 
00410     boost::shared_ptr<BoundaryConditionsContainer<DIM,DIM,1> > p_bcc(new BoundaryConditionsContainer<DIM,DIM,1>);
00411     p_bcc->DefineZeroNeumannOnMeshBoundary(mpElectricsMesh, 0);
00412     mpMonodomainProblem->SetBoundaryConditionsContainer(p_bcc);
00413 
00414     // get an electrics solver from the problem. Note that we don't call
00415     // Solve() on the CardiacProblem class, we do the looping here.
00416     AbstractDynamicLinearPdeSolver<DIM,DIM,1>* p_electrics_solver
00417        = mpMonodomainProblem->CreateSolver();
00418 
00419     // set up initial voltage etc
00420     Vec voltage=NULL; //This will be set and used later
00421     Vec initial_voltage = mpMonodomainProblem->CreateInitialCondition();
00422 
00423     unsigned num_quad_points = mpCardiacMechSolver->GetTotalNumQuadPoints();
00424     std::vector<double> interpolated_calcium_concs(num_quad_points, 0.0);
00425     std::vector<double> interpolated_voltages(num_quad_points, 0.0);
00426 
00427     // write the initial position
00428     unsigned counter = 0;
00429 
00430     TimeStepper stepper(0.0, mEndTime, mMechanicsTimeStep);
00431 
00432     CmguiDeformedSolutionsWriter<DIM>* p_cmgui_writer = NULL;
00433 
00434     unsigned mech_writer_counter = 0;
00435     if (mWriteOutput)
00436     {
00437         mpCardiacMechSolver->SetWriteOutput();
00438         mpCardiacMechSolver->WriteOutput(mech_writer_counter);
00439 
00440         p_cmgui_writer = new CmguiDeformedSolutionsWriter<DIM>(mOutputDirectory+"/deformation/cmgui", 
00441                                                                "solution", 
00442                                                                *(this->mpMechanicsMesh),
00443                                                                WRITE_QUADRATIC_MESH);
00444         std::vector<std::string> fields;
00445         fields.push_back("V");
00446         p_cmgui_writer->SetAdditionalFieldNames(fields);
00447         p_cmgui_writer->WriteInitialMesh();
00448 
00449 
00450         if(!mNoElectricsOutput)
00451         {
00452             mpMonodomainProblem->InitialiseWriter();
00453             mpMonodomainProblem->WriteOneStep(stepper.GetTime(), initial_voltage);
00454         }
00455 
00456         if(mIsWatchedLocation)
00457         {
00458             WriteWatchedLocationData(stepper.GetTime(), initial_voltage);
00459         }
00460     }
00461 
00462     while (!stepper.IsTimeAtEnd())
00463     {
00464         LOG(2, "\nCurrent time = " << stepper.GetTime());
00465         #ifdef MECH_VERBOSE // defined in AbstractNonlinearElasticitySolver
00466         // also output time to screen as newton solve information will be output
00467         std::cout << "\n\n ** Current time = " << stepper.GetTime() << "\n";
00468         #endif
00469 
00476         if(!mNoMechanoElectricFeedback)
00477         {
00478             //  Determine the stretch in each mechanics element (later: determine stretch, and 
00479             //  deformation gradient)
00480             mpCardiacMechSolver->ComputeDeformationGradientAndStretchInEachElement(mDeformationGradientsForEachMechanicsElement, mStretchesForEachMechanicsElement);
00481 
00482             //  Set the stretches on each of the cell models
00483             for(unsigned global_index = mpElectricsMesh->GetDistributedVectorFactory()->GetLow(); 
00484                          global_index < mpElectricsMesh->GetDistributedVectorFactory()->GetHigh();
00485                          global_index++)
00486             {
00487                 unsigned containing_elem = mpMeshPair->rGetCoarseElementsForFineNodes()[global_index];
00488                 double stretch = mStretchesForEachMechanicsElement[containing_elem];
00489                 mpMonodomainProblem->GetTissue()->GetCardiacCell(global_index)->SetStretch(stretch);
00490             }
00491             
00492             // finish #1244 (blocked on #1348)
00493             // NOW SET THE DEFORMATION GRADIENTS ON THE MONO/BI-DOMAIN SOLVER - do this once #1348 is done
00494             // something like:
00495             //   loop over electrics elements
00496             //   {
00497             //      unsigned containing_elem = mpMeshPair->rGetCoarseElementsForFineElementCentroids()[element_index];
00498             //      F = mDeformationGradientsForEachMechanicsElement[containing_element]
00499             //       ..
00500             //
00501         }
00502 
00503 
00509         LOG(2, "  Solving electrics");
00510         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::NON_MECH);
00511         for(unsigned i=0; i<mNumElecTimestepsPerMechTimestep; i++)
00512         {
00513             double current_time = stepper.GetTime() + i*mElectricsTimeStep;
00514             double next_time = stepper.GetTime() + (i+1)*mElectricsTimeStep;
00515 
00516             // solve the electrics
00517             p_electrics_solver->SetTimes(current_time, next_time, mElectricsTimeStep);
00518             p_electrics_solver->SetInitialCondition( initial_voltage );
00519 
00520             voltage = p_electrics_solver->Solve();
00521 
00522             PetscReal min_voltage, max_voltage;
00523             VecMax(voltage,PETSC_NULL,&max_voltage); //the second param is where the index would be returned
00524             VecMin(voltage,PETSC_NULL,&min_voltage);
00525             if(i==0)
00526             {
00527                 LOG(2, "  minimum and maximum voltage is " << min_voltage <<", "<<max_voltage);
00528             }
00529             else if(i==1)
00530             {
00531                 LOG(2, "  ..");
00532             }
00533 
00534             VecDestroy(initial_voltage);
00535             initial_voltage = voltage;
00536         }
00537 
00539 //        p_electrics_solver->SetMatrixIsNotAssembled();
00540 
00541 
00547 
00548         // compute Ca_I at each quad point (by interpolation, using the info on which
00549         // electrics element the quad point is in. Then set Ca_I on the mechanics solver
00550         LOG(2, "  Interpolating Ca_I and voltage");
00551 
00552         ReplicatableVector voltage_repl(voltage);
00553 
00554         // collect all the calcium concentrations (from the cells, which are
00555         // distributed) in one (replicated) vector
00556         ReplicatableVector calcium_repl(mpElectricsMesh->GetNumNodes());
00557         PetscInt lo, hi;
00558         VecGetOwnershipRange(voltage, &lo, &hi);
00559         for(int index=lo; index<hi; index++)
00560         {
00561             calcium_repl[index] = mpMonodomainProblem->GetTissue()->GetCardiacCell(index)->GetIntracellularCalciumConcentration();
00562         }
00563         PetscTools::Barrier();
00564         calcium_repl.Replicate(lo,hi);
00565 
00566 
00567         for(unsigned i=0; i<mpMeshPair->rGetElementsAndWeights().size(); i++)
00568         {
00569             double interpolated_CaI = 0;
00570             double interpolated_voltage = 0;
00571 
00572             Element<DIM,DIM>& element = *(mpElectricsMesh->GetElement(mpMeshPair->rGetElementsAndWeights()[i].ElementNum));
00573             for(unsigned node_index = 0; node_index<element.GetNumNodes(); node_index++)
00574             {
00575                 unsigned global_node_index = element.GetNodeGlobalIndex(node_index);
00576                 double CaI_at_node =  calcium_repl[global_node_index];
00577                 interpolated_CaI += CaI_at_node*mpMeshPair->rGetElementsAndWeights()[i].Weights(node_index);
00578                 interpolated_voltage += voltage_repl[global_node_index]*mpMeshPair->rGetElementsAndWeights()[i].Weights(node_index);
00579             }
00580 
00581             interpolated_calcium_concs[i] = interpolated_CaI;
00582             interpolated_voltages[i] = interpolated_voltage;
00583         }
00584 
00585         LOG(2, "  Setting Ca_I. max value = " << Max(interpolated_calcium_concs));
00586 
00587         // NOTE IF NHS: HERE WE SHOULD PERHAPS CHECK WHETHER THE CELL MODELS HAVE Ca_Trop
00588         // AND UPDATE FROM NHS TO CELL_MODEL, BUT NOT SURE HOW TO DO THIS.. (esp for implicit)
00589 
00590         // set [Ca], V, t
00591         mpCardiacMechSolver->SetCalciumAndVoltage(interpolated_calcium_concs, interpolated_voltages);
00592         MechanicsEventHandler::EndEvent(MechanicsEventHandler::NON_MECH);
00593 
00594 
00595 
00601         LOG(2, "  Solving mechanics ");
00602         mpCardiacMechSolver->SetWriteOutput(false);
00603 
00605 //        ApplyImpactTractions(stepper.GetTime());
00606 
00607         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ALL_MECH);
00608         mpCardiacMechSolver->Solve(stepper.GetTime(), stepper.GetNextTime(), mContractionModelOdeTimeStep);
00609         MechanicsEventHandler::EndEvent(MechanicsEventHandler::ALL_MECH);
00610 
00611         LOG(2, "    Number of newton iterations = " << mpCardiacMechSolver->GetNumNewtonIterations());
00612 
00613          // update the current time
00614         stepper.AdvanceOneTimeStep();
00615         counter++;
00616 
00617 
00623         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::OUTPUT);
00624         if(mWriteOutput && (counter%WRITE_EVERY_NTH_TIME==0))
00625         {
00626             LOG(2, "  Writing output");
00627             // write deformed position
00628             mech_writer_counter++;
00629             mpCardiacMechSolver->SetWriteOutput();
00630             mpCardiacMechSolver->WriteOutput(mech_writer_counter);
00631 
00632             p_cmgui_writer->WriteDeformationPositions(rGetDeformedPosition(), counter);
00633 
00634             if(!mNoElectricsOutput)
00635             {
00636                 mpMonodomainProblem->mpWriter->AdvanceAlongUnlimitedDimension();
00637                 mpMonodomainProblem->WriteOneStep(stepper.GetTime(), voltage);
00638             }
00639 
00640             if(mIsWatchedLocation)
00641             {
00642                 WriteWatchedLocationData(stepper.GetTime(), voltage);
00643             }
00644         }
00645         MechanicsEventHandler::EndEvent(MechanicsEventHandler::OUTPUT);
00646 
00647         // write the total elapsed time..
00648         LogFile::Instance()->WriteElapsedTime("  ");
00649     }
00650 
00651 
00652 
00653     if ((mWriteOutput) && (!mNoElectricsOutput))
00654     {
00655         HeartConfig::Instance()->Reset();
00656         mpMonodomainProblem->mpWriter->Close();
00657         delete mpMonodomainProblem->mpWriter;
00658 
00659         std::string input_dir = mOutputDirectory+"/electrics";
00660         
00661         // Convert simulation data to meshalyzer format - commented
00662         std::string config_directory = HeartConfig::Instance()->GetOutputDirectory();
00663         HeartConfig::Instance()->SetOutputDirectory(input_dir);
00664         
00665         //Hdf5ToMeshalyzerConverter<DIM,DIM> meshalyzer_converter(input_dir, "voltage", mpElectricsMesh);
00666         
00667         // convert output to CMGUI format        
00668         Hdf5ToCmguiConverter<DIM,DIM> cmgui_converter(input_dir,"voltage",mpElectricsMesh);
00669 
00670         // Write mesh in a suitable form for meshalyzer
00671         //std::string output_directory =  mOutputDirectory + "/electrics/output";
00672         // Write the mesh
00673         //MeshalyzerMeshWriter<DIM,DIM> mesh_writer(output_directory, "mesh", false);
00674         //mesh_writer.WriteFilesUsingMesh(*mpElectricsMesh);
00675         // Write the parameters out
00676         //HeartConfig::Instance()->Write();
00677     
00678         // interpolate the electrical data onto the mechanics mesh nodes and write CMGUI...
00679         // Note: this calculates the data on ALL nodes of the mechanics mesh (incl internal,
00680         // non-vertex ones), which won't be used if linear CMGUI visualisation 
00681         // of the mechanics solution is used.
00682         VoltageInterpolaterOntoMechanicsMesh<DIM> cnverter(*mpElectricsMesh,*mpMechanicsMesh,input_dir,"voltage");
00683         
00684         // reset to the default value
00685         HeartConfig::Instance()->SetOutputDirectory(config_directory);
00686     }
00687 
00688     if(p_cmgui_writer)
00689     {
00690         if(mNoElectricsOutput)
00691         {
00692             p_cmgui_writer->WriteCmguiScript();
00693         }
00694         else
00695         {
00696             p_cmgui_writer->WriteCmguiScript("../../electrics/cmgui_output/voltage_mechanics_mesh");
00697         }
00698         delete p_cmgui_writer;
00699     }
00700     VecDestroy(voltage);
00701     delete p_electrics_solver;
00702 
00703     MechanicsEventHandler::EndEvent(MechanicsEventHandler::ALL);
00704 }
00705 
00706 
00707 
00708 template<unsigned DIM>
00709 double CardiacElectroMechanicsProblem<DIM>::Max(std::vector<double>& vec)
00710 {
00711     double max = -1e200;
00712     for(unsigned i=0; i<vec.size(); i++)
00713     {
00714         if(vec[i]>max) max=vec[i];
00715     }
00716     return max;
00717 }
00718 
00719 template<unsigned DIM>
00720 void CardiacElectroMechanicsProblem<DIM>::SetNoElectricsOutput()
00721 {
00722     mNoElectricsOutput = true;
00723 }
00724 
00725 template<unsigned DIM>
00726 void CardiacElectroMechanicsProblem<DIM>::SetWatchedPosition(c_vector<double,DIM> watchedLocation)
00727 {
00728     mIsWatchedLocation = true;
00729     mWatchedLocation = watchedLocation;
00730 }
00731 
00732 template<unsigned DIM>
00733 void CardiacElectroMechanicsProblem<DIM>::SetVariableFibreSheetDirectionsFile(std::string fibreSheetDirectionsFile, bool definedPerQuadraturePoint)
00734 {
00735     mFibreSheetDirectionsFile = fibreSheetDirectionsFile;
00736     mFibreSheetDirectionsDefinedPerQuadraturePoint = definedPerQuadraturePoint;
00737 }
00738 
00739 template<unsigned DIM>
00740 std::vector<c_vector<double,DIM> >& CardiacElectroMechanicsProblem<DIM>::rGetDeformedPosition()
00741 {
00742     return mpCardiacMechSolver->rGetDeformedPosition();
00743 }
00744 
00745 
00747 // Explicit instantiation
00749 
00750 //note: 1d incompressible material doesn't make sense
00751 template class CardiacElectroMechanicsProblem<2>;
00752 template class CardiacElectroMechanicsProblem<3>;

Generated on Mon Nov 1 12:35:17 2010 for Chaste by  doxygen 1.5.5