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     BaseClassType::InitialiseForSolve(initialSolution);
00051     if (HeartConfig::Instance()->GetUseAbsoluteTolerance())
00052     {
00053         this->mpLinearSystem->SetAbsoluteTolerance(mpConfig->GetAbsoluteTolerance());
00054     }
00055     else
00056     {
00057         this->mpLinearSystem->SetRelativeTolerance(mpConfig->GetRelativeTolerance());
00058     }
00059     this->mpLinearSystem->SetKspType(HeartConfig::Instance()->GetKSPSolver());
00060     this->mpLinearSystem->SetPcType(HeartConfig::Instance()->GetKSPPreconditioner());
00061 }
00062 
00063 
00064 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00065 void BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::IncrementInterpolatedQuantities(double phi_i, const Node<SPACE_DIM>* pNode)
00066 {
00067     unsigned node_global_index = pNode->GetIndex();
00068 
00069     mIionic                 += phi_i * mpBidomainPde->rGetIionicCacheReplicated()[ node_global_index ];
00070     mIIntracellularStimulus += phi_i * mpBidomainPde->rGetIntracellularStimulusCacheReplicated()[ node_global_index ];
00071 }
00072 
00073 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00074 c_matrix<double,2*(ELEMENT_DIM+1),2*(ELEMENT_DIM+1)>
00075     BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeMatrixTerm(
00076             c_vector<double, ELEMENT_DIM+1> &rPhi,
00077             c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> &rGradPhi,
00078             ChastePoint<SPACE_DIM> &rX,
00079             c_vector<double,2> &u,
00080             c_matrix<double, 2, SPACE_DIM> &rGradU /* not used */,
00081             Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00082 {
00083     // get bidomain parameters
00084     double Am = mpConfig->GetSurfaceAreaToVolumeRatio();
00085     double Cm = mpConfig->GetCapacitance();
00086 
00087     const c_matrix<double, SPACE_DIM, SPACE_DIM>& sigma_i = mpBidomainPde->rGetIntracellularConductivityTensor(pElement->GetIndex());
00088     const c_matrix<double, SPACE_DIM, SPACE_DIM>& sigma_e = mpBidomainPde->rGetExtracellularConductivityTensor(pElement->GetIndex());
00089 
00090 
00091     c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> temp = prod(sigma_i, rGradPhi);
00092     c_matrix<double, ELEMENT_DIM+1, ELEMENT_DIM+1> grad_phi_sigma_i_grad_phi =
00093         prod(trans(rGradPhi), temp);
00094 
00095     c_matrix<double, ELEMENT_DIM+1, ELEMENT_DIM+1> basis_outer_prod =
00096         outer_prod(rPhi, rPhi);
00097 
00098     c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> temp2 = prod(sigma_e, rGradPhi);
00099     c_matrix<double, ELEMENT_DIM+1, ELEMENT_DIM+1> grad_phi_sigma_e_grad_phi =
00100         prod(trans(rGradPhi), temp2);
00101 
00102 
00103     c_matrix<double,2*(ELEMENT_DIM+1),2*(ELEMENT_DIM+1)> ret;
00104 
00105     // even rows, even columns
00106     matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00107     slice00(ret, slice (0, 2, ELEMENT_DIM+1), slice (0, 2, ELEMENT_DIM+1));
00108     slice00 = (Am*Cm/this->mDt)*basis_outer_prod + grad_phi_sigma_i_grad_phi;
00109 
00110     // odd rows, even columns
00111     matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00112     slice10(ret, slice (1, 2, ELEMENT_DIM+1), slice (0, 2, ELEMENT_DIM+1));
00113     slice10 = grad_phi_sigma_i_grad_phi;
00114 
00115     // even rows, odd columns
00116     matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00117     slice01(ret, slice (0, 2, ELEMENT_DIM+1), slice (1, 2, ELEMENT_DIM+1));
00118     slice01 = grad_phi_sigma_i_grad_phi;
00119 
00120     // odd rows, odd columns
00121     matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00122     slice11(ret, slice (1, 2, ELEMENT_DIM+1), slice (1, 2, ELEMENT_DIM+1));
00123     slice11 = grad_phi_sigma_i_grad_phi + grad_phi_sigma_e_grad_phi;
00124 
00125     return ret;
00126 }
00127 
00128 
00129 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00130 c_vector<double,2*(ELEMENT_DIM+1)>
00131     BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeVectorTerm(
00132             c_vector<double, ELEMENT_DIM+1> &rPhi,
00133             c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> &rGradPhi,
00134             ChastePoint<SPACE_DIM> &rX,
00135             c_vector<double,2> &u,
00136             c_matrix<double, 2, SPACE_DIM> &rGradU /* not used */,
00137             Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00138 {
00139     // get bidomain parameters
00140     double Am = mpConfig->GetSurfaceAreaToVolumeRatio();
00141     double Cm = mpConfig->GetCapacitance();
00142 
00143     c_vector<double,2*(ELEMENT_DIM+1)> ret;
00144 
00145     vector_slice<c_vector<double, 2*(ELEMENT_DIM+1)> > slice_V  (ret, slice (0, 2, ELEMENT_DIM+1));
00146     vector_slice<c_vector<double, 2*(ELEMENT_DIM+1)> > slice_Phi(ret, slice (1, 2, ELEMENT_DIM+1));
00147 
00148     // u(0) = voltage
00149     noalias(slice_V)   = (Am*Cm*u(0)/this->mDt - Am*mIionic - mIIntracellularStimulus) * rPhi;
00150     noalias(slice_Phi) = zero_vector<double>(ELEMENT_DIM+1);
00151 
00152     return ret;
00153 }
00154 
00155 
00156 
00157 #define COVERAGE_IGNORE
00158 // never called, since no non-zero neumann conditions
00159 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00160 c_vector<double, 2*ELEMENT_DIM> BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeVectorSurfaceTerm(
00161     const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM> &rSurfaceElement,
00162     c_vector<double,ELEMENT_DIM> &rPhi,
00163     ChastePoint<SPACE_DIM> &rX)
00164 {
00165     // D_times_gradu_dot_n = [D grad(u)].n, D=diffusion matrix
00166     double D_times_grad_v_dot_n     = this->mpBoundaryConditions->GetNeumannBCValue(&rSurfaceElement, rX, 0);
00167     double D_times_grad_phi_e_dot_n = this->mpBoundaryConditions->GetNeumannBCValue(&rSurfaceElement, rX, 1);
00168 
00169     c_vector<double, 2*ELEMENT_DIM> ret;
00170     for (unsigned i=0; i<ELEMENT_DIM; i++)
00171     {
00172         ret(2*i)   = rPhi(i)*D_times_grad_v_dot_n;
00173         ret(2*i+1) = rPhi(i)*D_times_grad_phi_e_dot_n;
00174     }
00175 
00176     return ret;
00177 }
00178 #undef COVERAGE_IGNORE
00179 
00180 
00181 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00182 void BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::PrepareForAssembleSystem(Vec currentSolution, double time)
00183 {
00184     mpBidomainPde->SolveCellSystems(currentSolution, time, time+this->mDt);
00185 }
00186 
00187 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00188 void BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::FinaliseAssembleSystem(Vec currentSolution, double currentTime)
00189 {
00190     if (mFixedExtracellularPotentialNodes.empty())
00191     {
00192         // We're not pinning any nodes.
00193         if (mRowForAverageOfPhiZeroed==INT_MAX)
00194         {
00195             // We're not using the 'Average phi_e = 0' method, hence use a null space
00196             if (!mNullSpaceCreated)
00197             {
00198                 // No null space set up, so create one and pass it to the linear system
00199                 Vec nullbasis[1];
00200                 nullbasis[0]=DistributedVector::CreateVec(2);
00201                 DistributedVector dist_null_basis(nullbasis[0]);
00202                 DistributedVector::Stripe null_basis_stripe_0(dist_null_basis,0);
00203                 DistributedVector::Stripe null_basis_stripe_1(dist_null_basis,1);
00204                 for (DistributedVector::Iterator index = DistributedVector::Begin();
00205                      index != DistributedVector::End();
00206                      ++index)
00207                 {
00208                     null_basis_stripe_0[index] = 0;
00209                     null_basis_stripe_1[index] = 1;
00210                 }
00211                 dist_null_basis.Restore();
00212 
00213                 this->mpLinearSystem->SetNullBasis(nullbasis, 1);
00214 
00215                 VecDestroy(nullbasis[0]);
00216                 mNullSpaceCreated = true;
00217 
00218                 //Make a mask to use if we need to shift the external voltage
00219                 VecDuplicate(currentSolution, &mExternalVoltageMask);
00220                 DistributedVector mask(mExternalVoltageMask);
00221                 DistributedVector::Stripe v_m(mask,0);
00222                 DistributedVector::Stripe phi_e(mask,1);
00223                 for (DistributedVector::Iterator index = DistributedVector::Begin();
00224                      index != DistributedVector::End();
00225                      ++index)
00226                 {
00227                     v_m[index] = 0.0;
00228                     phi_e[index] = 1.0;
00229                 }
00230                 mask.Restore();
00231             }
00232             //Try to fudge the solution vector with respect to the external voltage
00233             //Find the largest absolute value
00234             double min, max;
00235 
00236 #if (PETSC_VERSION_MINOR == 2) //Old API
00237             PetscInt position;
00238             VecMax(currentSolution, &position, &max);
00239             VecMin(currentSolution, &position, &min);
00240 #else
00241             VecMax(currentSolution, PETSC_NULL, &max);
00242             VecMin(currentSolution, PETSC_NULL, &min);
00243 #endif
00244             if ( -min > max )
00245             {
00246                 //Largest value is negative
00247                 max=min;
00248             }
00249             //Standard transmembrane potentials are within +-100 mV
00250             if (fabs(max) > 500)
00251             {
00252 #define COVERAGE_IGNORE
00253                 // std::cout<<"warning: shifting phi_e by "<<-max<<"\n";
00254                 //Use mask currentSolution=currentSolution - max*mExternalVoltageMask
00255 #if (PETSC_VERSION_MINOR == 2) //Old API
00256                 max *= -1;
00257                 VecAXPY(&max, mExternalVoltageMask, currentSolution);
00258 #else
00259                 VecAXPY(currentSolution, -max, mExternalVoltageMask);
00260 #endif
00261 #undef COVERAGE_IGNORE
00262             }
00263         }
00264         else
00265         {
00266             // mRowForAverageOfPhiZeroed!=INT_MAX, i.e. we're using the 'Average phi_e = 0; method
00267             // Set average phi_e to zero
00268             unsigned matrix_size = this->mpLinearSystem->GetSize();
00269             if (!this->mMatrixIsAssembled)
00270             {
00271 
00272                 // Set the mRowForAverageOfPhiZeroed-th matrix row to 0 1 0 1 ...
00273                 this->mpLinearSystem->ZeroMatrixRow(mRowForAverageOfPhiZeroed);
00274                 for (unsigned col_index=0; col_index<matrix_size; col_index++)
00275                 {
00276                     if (col_index%2 == 1)
00277                     {
00278                         this->mpLinearSystem->SetMatrixElement(mRowForAverageOfPhiZeroed, col_index, 1);
00279                     }
00280                 }
00281                 this->mpLinearSystem->AssembleFinalLhsMatrix();
00282 
00283             }
00284             // Set the mRowForAverageOfPhiZeroed-th rhs vector row to 0
00285             this->mpLinearSystem->SetRhsVectorElement(mRowForAverageOfPhiZeroed, 0);
00286 
00287             this->mpLinearSystem->AssembleRhsVector();
00288         }
00289     }
00290 }
00291 
00292 
00293 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00294 BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::BidomainDg0Assembler(
00295             AbstractMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00296             BidomainPde<SPACE_DIM>* pPde,
00297             BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, 2>* pBcc,
00298             unsigned numQuadPoints)
00299     : AbstractAssembler<ELEMENT_DIM,SPACE_DIM,2>(),
00300       BaseClassType(numQuadPoints),
00301       AbstractDynamicAssemblerMixin<ELEMENT_DIM,SPACE_DIM,2>()
00302 {
00303     assert(pPde != NULL);
00304     assert(pMesh != NULL);
00305     assert(pBcc != NULL);
00306 
00307     mpBidomainPde = pPde;
00308     this->SetMesh(pMesh);
00309 
00310     this->mpBoundaryConditions = pBcc;
00311 
00312     mNullSpaceCreated = false;
00313 
00314     this->SetMatrixIsConstant();
00315 
00316     mRowForAverageOfPhiZeroed = INT_MAX; //this->mpLinearSystem->GetSize() - 1;
00317     
00318     mpConfig = HeartConfig::Instance();
00319 }
00320 
00321 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00322 BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::~BidomainDg0Assembler()
00323 {
00324     if (mNullSpaceCreated)
00325     {
00326         VecDestroy(mExternalVoltageMask);
00327     }
00328 }
00329 
00330 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00331 void BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::SetFixedExtracellularPotentialNodes(
00332             std::vector<unsigned> fixedExtracellularPotentialNodes)
00333 {
00334     for (unsigned i=0; i<fixedExtracellularPotentialNodes.size(); i++)
00335     {
00336         if (fixedExtracellularPotentialNodes[i] >= this->mpMesh->GetNumNodes() )
00337         {
00338             EXCEPTION("Fixed node number must be less than total number nodes");
00339         }
00340     }
00341 
00342     mFixedExtracellularPotentialNodes = fixedExtracellularPotentialNodes;
00343 
00344     for (unsigned i=0; i<mFixedExtracellularPotentialNodes.size(); i++)
00345     {
00346         ConstBoundaryCondition<SPACE_DIM>* p_boundary_condition
00347          = new ConstBoundaryCondition<SPACE_DIM>(0.0);
00348 
00349         Node<SPACE_DIM>* p_node = this->mpMesh->GetNode(mFixedExtracellularPotentialNodes[i]);
00350 
00351         this->mpBoundaryConditions->AddDirichletBoundaryCondition(p_node, p_boundary_condition, 1);
00352     }
00353 }
00354 
00355 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00356 void BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::SetRowForAverageOfPhiZeroed(unsigned row)
00357 {
00358     // Row should be odd in C++-like indexing
00359     if (row % 2 == 0)
00360     {
00361         EXCEPTION("Row for applying the constraint 'Average of phi_e = zero' should be odd in C++ like indexing");
00362     }
00363 
00364     mRowForAverageOfPhiZeroed = row;
00365 }
00366 
00368 // Explicit instantiation
00370 
00371 template class BidomainDg0Assembler<1,1>;
00372 template class BidomainDg0Assembler<2,2>;
00373 template class BidomainDg0Assembler<3,3>;

Generated on Wed Mar 18 12:51:52 2009 for Chaste by  doxygen 1.5.5