BidomainWithBathAssembler.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 "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, ELEMENT_DIM, ELEMENT_DIM+1> &rGradPhi,
00042             ChastePoint<SPACE_DIM> &rX,
00043             c_vector<double,2> &u,
00044             c_matrix<double, 2, SPACE_DIM> &rGradU /* not used */,
00045             Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00046 {
00047     if (pElement->GetRegion() == HeartRegionCode::TISSUE) // ie if a tissue element
00048     {
00049         return BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeMatrixTerm(rPhi,rGradPhi,rX,u,rGradU,pElement);
00050     }
00051     else // bath element
00052     {
00053          
00055         double bath_cond=HeartConfig::Instance()->GetBathConductivity(); 
00056         const c_matrix<double, SPACE_DIM, SPACE_DIM>& sigma_b = bath_cond*identity_matrix<double>(SPACE_DIM);
00057 
00058         c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> temp = prod(sigma_b, rGradPhi);
00059         c_matrix<double, ELEMENT_DIM+1, ELEMENT_DIM+1> grad_phi_sigma_b_grad_phi =
00060             prod(trans(rGradPhi), temp);
00061 
00062         c_matrix<double,2*(ELEMENT_DIM+1),2*(ELEMENT_DIM+1)> ret = zero_matrix<double>(2*(ELEMENT_DIM+1));
00063 
00064         // even rows, even columns
00065         //matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00066         //slice00(ret, slice (0, 2, ELEMENT_DIM+1), slice (0, 2, ELEMENT_DIM+1));
00067         //slice00 = 0;
00068 
00069         // odd rows, even columns
00070         //matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00071         //slice10(ret, slice (1, 2, ELEMENT_DIM+1), slice (0, 2, ELEMENT_DIM+1));
00072         //slice10 = 0
00073 
00074         // even rows, odd columns
00075         //matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00076         //slice01(ret, slice (0, 2, ELEMENT_DIM+1), slice (1, 2, ELEMENT_DIM+1));
00077         //slice01 = 0;
00078 
00079         // odd rows, odd columns
00080         matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00081         slice11(ret, slice (1, 2, ELEMENT_DIM+1), slice (1, 2, ELEMENT_DIM+1));
00082         slice11 = grad_phi_sigma_b_grad_phi;
00083 
00084         return ret;
00085     }
00086 }
00087 
00088 
00089 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00090 c_vector<double,2*(ELEMENT_DIM+1)>
00091     BidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>::ComputeVectorTerm(
00092             c_vector<double, ELEMENT_DIM+1> &rPhi,
00093             c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> &rGradPhi,
00094             ChastePoint<SPACE_DIM> &rX,
00095             c_vector<double,2> &u,
00096             c_matrix<double, 2, SPACE_DIM> &rGradU /* not used */,
00097             Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00098 {
00099     if (pElement->GetRegion() == HeartRegionCode::TISSUE) // ie if a tissue element
00100     {
00101         return BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeVectorTerm(rPhi,rGradPhi,rX,u,rGradU,pElement);
00102     }
00103     else // bath element
00104     {
00105         c_vector<double,2*(ELEMENT_DIM+1)> ret = zero_vector<double>(2*(ELEMENT_DIM+1));
00106 
00107         //vector_slice<c_vector<double, 2*(ELEMENT_DIM+1)> > slice_V  (ret, slice (0, 2, ELEMENT_DIM+1));
00108         //vector_slice<c_vector<double, 2*(ELEMENT_DIM+1)> > slice_Phi(ret, slice (1, 2, ELEMENT_DIM+1));
00109 
00110         // u(0) = voltage
00111         //noalias(slice_V) = zero_vector<double>(ELEMENT_DIM+1); 
00112         //noalias(slice_Phi) = zero_vector<double>(ELEMENT_DIM+1);
00113 
00114         return ret;
00115     }
00116 }
00117 
00118 
00119 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00120 void BidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>::FinaliseLinearSystem(
00121             Vec currentSolutionOrGuess,
00122             double currentTime,
00123             bool assembleVector, bool assembleMatrix)
00124 {
00125     for(unsigned i=0; i<this->mpMesh->GetNumNodes(); i++)
00126     {
00127         if(this->mpMesh->GetNode(i)->GetRegion() == HeartRegionCode::BATH) // ie is a bath node
00128         {
00129             PetscInt index[1];
00130             index[0] = 2*i;
00131 
00132             if(assembleMatrix)
00133             {
00134                 // zero the row corresponding to V for this bath node
00135                 (*(this->GetLinearSystem()))->ZeroMatrixRow(2*i);
00136                 // zero the column corresponding to V for this bath node.
00137                 (*(this->GetLinearSystem()))->ZeroMatrixColumn(2*i);
00138 
00139                 // put 1.0 on the diagonal
00140                 Mat& r_matrix = (*(this->GetLinearSystem()))->rGetLhsMatrix();
00141                 MatSetValue(r_matrix,index[0],index[0],1.0,INSERT_VALUES);
00142             }
00143             
00144             if(assembleVector)
00145             {
00146                 // zero rhs vector entry
00147                 VecSetValue((*(this->GetLinearSystem()))->rGetRhsVector(), index[0], 0.0, INSERT_VALUES);
00148             }
00149         }
00150     }
00151 }
00152 
00153 
00154 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00155 BidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>::BidomainWithBathAssembler(
00156             AbstractMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00157             BidomainPde<SPACE_DIM>* pPde,
00158             BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, 2>* pBcc,
00159             unsigned numQuadPoints)
00160     : BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>(pMesh, pPde, pBcc, numQuadPoints)
00161 {        
00162 }
00163 
00165 // Explicit instantiation
00167 
00168 template class BidomainWithBathAssembler<1,1>;
00169 template class BidomainWithBathAssembler<2,2>;
00170 template class BidomainWithBathAssembler<3,3>;

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