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 "AbstractExtendedBidomainSolver.hpp"
00030
00031 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00032 void AbstractExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00033 {
00034 if (this->mpLinearSystem != NULL)
00035 {
00036 return;
00037 }
00038
00039
00040 AbstractDynamicLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,3>::InitialiseForSolve(initialSolution);
00041
00042 if (HeartConfig::Instance()->GetUseAbsoluteTolerance())
00043 {
00044 #ifdef TRACE_KSP
00045 std::cout << "Using absolute tolerance: " << mpConfig->GetAbsoluteTolerance() <<"\n";
00046 #endif
00047 this->mpLinearSystem->SetAbsoluteTolerance(mpConfig->GetAbsoluteTolerance());
00048 }
00049 else
00050 {
00051 #ifdef TRACE_KSP
00052 std::cout << "Using relative tolerance: " << mpConfig->GetRelativeTolerance() <<"\n";
00053 #endif
00054 this->mpLinearSystem->SetRelativeTolerance(mpConfig->GetRelativeTolerance());
00055 }
00056
00057 this->mpLinearSystem->SetKspType(HeartConfig::Instance()->GetKSPSolver());
00058 this->mpLinearSystem->SetPcType(HeartConfig::Instance()->GetKSPPreconditioner());
00059
00060 if (mRowForAverageOfPhiZeroed==INT_MAX)
00061 {
00062
00063 this->mpLinearSystem->SetMatrixIsSymmetric(true);
00064 }
00065 else
00066 {
00067
00068 this->mpLinearSystem->SetMatrixIsSymmetric(false);
00069 }
00070 }
00071
00072
00073 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00074 void AbstractExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::PrepareForSetupLinearSystem(Vec existingSolution)
00075 {
00076 double time = PdeSimulationTime::GetTime();
00077 double delta_t = PdeSimulationTime::GetPdeTimeStep();
00078 mpExtendedBidomainTissue->SolveCellSystems(existingSolution, time, time + delta_t);
00079 }
00080
00081 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00082 Vec AbstractExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::GenerateNullBasis() const
00083 {
00084 double sqrrt_num_nodes = sqrt( this->mpMesh->GetNumNodes());
00085 double normalised_basis_value = 1.0 / sqrrt_num_nodes;
00086
00087 Vec nullbasis;
00088 DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00089 nullbasis=p_factory->CreateVec(3);
00090 DistributedVector dist_null_basis = p_factory->CreateDistributedVector(nullbasis);
00091 DistributedVector::Stripe null_basis_stripe_0(dist_null_basis,0);
00092 DistributedVector::Stripe null_basis_stripe_1(dist_null_basis,1);
00093 DistributedVector::Stripe null_basis_stripe_2(dist_null_basis,2);
00094 for (DistributedVector::Iterator index = dist_null_basis.Begin();
00095 index != dist_null_basis.End();
00096 ++index)
00097 {
00098 null_basis_stripe_0[index] = 0.0;
00099 null_basis_stripe_1[index] = 0.0;
00100 null_basis_stripe_2[index] = normalised_basis_value;
00101 }
00102 dist_null_basis.Restore();
00103 return nullbasis;
00104 }
00105
00106
00107 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00108 void AbstractExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::FinaliseLinearSystem(Vec existingSolution)
00109 {
00110 if (!(this->mpBoundaryConditions->HasDirichletBoundaryConditions()))
00111 {
00112
00113 if (mRowForAverageOfPhiZeroed==INT_MAX)
00114 {
00115
00116 if (!mNullSpaceCreated)
00117 {
00118
00119 Vec nullbasis[] = {GenerateNullBasis()};
00120
00121 this->mpLinearSystem->SetNullBasis(nullbasis, 1);
00122
00123 VecDestroy(nullbasis[0]);
00124 mNullSpaceCreated = true;
00125 }
00126 }
00127 else
00128 {
00129
00130 this->mpLinearSystem->SetKspType("gmres");
00131 mpConfig->SetKSPSolver("gmres");
00132
00133
00134 unsigned matrix_size = this->mpLinearSystem->GetSize();
00135 if (!this->mMatrixIsAssembled)
00136 {
00137
00138
00139 std::vector<unsigned> row_for_average;
00140 row_for_average.push_back(mRowForAverageOfPhiZeroed);
00141 this->mpLinearSystem->ZeroMatrixRowsWithValueOnDiagonal(row_for_average, 0.0);
00142 for (unsigned col_index=0; col_index<matrix_size; col_index++)
00143 {
00144 if (((col_index-2)%3 == 0) && (col_index>1))
00145 {
00146 this->mpLinearSystem->SetMatrixElement(mRowForAverageOfPhiZeroed, col_index, 1);
00147 }
00148
00149 }
00150 this->mpLinearSystem->FinaliseLhsMatrix();
00151 }
00152
00153 this->mpLinearSystem->SetRhsVectorElement(mRowForAverageOfPhiZeroed, 0);
00154
00155 this->mpLinearSystem->FinaliseRhsVector();
00156 }
00157 }
00158 CheckCompatibilityCondition();
00159 }
00160
00161 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00162 void AbstractExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::CheckCompatibilityCondition()
00163 {
00164 if (this->mpBoundaryConditions->HasDirichletBoundaryConditions() || mRowForAverageOfPhiZeroed!=INT_MAX )
00165 {
00166
00167 return;
00168 }
00169
00170 #ifndef NDEBUG
00171 DistributedVector distributed_rhs=this->mpMesh->GetDistributedVectorFactory()->CreateDistributedVector(this->mpLinearSystem->rGetRhsVector());
00172 DistributedVector::Stripe distributed_rhs_phi_entries(distributed_rhs,2);
00173
00174 double local_sum=0;
00175 for (DistributedVector::Iterator index = distributed_rhs.Begin();
00176 index!= distributed_rhs.End();
00177 ++index)
00178 {
00179 local_sum += distributed_rhs_phi_entries[index];
00180 }
00181
00182 double global_sum;
00183 int mpi_ret = MPI_Allreduce(&local_sum, &global_sum, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
00184 assert(mpi_ret == MPI_SUCCESS);
00185
00186 if(fabs(global_sum)>1e-6)
00187 {
00188 #define COVERAGE_IGNORE
00189
00190 std::cout << "Sum of b_{every 3 items} = " << global_sum << " (should be zero for compatibility)\n";
00191 EXCEPTION("Linear system does not satisfy compatibility constraint!");
00192 #undef COVERAGE_IGNORE
00193 }
00194 #endif
00195 }
00196
00197 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00198 AbstractExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::AbstractExtendedBidomainSolver(
00199 bool bathSimulation,
00200 AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00201 ExtendedBidomainTissue<SPACE_DIM>* pTissue,
00202 BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, 3>* pBcc,
00203 unsigned numQuadPoints)
00204 : AbstractDynamicLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,3>(pMesh),
00205 mBathSimulation(bathSimulation),
00206 mpExtendedBidomainTissue(pTissue),
00207 mpBoundaryConditions(pBcc),
00208 mNumQuadPoints(numQuadPoints)
00209 {
00210 assert(pTissue != NULL);
00211 assert(pBcc != NULL);
00212
00213 mNullSpaceCreated = false;
00214
00215
00216 this->mMatrixIsConstant = true;
00217
00218 mRowForAverageOfPhiZeroed = INT_MAX;
00219 mpConfig = HeartConfig::Instance();
00220
00221 mpExtendedBidomainAssembler = NULL;
00222 }
00223
00224 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00225 AbstractExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::~AbstractExtendedBidomainSolver()
00226 {
00227 }
00228
00229 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00230 void AbstractExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::SetFixedExtracellularPotentialNodes(
00231 std::vector<unsigned> fixedExtracellularPotentialNodes)
00232 {
00233 for (unsigned i=0; i<fixedExtracellularPotentialNodes.size(); i++)
00234 {
00235 if (fixedExtracellularPotentialNodes[i] >= this->mpMesh->GetNumNodes() )
00236 {
00237 EXCEPTION("Fixed node number must be less than total number nodes");
00238 }
00239 }
00240
00241 mFixedExtracellularPotentialNodes = fixedExtracellularPotentialNodes;
00242
00243
00244 this->mpBoundaryConditions->ResetDirichletCommunication();
00245
00246 for (unsigned i=0; i<mFixedExtracellularPotentialNodes.size(); i++)
00247 {
00248 if (this->mpMesh->GetDistributedVectorFactory()->IsGlobalIndexLocal(mFixedExtracellularPotentialNodes[i]))
00249 {
00250 ConstBoundaryCondition<SPACE_DIM>* p_boundary_condition = new ConstBoundaryCondition<SPACE_DIM>(0.0);
00251
00252
00253 Node<SPACE_DIM>* p_node = this->mpMesh->GetNode(mFixedExtracellularPotentialNodes[i]);
00254
00255
00256 GetBoundaryConditions()->AddDirichletBoundaryCondition(p_node, p_boundary_condition, 2, false);
00257 }
00258 }
00259 }
00260
00261 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00262 void AbstractExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::SetRowForAverageOfPhiZeroed(unsigned row)
00263 {
00264
00265 if ( ((row-2)%3 != 0) && row!=INT_MAX)
00266 {
00267 EXCEPTION("Row for applying the constraint 'Average of phi_e = zero' should be every 3 rows");
00268 }
00269
00270 mRowForAverageOfPhiZeroed = row;
00271 }
00272
00274
00276
00277 template class AbstractExtendedBidomainSolver<1,1>;
00278 template class AbstractExtendedBidomainSolver<2,2>;
00279 template class AbstractExtendedBidomainSolver<3,3>;