BidomainWithBathAssembler.cpp

00001 
00002 /*
00003 
00004 Copyright (C) University of Oxford, 2005-2011
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Chaste is free software: you can redistribute it and/or modify it
00013 under the terms of the GNU Lesser General Public License as published
00014 by the Free Software Foundation, either version 2.1 of the License, or
00015 (at your option) any later version.
00016 
00017 Chaste is distributed in the hope that it will be useful, but WITHOUT
00018 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00019 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00020 License for more details. The offer of Chaste under the terms of the
00021 License is subject to the License being interpreted in accordance with
00022 English Law and subject to any action against the University of Oxford
00023 being under the jurisdiction of the English Courts.
00024 
00025 You should have received a copy of the GNU Lesser General Public License
00026 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00027 
00028 */
00029 
00030 #include "BidomainWithBathAssembler.hpp"
00031 
00032 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00033 c_matrix<double,2*(ELEMENT_DIM+1),2*(ELEMENT_DIM+1)>
00034     BidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>::ComputeMatrixTerm(
00035             c_vector<double, ELEMENT_DIM+1> &rPhi,
00036             c_matrix<double, SPACE_DIM, ELEMENT_DIM+1> &rGradPhi,
00037             ChastePoint<SPACE_DIM> &rX,
00038             c_vector<double,2> &rU,
00039             c_matrix<double, 2, SPACE_DIM> &rGradU /* not used */,
00040             Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00041 {
00042     if (!HeartRegionCode::IsRegionBath( pElement->GetRegion() )) // ie if a tissue element
00043     {
00044         return BidomainAssembler<ELEMENT_DIM,SPACE_DIM>::ComputeMatrixTerm(rPhi,rGradPhi,rX,rU,rGradU,pElement);
00045     }
00046     else // bath element
00047     {
00048         double bath_cond=HeartConfig::Instance()->GetBathConductivity(pElement->GetRegion());
00049 
00050         c_matrix<double, ELEMENT_DIM+1, ELEMENT_DIM+1> grad_phi_sigma_b_grad_phi =
00051             bath_cond * prod(trans(rGradPhi), rGradPhi);
00052 
00053         c_matrix<double,2*(ELEMENT_DIM+1),2*(ELEMENT_DIM+1)> ret = zero_matrix<double>(2*(ELEMENT_DIM+1));
00054 
00055         // even rows, even columns
00056         //matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00057         //slice00(ret, slice (0, 2, ELEMENT_DIM+1), slice (0, 2, ELEMENT_DIM+1));
00058         //slice00 = 0;
00059 
00060         // odd rows, even columns
00061         //matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00062         //slice10(ret, slice (1, 2, ELEMENT_DIM+1), slice (0, 2, ELEMENT_DIM+1));
00063         //slice10 = 0
00064 
00065         // even rows, odd columns
00066         //matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00067         //slice01(ret, slice (0, 2, ELEMENT_DIM+1), slice (1, 2, ELEMENT_DIM+1));
00068         //slice01 = 0;
00069 
00070         // odd rows, odd columns
00071         matrix_slice<c_matrix<double, 2*ELEMENT_DIM+2, 2*ELEMENT_DIM+2> >
00072         slice11(ret, slice (1, 2, ELEMENT_DIM+1), slice (1, 2, ELEMENT_DIM+1));
00073         slice11 = grad_phi_sigma_b_grad_phi;
00074 
00075         return ret;
00076     }
00077 }
00078 
00079 
00080 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00081 c_vector<double,2*(ELEMENT_DIM+1)>
00082     BidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>::ComputeVectorTerm(
00083             c_vector<double, ELEMENT_DIM+1> &rPhi,
00084             c_matrix<double, SPACE_DIM, ELEMENT_DIM+1> &rGradPhi,
00085             ChastePoint<SPACE_DIM> &rX,
00086             c_vector<double,2> &rU,
00087             c_matrix<double, 2, SPACE_DIM> &rGradU /* not used */,
00088             Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00089 {
00090     if (HeartRegionCode::IsRegionTissue( pElement->GetRegion() )) // ie if a tissue element
00091     {
00092         return BidomainAssembler<ELEMENT_DIM,SPACE_DIM>::ComputeVectorTerm(rPhi,rGradPhi,rX,rU,rGradU,pElement);
00093     }
00094     else // bath element
00095     {
00096         return zero_vector<double>(2*(ELEMENT_DIM+1));
00097     }
00098 }
00099 
00100 
00102 // explicit instantiation
00104 
00105 template class BidomainWithBathAssembler<1,1>;
00106 template class BidomainWithBathAssembler<2,2>;
00107 template class BidomainWithBathAssembler<3,3>;

Generated on Tue May 31 14:31:41 2011 for Chaste by  doxygen 1.5.5