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, 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 ,
00045 Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00046 {
00047 if (pElement->GetRegion() != HeartRegionCode::BATH)
00048 {
00049 return BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeMatrixTerm(rPhi,rGradPhi,rX,rU,rGradU,pElement);
00050 }
00051 else
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
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
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 ,
00095 Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00096 {
00097 if (pElement->GetRegion() != HeartRegionCode::BATH)
00098 {
00099 return BidomainDg0Assembler<ELEMENT_DIM,SPACE_DIM>::ComputeVectorTerm(rPhi,rGradPhi,rX,rU,rGradU,pElement);
00100 }
00101 else
00102 {
00103 c_vector<double,2*(ELEMENT_DIM+1)> ret = zero_vector<double>(2*(ELEMENT_DIM+1));
00104
00105
00106
00107
00108
00109
00110
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
00125
00126
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)
00155 {
00156 PetscInt index[1];
00157 index[0] = 2*i;
00158
00159 if(assembleMatrix)
00160 {
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174
00175
00176 #ifndef NDEBUG
00177 int num_equation = 2*i;
00178
00179
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
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
00196 for (int row=local_lo; row<local_hi; row++)
00197 {
00198 assert((*this->GetLinearSystem())->GetMatrixElement(row, num_equation)==0);
00199 }
00200 #endif
00201
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
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
00232
00233 template class BidomainWithBathAssembler<1,1>;
00234 template class BidomainWithBathAssembler<2,2>;
00235 template class BidomainWithBathAssembler<3,3>;