Chaste Release::3.1
|
00001 /* 00002 00003 Copyright (c) 2005-2012, University of Oxford. 00004 All rights reserved. 00005 00006 University of Oxford means the Chancellor, Masters and Scholars of the 00007 University of Oxford, having an administrative office at Wellington 00008 Square, Oxford OX1 2JD, UK. 00009 00010 This file is part of Chaste. 00011 00012 Redistribution and use in source and binary forms, with or without 00013 modification, are permitted provided that the following conditions are met: 00014 * Redistributions of source code must retain the above copyright notice, 00015 this list of conditions and the following disclaimer. 00016 * Redistributions in binary form must reproduce the above copyright notice, 00017 this list of conditions and the following disclaimer in the documentation 00018 and/or other materials provided with the distribution. 00019 * Neither the name of the University of Oxford nor the names of its 00020 contributors may be used to endorse or promote products derived from this 00021 software without specific prior written permission. 00022 00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" 00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE 00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE 00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR 00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE 00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) 00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT 00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT 00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00033 00034 */ 00035 00036 #include "OperatorSplittingMonodomainSolver.hpp" 00037 00038 00039 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM> 00040 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(Vec currentSolution, bool computeMatrix) 00041 { 00042 assert(this->mpLinearSystem->rGetLhsMatrix() != NULL); 00043 assert(this->mpLinearSystem->rGetRhsVector() != NULL); 00044 00046 // set up LHS matrix (and mass matrix) 00048 if(computeMatrix) 00049 { 00050 mpMonodomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix()); 00051 mpMonodomainAssembler->AssembleMatrix(); 00052 00053 MassMatrixAssembler<ELEMENT_DIM,SPACE_DIM> mass_matrix_assembler(this->mpMesh, HeartConfig::Instance()->GetUseMassLumping()); 00054 mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix); 00055 mass_matrix_assembler.Assemble(); 00056 00057 this->mpLinearSystem->FinaliseLhsMatrix(); 00058 PetscMatTools::Finalise(mMassMatrix); 00059 } 00060 00061 HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS); 00062 00064 // Set up z in b=Mz 00066 DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory(); 00067 // dist stripe for the current Voltage 00068 DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution); 00069 // dist stripe for z (return value) 00070 DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs); 00071 00072 double Am = HeartConfig::Instance()->GetSurfaceAreaToVolumeRatio(); 00073 double Cm = HeartConfig::Instance()->GetCapacitance(); 00074 00075 for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin(); 00076 index!= dist_vec_matrix_based.End(); 00077 ++index) 00078 { 00079 double V = distributed_current_solution[index]; 00080 // in the main solver, the nodal ionic current and stimuli is computed and used. 00081 // However in operator splitting, this part of the solve is diffusion only, no reaction terms 00082 //double F = - Am*this->mpMonodomainTissue->rGetIionicCacheReplicated()[index.Global] 00083 // - this->mpMonodomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global]; 00084 00085 dist_vec_matrix_based[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse(); 00086 } 00087 dist_vec_matrix_based.Restore(); 00088 00090 // b = Mz 00092 MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector()); 00093 00094 // assembling RHS is not finished yet, as Neumann bcs are added below, but 00095 // the event will be begun again inside mpMonodomainAssembler->AssembleVector(); 00096 HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS); 00097 00099 // apply Neumann boundary conditions 00101 mpNeumannSurfaceTermsAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/); 00102 mpNeumannSurfaceTermsAssembler->AssembleVector(); 00103 00104 // finalise 00105 this->mpLinearSystem->FinaliseRhsVector(); 00106 } 00107 00108 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM> 00109 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::PrepareForSetupLinearSystem(Vec currentSolution) 00110 { 00111 double time = PdeSimulationTime::GetTime(); 00112 double dt = PdeSimulationTime::GetPdeTimeStep(); 00113 mpMonodomainTissue->SolveCellSystems(currentSolution, time, time+dt/2.0, true); 00114 } 00115 00116 00117 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM> 00118 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::FollowingSolveLinearSystem(Vec currentSolution) 00119 { 00120 // solve cell models for second half timestep 00121 double time = PdeSimulationTime::GetTime(); 00122 double dt = PdeSimulationTime::GetPdeTimeStep(); 00123 mpMonodomainTissue->SolveCellSystems(currentSolution, time + dt/2, time+dt, true); 00124 } 00125 00126 00127 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM> 00128 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution) 00129 { 00130 if (this->mpLinearSystem != NULL) 00131 { 00132 return; 00133 } 00134 00135 // call base class version... 00136 AbstractLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>::InitialiseForSolve(initialSolution); 00137 00138 //..then do a bit extra 00139 if(HeartConfig::Instance()->GetUseAbsoluteTolerance()) 00140 { 00141 this->mpLinearSystem->SetAbsoluteTolerance(HeartConfig::Instance()->GetAbsoluteTolerance()); 00142 } 00143 else 00144 { 00145 NEVER_REACHED; 00146 // re-implement when needed 00147 //this->mpLinearSystem->SetRelativeTolerance(HeartConfig::Instance()->GetRelativeTolerance()); 00148 } 00149 00150 this->mpLinearSystem->SetKspType(HeartConfig::Instance()->GetKSPSolver()); 00151 this->mpLinearSystem->SetPcType(HeartConfig::Instance()->GetKSPPreconditioner()); 00152 this->mpLinearSystem->SetMatrixIsSymmetric(true); 00153 this->mpLinearSystem->SetUseFixedNumberIterations(HeartConfig::Instance()->GetUseFixedNumberIterationsLinearSolver(), HeartConfig::Instance()->GetEvaluateNumItsEveryNSolves()); 00154 00155 // initialise matrix-based RHS vector and matrix, and use the linear 00156 // system rhs as a template 00157 Vec& r_template = this->mpLinearSystem->rGetRhsVector(); 00158 VecDuplicate(r_template, &mVecForConstructingRhs); 00159 PetscInt ownership_range_lo; 00160 PetscInt ownership_range_hi; 00161 VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi); 00162 PetscInt local_size = ownership_range_hi - ownership_range_lo; 00163 PetscTools::SetupMat(mMassMatrix, this->mpMesh->GetNumNodes(), this->mpMesh->GetNumNodes(), 00164 this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(), 00165 local_size, local_size); 00166 } 00167 00168 00169 00170 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM> 00171 OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::OperatorSplittingMonodomainSolver( 00172 AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh, 00173 MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* pTissue, 00174 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,1>* pBoundaryConditions, 00175 unsigned numQuadPoints) 00176 : AbstractDynamicLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>(pMesh), 00177 mpBoundaryConditions(pBoundaryConditions), 00178 mpMonodomainTissue(pTissue), 00179 mNumQuadPoints(numQuadPoints) 00180 { 00181 assert(pTissue); 00182 assert(pBoundaryConditions); 00183 this->mMatrixIsConstant = true; 00184 00185 mpMonodomainAssembler = new MonodomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints); 00186 mpNeumannSurfaceTermsAssembler = new NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,1>(pMesh,pBoundaryConditions); 00187 00188 // Tell tissue there's no need to replicate ionic caches 00189 pTissue->SetCacheReplication(false); 00190 mVecForConstructingRhs = NULL; 00191 } 00192 00193 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM> 00194 OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::~OperatorSplittingMonodomainSolver() 00195 { 00196 delete mpMonodomainAssembler; 00197 delete mpNeumannSurfaceTermsAssembler; 00198 00199 if(mVecForConstructingRhs) 00200 { 00201 PetscTools::Destroy(mVecForConstructingRhs); 00202 PetscTools::Destroy(mMassMatrix); 00203 } 00204 } 00205 00206 00207 00209 // explicit instantiation 00211 00212 template class OperatorSplittingMonodomainSolver<1,1>; 00213 template class OperatorSplittingMonodomainSolver<1,2>; 00214 template class OperatorSplittingMonodomainSolver<1,3>; 00215 template class OperatorSplittingMonodomainSolver<2,2>; 00216 template class OperatorSplittingMonodomainSolver<3,3>; 00217