00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
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 ,
00045 Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00046 {
00047 if (pElement->GetRegion() == HeartRegionCode::TISSUE)
00048 {
00049 return BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeMatrixTerm(rPhi,rGradPhi,rX,u,rGradU,pElement);
00050 }
00051 else
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
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
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 ,
00097 Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00098 {
00099 if (pElement->GetRegion() == HeartRegionCode::TISSUE)
00100 {
00101 return BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeVectorTerm(rPhi,rGradPhi,rX,u,rGradU,pElement);
00102 }
00103 else
00104 {
00105 c_vector<double,2*(ELEMENT_DIM+1)> ret = zero_vector<double>(2*(ELEMENT_DIM+1));
00106
00107
00108
00109
00110
00111
00112
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)
00128 {
00129 PetscInt index[1];
00130 index[0] = 2*i;
00131
00132 if(assembleMatrix)
00133 {
00134
00135 (*(this->GetLinearSystem()))->ZeroMatrixRow(2*i);
00136
00137 (*(this->GetLinearSystem()))->ZeroMatrixColumn(2*i);
00138
00139
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
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
00167
00168 template class BidomainWithBathAssembler<1,1>;
00169 template class BidomainWithBathAssembler<2,2>;
00170 template class BidomainWithBathAssembler<3,3>;