BidomainDg0Assembler.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2009
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 "BidomainDg0Assembler.hpp"
00030 #include <boost/numeric/ublas/vector_proxy.hpp>
00031 
00032 #include "Exception.hpp"
00033 #include "DistributedVector.hpp"
00034 #include "ConstBoundaryCondition.hpp"
00035 
00036 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00037 void BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ResetInterpolatedQuantities( void )
00038 {
00039     mIionic=0;
00040     mIIntracellularStimulus=0;
00041 }
00042 
00043 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00044 void BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00045 {
00046     if (this->mpLinearSystem != NULL)
00047     {
00048         return;
00049     }
00050 
00051     // linear system created here
00052     BaseClassType::InitialiseForSolve(initialSolution);
00053 
00054     if (HeartConfig::Instance()->GetUseAbsoluteTolerance())
00055     {
00056 #ifdef TRACE_KSP        
00057         std::cout << "Using absolute tolerance: " << mpConfig->GetAbsoluteTolerance() <<"\n";
00058 #endif
00059         this->mpLinearSystem->SetAbsoluteTolerance(mpConfig->GetAbsoluteTolerance());
00060     }
00061     else
00062     {
00063 #ifdef TRACE_KSP        
00064         std::cout << "Using relative tolerance: " << mpConfig->GetRelativeTolerance() <<"\n";
00065 #endif
00066         this->mpLinearSystem->SetRelativeTolerance(mpConfig->GetRelativeTolerance());
00067     }
00068 
00069     this->mpLinearSystem->SetKspType(HeartConfig::Instance()->GetKSPSolver());
00070     this->mpLinearSystem->SetPcType(HeartConfig::Instance()->GetKSPPreconditioner());
00071 
00072     if (mRowForAverageOfPhiZeroed==INT_MAX)
00073     {   
00074         // not applying average(phi)=0 constraint, so matrix is symmetric     
00075         this->mpLinearSystem->SetMatrixIsSymmetric(true);
00076     }
00077     else
00078     {    
00079         // applying average(phi)=0 constraint, so matrix is not symmetric      
00080         this->mpLinearSystem->SetMatrixIsSymmetric(false);
00081     }
00082 }
00083 
00084 
00085 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00086 void BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::IncrementInterpolatedQuantities(double phiI, const Node<SPACE_DIM>* pNode)
00087 {
00088     unsigned node_global_index = pNode->GetIndex();
00089 
00090     mIionic                 += phiI * mpBidomainPde->rGetIionicCacheReplicated()[ node_global_index ];
00091     mIIntracellularStimulus += phiI * mpBidomainPde->rGetIntracellularStimulusCacheReplicated()[ node_global_index ];
00092 }
00093 
00094 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00095 c_matrix<double,2*(ELEMENT_DIM+1),2*(ELEMENT_DIM+1)>
00096     BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeMatrixTerm(
00097             c_vector<double, ELEMENT_DIM+1> &rPhi,
00098             c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> &rGradPhi,
00099             ChastePoint<SPACE_DIM> &rX,
00100             c_vector<double,2> &rU,
00101             c_matrix<double, 2, SPACE_DIM> &rGradU /* not used */,
00102             Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00103 {
00104     // get bidomain parameters
00105     double Am = mpConfig->GetSurfaceAreaToVolumeRatio();
00106     double Cm = mpConfig->GetCapacitance();
00107 
00108     const c_matrix<double, SPACE_DIM, SPACE_DIM>& sigma_i = mpBidomainPde->rGetIntracellularConductivityTensor(pElement->GetIndex());
00109     const c_matrix<double, SPACE_DIM, SPACE_DIM>& sigma_e = mpBidomainPde->rGetExtracellularConductivityTensor(pElement->GetIndex());
00110 
00111 
00112     c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> temp = prod(sigma_i, rGradPhi);
00113     c_matrix<double, ELEMENT_DIM+1, ELEMENT_DIM+1> grad_phi_sigma_i_grad_phi =
00114         prod(trans(rGradPhi), temp);
00115 
00116     c_matrix<double, ELEMENT_DIM+1, ELEMENT_DIM+1> basis_outer_prod =
00117         outer_prod(rPhi, rPhi);
00118 
00119     c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> temp2 = prod(sigma_e, rGradPhi);
00120     c_matrix<double, ELEMENT_DIM+1, ELEMENT_DIM+1> grad_phi_sigma_e_grad_phi =
00121         prod(trans(rGradPhi), temp2);
00122 
00123 
00124     c_matrix<double,2*(ELEMENT_DIM+1),2*(ELEMENT_DIM+1)> ret;
00125 
00126     // even rows, even columns
00127     matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00128     slice00(ret, slice (0, 2, ELEMENT_DIM+1), slice (0, 2, ELEMENT_DIM+1));
00129     slice00 = (Am*Cm/this->mDt)*basis_outer_prod + grad_phi_sigma_i_grad_phi;
00130 
00131     // odd rows, even columns
00132     matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00133     slice10(ret, slice (1, 2, ELEMENT_DIM+1), slice (0, 2, ELEMENT_DIM+1));
00134     slice10 = grad_phi_sigma_i_grad_phi;
00135 
00136     // even rows, odd columns
00137     matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00138     slice01(ret, slice (0, 2, ELEMENT_DIM+1), slice (1, 2, ELEMENT_DIM+1));
00139     slice01 = grad_phi_sigma_i_grad_phi;
00140 
00141     // odd rows, odd columns
00142     matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00143     slice11(ret, slice (1, 2, ELEMENT_DIM+1), slice (1, 2, ELEMENT_DIM+1));
00144     slice11 = grad_phi_sigma_i_grad_phi + grad_phi_sigma_e_grad_phi;
00145 
00146     return ret;
00147 }
00148 
00149 
00150 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00151 c_vector<double,2*(ELEMENT_DIM+1)>
00152     BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeVectorTerm(
00153             c_vector<double, ELEMENT_DIM+1> &rPhi,
00154             c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> &rGradPhi,
00155             ChastePoint<SPACE_DIM> &rX,
00156             c_vector<double,2> &u,
00157             c_matrix<double, 2, SPACE_DIM> &rGradU /* not used */,
00158             Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00159 {
00160     // get bidomain parameters
00161     double Am = mpConfig->GetSurfaceAreaToVolumeRatio();
00162     double Cm = mpConfig->GetCapacitance();
00163 
00164     c_vector<double,2*(ELEMENT_DIM+1)> ret;
00165 
00166     vector_slice<c_vector<double, 2*(ELEMENT_DIM+1)> > slice_V  (ret, slice (0, 2, ELEMENT_DIM+1));
00167     vector_slice<c_vector<double, 2*(ELEMENT_DIM+1)> > slice_Phi(ret, slice (1, 2, ELEMENT_DIM+1));
00168 
00169     // u(0) = voltage
00170     noalias(slice_V)   = (Am*Cm*u(0)/this->mDt - Am*mIionic - mIIntracellularStimulus) * rPhi;
00171     noalias(slice_Phi) = zero_vector<double>(ELEMENT_DIM+1);
00172 
00173     return ret;
00174 }
00175 
00176 
00177 
00178 #define COVERAGE_IGNORE
00179 // never called, since no non-zero neumann conditions
00180 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00181 c_vector<double, 2*ELEMENT_DIM> BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeVectorSurfaceTerm(
00182     const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM> &rSurfaceElement,
00183     c_vector<double,ELEMENT_DIM> &rPhi,
00184     ChastePoint<SPACE_DIM> &rX)
00185 {
00186     // D_times_gradu_dot_n = [D grad(u)].n, D=diffusion matrix
00187     double sigma_i_times_grad_phi_i_dot_n = this->mpBoundaryConditions->GetNeumannBCValue(&rSurfaceElement, rX, 0);
00188     double sigma_e_times_grad_phi_e_dot_n = this->mpBoundaryConditions->GetNeumannBCValue(&rSurfaceElement, rX, 1);
00189 
00190     c_vector<double, 2*ELEMENT_DIM> ret;
00191     for (unsigned i=0; i<ELEMENT_DIM; i++)
00192     {
00193         ret(2*i)   = rPhi(i)*sigma_i_times_grad_phi_i_dot_n;
00194         ret(2*i+1) = rPhi(i)*(sigma_i_times_grad_phi_i_dot_n + sigma_e_times_grad_phi_e_dot_n);
00195     }
00196 
00197     return ret;
00198 }
00199 #undef COVERAGE_IGNORE
00200 
00201 
00202 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00203 void BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::PrepareForAssembleSystem(Vec existingSolution, double time)
00204 {
00205     mpBidomainPde->SolveCellSystems(existingSolution, time, time+this->mDt);
00206 }
00207 
00208 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00209 Vec BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::GenerateNullBasis() const
00210 {
00211     double sqrt_num_nodes = sqrt((double) this->mpMesh->GetNumNodes());
00212     
00213     Vec nullbasis;
00214     DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00215     nullbasis=p_factory->CreateVec(2);
00216     DistributedVector dist_null_basis = p_factory->CreateDistributedVector(nullbasis);
00217     DistributedVector::Stripe null_basis_stripe_0(dist_null_basis,0);
00218     DistributedVector::Stripe null_basis_stripe_1(dist_null_basis,1);
00219     for (DistributedVector::Iterator index = dist_null_basis.Begin();
00220          index != dist_null_basis.End();
00221          ++index)
00222     {
00223         null_basis_stripe_0[index] = 0.0;
00224         null_basis_stripe_1[index] = 1.0/sqrt_num_nodes; // normalised vector        
00225     }
00226     dist_null_basis.Restore();
00227     
00228     return nullbasis;    
00229 }
00230 
00231 
00232 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00233 void BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::FinaliseAssembleSystem(Vec existingSolution, double time)
00234 {
00235     if (!(this->mpBoundaryConditions->HasDirichletBoundaryConditions()))
00236     {
00237         // We're not pinning any nodes.
00238         if (mRowForAverageOfPhiZeroed==INT_MAX)
00239         {
00240             // We're not using the 'Average phi_e = 0' method, hence use a null space
00241             if (!mNullSpaceCreated)
00242             {
00243                 // No null space set up, so create one and pass it to the linear system                
00244                 Vec nullbasis[] = {GenerateNullBasis()};
00245                 
00246                 this->mpLinearSystem->SetNullBasis(nullbasis, 1);
00247 
00248                 VecDestroy(nullbasis[0]);
00249                 mNullSpaceCreated = true;
00250             }
00251         }
00252         else
00253         {
00254             // mRowForAverageOfPhiZeroed!=INT_MAX, i.e. we're using the 'Average phi_e = 0' method
00255             
00256             // CG (default solver) won't work since the system isn't symmetric anymore. Switch to GMRES
00257             this->mpLinearSystem->SetKspType("gmres"); // Switches the solver
00258             mpConfig->SetKSPSolver("gmres"); // Makes sure this change will be reflected in the XML file written to disk at the end of the simulation.            
00259             
00260             // Set average phi_e to zero
00261             unsigned matrix_size = this->mpLinearSystem->GetSize();
00262             if (!this->mMatrixIsAssembled)
00263             {
00264 
00265                 // Set the mRowForAverageOfPhiZeroed-th matrix row to 0 1 0 1 ...
00266                 this->mpLinearSystem->ZeroMatrixRow(mRowForAverageOfPhiZeroed);
00267                 for (unsigned col_index=0; col_index<matrix_size; col_index++)
00268                 {
00269                     if (col_index%2 == 1)
00270                     {
00271                         this->mpLinearSystem->SetMatrixElement(mRowForAverageOfPhiZeroed, col_index, 1);
00272                     }
00273                 }
00274                 this->mpLinearSystem->AssembleFinalLhsMatrix();
00275 
00276             }
00277             // Set the mRowForAverageOfPhiZeroed-th rhs vector row to 0
00278             this->mpLinearSystem->SetRhsVectorElement(mRowForAverageOfPhiZeroed, 0);
00279 
00280             this->mpLinearSystem->AssembleRhsVector();
00281         }
00282     }
00283 
00284     CheckCompatibilityCondition();
00285 }
00286 
00287 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00288 void BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::CheckCompatibilityCondition()
00289 {
00290     if(!PetscTools::GetNumProcs()>1)
00291     {
00292         // don't do test in parallel
00293         return;
00294     }  
00295     
00296     if (this->mpBoundaryConditions->HasDirichletBoundaryConditions() || mRowForAverageOfPhiZeroed!=INT_MAX )
00297     {
00298         // not a singular system, no compability condition
00299         return;
00300     }
00301 
00302 #ifndef NDEBUG
00303     ReplicatableVector rep(this->mpLinearSystem->rGetRhsVector());
00304     double sum = 0;
00305     for(unsigned i=1; i<rep.size(); i+=2) // i=1,3,5,..  ie all the phi_e components
00306     {
00307         sum += rep[i];
00308     }
00309     
00310     if(fabs(sum)>1e-6) // magic number! sum should really be a sum of zeros and exactly zero though anyway (or a-a+b-b+c-c.. etc in the case of electrodes)
00311     {
00312         #define COVERAGE_IGNORE
00313         // shouldn't ever reach this line but useful to have the error printed out if you do
00314         //std::cout << "Sum of b_{2i+1} = " << sum << " (should be zero for compatibility)\n";
00315         EXCEPTION("Linear system does not satisfy compatibility constraint!");
00316         #undef COVERAGE_IGNORE
00317     }
00318 #endif
00319 }
00320 
00321 
00322 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00323 BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::BidomainDg0Assembler(
00324             AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00325             BidomainPde<SPACE_DIM>* pPde,
00326             BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, 2>* pBcc,
00327             unsigned numQuadPoints)
00328     : AbstractAssembler<ELEMENT_DIM,SPACE_DIM,2>(),
00329       BaseClassType(numQuadPoints),
00330       AbstractDynamicAssemblerMixin<ELEMENT_DIM,SPACE_DIM,2>()
00331 {
00332     assert(pPde != NULL);
00333     assert(pMesh != NULL);
00334     assert(pBcc != NULL);
00335 
00336     mpBidomainPde = pPde;
00337     this->SetMesh(pMesh);
00338 
00339     this->mpBoundaryConditions = pBcc;
00340 
00341     mNullSpaceCreated = false;
00342 
00343     this->SetMatrixIsConstant();
00344 
00345     mRowForAverageOfPhiZeroed = INT_MAX; //this->mpLinearSystem->GetSize() - 1;
00346 
00347     mpConfig = HeartConfig::Instance();
00348 }
00349 
00350 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00351 BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::~BidomainDg0Assembler()
00352 {
00353 }
00354 
00355 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00356 void BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::SetFixedExtracellularPotentialNodes(
00357             std::vector<unsigned> fixedExtracellularPotentialNodes)
00358 {
00359     for (unsigned i=0; i<fixedExtracellularPotentialNodes.size(); i++)
00360     {
00361         if (fixedExtracellularPotentialNodes[i] >= this->mpMesh->GetNumNodes() )
00362         {
00363             EXCEPTION("Fixed node number must be less than total number nodes");
00364         }
00365     }
00366 
00367     mFixedExtracellularPotentialNodes = fixedExtracellularPotentialNodes;
00368 
00369     for (unsigned i=0; i<mFixedExtracellularPotentialNodes.size(); i++)
00370     {
00371         ConstBoundaryCondition<SPACE_DIM>* p_boundary_condition
00372          = new ConstBoundaryCondition<SPACE_DIM>(0.0);
00373 
00374         Node<SPACE_DIM>* p_node = this->mpMesh->GetNode(mFixedExtracellularPotentialNodes[i]);
00375 
00376         this->mpBoundaryConditions->AddDirichletBoundaryCondition(p_node, p_boundary_condition, 1);
00377     }
00378 }
00379 
00380 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00381 void BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::SetRowForAverageOfPhiZeroed(unsigned row)
00382 {
00383     // Row should be odd in C++-like indexing
00384     if (row % 2 == 0)
00385     {
00386         EXCEPTION("Row for applying the constraint 'Average of phi_e = zero' should be odd in C++ like indexing");
00387     }
00388 
00389     mRowForAverageOfPhiZeroed = row;
00390 }
00391 
00393 // Explicit instantiation
00395 
00396 template class BidomainDg0Assembler<1,1>;
00397 template class BidomainDg0Assembler<2,2>;
00398 template class BidomainDg0Assembler<3,3>;

Generated on Tue Aug 4 16:10:22 2009 for Chaste by  doxygen 1.5.5