BidomainWithBathAssembler.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 "BidomainWithBathAssembler.hpp"
00030 
00031 #include <boost/numeric/ublas/vector_proxy.hpp>
00032 
00033 #include "ConstBoundaryCondition.hpp"
00034 #include "HeartRegionCodes.hpp"
00035 
00036 
00037 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00038 c_matrix<double,2*(ELEMENT_DIM+1),2*(ELEMENT_DIM+1)>
00039     BidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>::ComputeMatrixTerm(
00040             c_vector<double, ELEMENT_DIM+1> &rPhi,
00041             c_matrix<double, SPACE_DIM, ELEMENT_DIM+1> &rGradPhi,
00042             ChastePoint<SPACE_DIM> &rX,
00043             c_vector<double,2> &rU,
00044             c_matrix<double, 2, SPACE_DIM> &rGradU /* not used */,
00045             Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00046 {
00047     if (pElement->GetRegion() != HeartRegionCode::BATH) // ie if a tissue element
00048     {
00049         return BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeMatrixTerm(rPhi,rGradPhi,rX,rU,rGradU,pElement);
00050     }
00051     else // bath element
00052     {
00053         double bath_cond=HeartConfig::Instance()->GetBathConductivity();
00054         const c_matrix<double, SPACE_DIM, SPACE_DIM>& sigma_b = bath_cond*identity_matrix<double>(SPACE_DIM);
00055 
00056         c_matrix<double, SPACE_DIM, ELEMENT_DIM+1> temp = prod(sigma_b, rGradPhi);
00057         c_matrix<double, ELEMENT_DIM+1, ELEMENT_DIM+1> grad_phi_sigma_b_grad_phi =
00058             prod(trans(rGradPhi), temp);
00059 
00060         c_matrix<double,2*(ELEMENT_DIM+1),2*(ELEMENT_DIM+1)> ret = zero_matrix<double>(2*(ELEMENT_DIM+1));
00061 
00062         // even rows, even columns
00063         //matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00064         //slice00(ret, slice (0, 2, ELEMENT_DIM+1), slice (0, 2, ELEMENT_DIM+1));
00065         //slice00 = 0;
00066 
00067         // odd rows, even columns
00068         //matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00069         //slice10(ret, slice (1, 2, ELEMENT_DIM+1), slice (0, 2, ELEMENT_DIM+1));
00070         //slice10 = 0
00071 
00072         // even rows, odd columns
00073         //matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00074         //slice01(ret, slice (0, 2, ELEMENT_DIM+1), slice (1, 2, ELEMENT_DIM+1));
00075         //slice01 = 0;
00076 
00077         // odd rows, odd columns
00078         matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00079         slice11(ret, slice (1, 2, ELEMENT_DIM+1), slice (1, 2, ELEMENT_DIM+1));
00080         slice11 = grad_phi_sigma_b_grad_phi;
00081 
00082         return ret;
00083     }
00084 }
00085 
00086 
00087 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00088 c_vector<double,2*(ELEMENT_DIM+1)>
00089     BidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>::ComputeVectorTerm(
00090             c_vector<double, ELEMENT_DIM+1> &rPhi,
00091             c_matrix<double, SPACE_DIM, ELEMENT_DIM+1> &rGradPhi,
00092             ChastePoint<SPACE_DIM> &rX,
00093             c_vector<double,2> &rU,
00094             c_matrix<double, 2, SPACE_DIM> &rGradU /* not used */,
00095             Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00096 {
00097     if (pElement->GetRegion() != HeartRegionCode::BATH) // ie if a tissue element
00098     {
00099         return BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeVectorTerm(rPhi,rGradPhi,rX,rU,rGradU,pElement);
00100     }
00101     else // bath element
00102     {
00103         c_vector<double,2*(ELEMENT_DIM+1)> ret = zero_vector<double>(2*(ELEMENT_DIM+1));
00104 
00105         //vector_slice<c_vector<double, 2*(ELEMENT_DIM+1)> > slice_V  (ret, slice (0, 2, ELEMENT_DIM+1));
00106         //vector_slice<c_vector<double, 2*(ELEMENT_DIM+1)> > slice_Phi(ret, slice (1, 2, ELEMENT_DIM+1));
00107 
00108         // u(0) = voltage
00109         //noalias(slice_V) = zero_vector<double>(ELEMENT_DIM+1);
00110         //noalias(slice_Phi) = zero_vector<double>(ELEMENT_DIM+1);
00111 
00112         return ret;
00113     }
00114 }
00115 
00116 
00117 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00118 void BidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>::FinaliseLinearSystem(
00119             Vec existingSolutionOrGuess,
00120             double time,
00121             bool assembleVector, bool assembleMatrix)
00122 {
00124     // CG (default solver) won't work since the system is indefinite. Switch to SYMMLQ
00125 //    this->mpLinearSystem->SetKspType("symmlq"); // Switches the solver
00126 //    this->mpConfig->SetKSPSolver("symmlq"); // Makes sure this change will be reflected in the XML file written to disk at the end of the simulation.
00127 
00128     unsigned* is_node_bath = new unsigned[this->mpMesh->GetNumNodes()];
00129     for(unsigned i = 0; i < this->mpMesh->GetNumNodes(); ++i)
00130     {
00131         is_node_bath[i] = 0;
00132     }
00133 
00134     for (typename AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>::NodeIterator iter=this->mpMesh->GetNodeIteratorBegin();
00135          iter != this->mpMesh->GetNodeIteratorEnd();
00136          ++iter)
00137     {
00143         if ((*iter).GetRegion() == HeartRegionCode::BATH)
00144         {
00145             is_node_bath[(*iter).GetIndex()] = 1;
00146         }
00147     }
00148 
00149     unsigned* is_node_bath_reduced = new unsigned[this->mpMesh->GetNumNodes()];
00150     MPI_Allreduce(is_node_bath, is_node_bath_reduced, this->mpMesh->GetNumNodes(), MPI_UNSIGNED, MPI_SUM, PETSC_COMM_WORLD);
00151 
00152     for(unsigned i=0; i<this->mpMesh->GetNumNodes(); i++)
00153     {
00154         if(is_node_bath_reduced[i] > 0) // ie is a bath node
00155         {
00156             PetscInt index[1];
00157             index[0] = 2*i; // assumes Vm and Phie are interleaved
00158 
00159             if(assembleMatrix)
00160             {
00161                 /*
00162                  *  Before revision 6516, we used to zero out i-th row and column here. It seems to be redundant because they are already zero after assembly.
00163                  *  When assembling a bath element you get a matrix subblock that looks like (2D example):
00164                  *
00165                  *     Vm   0 0 0 0 0 0
00166                  *     Vm   0 0 0 0 0 0
00167                  *     Vm   0 0 0 0 0 0
00168                  *     Phie 0 0 0 x x x
00169                  *     Phie 0 0 0 x x x  -> the x subblock is assembled from div(grad_phi) = 0
00170                  *     Phie 0 0 0 x x x
00171                  *
00172                  *  Therefore, all the Vm entries of this node are already 0.
00173                  *
00174                  *  Explicitly checking it in non-production builds.
00175                  */
00176 #ifndef NDEBUG
00177                 int num_equation = 2*i; // assumes Vm and Phie are interleaved
00178 
00179                 // Matrix need to be assembled in order to use GetMatrixElement()
00180                 MatAssemblyBegin((*this->GetLinearSystem())->rGetLhsMatrix(), MAT_FINAL_ASSEMBLY);
00181                 MatAssemblyEnd((*this->GetLinearSystem())->rGetLhsMatrix(), MAT_FINAL_ASSEMBLY);
00182 
00183                 PetscInt local_lo, local_hi;
00184                 (*this->GetLinearSystem())->GetOwnershipRange(local_lo, local_hi);
00185 
00186                 // If this processor owns i-th row, check it.
00187                 if ((local_lo <= (int)num_equation) && ((int)num_equation < local_hi))
00188                 {
00189                     for (unsigned column=0; column < (*this->GetLinearSystem())->GetSize(); column++)
00190                     {
00191                         assert((*this->GetLinearSystem())->GetMatrixElement(num_equation, column)==0.0);
00192                     }
00193                 }
00194 
00195                 // Check the local entries of the i-th column
00196                 for (int row=local_lo; row<local_hi; row++)
00197                 {
00198                     assert((*this->GetLinearSystem())->GetMatrixElement(row, num_equation)==0);
00199                 }
00200 #endif
00201                 // put 1.0 on the diagonal
00202                 Mat& r_matrix = (*(this->GetLinearSystem()))->rGetLhsMatrix();
00203                 MatSetValue(r_matrix,index[0],index[0],1.0,INSERT_VALUES);
00204             }
00205 
00206             if(assembleVector)
00207             {
00208                 // zero rhs vector entry
00209                 VecSetValue((*(this->GetLinearSystem()))->rGetRhsVector(), index[0], 0.0, INSERT_VALUES);
00210             }
00211         }
00212     }
00213 
00214     delete[] is_node_bath;
00215     delete[] is_node_bath_reduced;
00216 }
00217 
00218 
00219 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00220 BidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>::BidomainWithBathAssembler(
00221             AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00222             BidomainPde<SPACE_DIM>* pPde,
00223             BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, 2>* pBcc,
00224             unsigned numQuadPoints)
00225     : BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>(pMesh, pPde, pBcc, numQuadPoints)
00226 {
00227 }
00228 
00230 // Explicit instantiation
00232 
00233 template class BidomainWithBathAssembler<1,1>;
00234 template class BidomainWithBathAssembler<2,2>;
00235 template class BidomainWithBathAssembler<3,3>;

Generated by  doxygen 1.6.2