Chaste Release::3.1
|
00001 /* 00002 00003 Copyright (c) 2005-2012, University of Oxford. 00004 All rights reserved. 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 Redistribution and use in source and binary forms, with or without 00013 modification, are permitted provided that the following conditions are met: 00014 * Redistributions of source code must retain the above copyright notice, 00015 this list of conditions and the following disclaimer. 00016 * Redistributions in binary form must reproduce the above copyright notice, 00017 this list of conditions and the following disclaimer in the documentation 00018 and/or other materials provided with the distribution. 00019 * Neither the name of the University of Oxford nor the names of its 00020 contributors may be used to endorse or promote products derived from this 00021 software without specific prior written permission. 00022 00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" 00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE 00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE 00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR 00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE 00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) 00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT 00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT 00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00033 00034 */ 00035 00036 #ifndef ABSTRACTFECABLEINTEGRALASSEMBLER_HPP_ 00037 #define ABSTRACTFECABLEINTEGRALASSEMBLER_HPP_ 00038 00039 #include "AbstractFeAssemblerCommon.hpp" 00040 #include "GaussianQuadratureRule.hpp" 00041 #include "PetscVecTools.hpp" 00042 #include "PetscMatTools.hpp" 00043 00050 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL> 00051 class AbstractFeCableIntegralAssembler : public AbstractFeAssemblerCommon<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,CAN_ASSEMBLE_VECTOR,CAN_ASSEMBLE_MATRIX,INTERPOLATION_LEVEL> 00052 { 00053 protected: 00055 static const unsigned CABLE_ELEMENT_DIM = 1; 00056 00058 static const unsigned NUM_CABLE_ELEMENT_NODES = 2; 00059 00061 MixedDimensionMesh<ELEMENT_DIM, SPACE_DIM>* mpMesh; 00062 00064 GaussianQuadratureRule<1>* mpCableQuadRule; 00065 00067 typedef LinearBasisFunction<1> CableBasisFunction; 00068 00083 void ComputeTransformedBasisFunctionDerivatives(const ChastePoint<CABLE_ELEMENT_DIM>& rPoint, 00084 const c_matrix<double, CABLE_ELEMENT_DIM, SPACE_DIM>& rInverseJacobian, 00085 c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rReturnValue); 00086 00092 void DoAssemble(); 00093 00116 virtual c_matrix<double,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES> ComputeCableMatrixTerm( 00117 c_vector<double, NUM_CABLE_ELEMENT_NODES>& rPhi, 00118 c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rGradPhi, 00119 ChastePoint<SPACE_DIM>& rX, 00120 c_vector<double,PROBLEM_DIM>& rU, 00121 c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU, 00122 Element<CABLE_ELEMENT_DIM,SPACE_DIM>* pElement) 00123 { 00124 // If this line is reached this means this method probably hasn't been over-ridden correctly in 00125 // the concrete class 00126 NEVER_REACHED; 00127 return zero_matrix<double>(PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES); 00128 } 00129 00150 virtual c_vector<double,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES> ComputeCableVectorTerm( 00151 c_vector<double, NUM_CABLE_ELEMENT_NODES>& rPhi, 00152 c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rGradPhi, 00153 ChastePoint<SPACE_DIM>& rX, 00154 c_vector<double,PROBLEM_DIM>& rU, 00155 c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU, 00156 Element<CABLE_ELEMENT_DIM,SPACE_DIM>* pElement) 00157 { 00158 // If this line is reached this means this method probably hasn't been over-ridden correctly in 00159 // the concrete class 00160 NEVER_REACHED; 00161 return zero_vector<double>(PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES); 00162 } 00163 00178 virtual void AssembleOnCableElement(Element<CABLE_ELEMENT_DIM,SPACE_DIM>& rElement, 00179 c_matrix<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES >& rAElem, 00180 c_vector<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES>& rBElem); 00181 00182 00190 virtual bool ElementAssemblyCriterion(Element<CABLE_ELEMENT_DIM,SPACE_DIM>& rElement) 00191 { 00192 return true; 00193 } 00194 00195 public: 00196 00204 AbstractFeCableIntegralAssembler(MixedDimensionMesh<ELEMENT_DIM,SPACE_DIM>* pMesh, unsigned numQuadPoints=2); 00205 00209 virtual ~AbstractFeCableIntegralAssembler() 00210 { 00211 delete mpCableQuadRule; 00212 } 00213 }; 00214 00215 00216 00217 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL> 00218 AbstractFeCableIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::AbstractFeCableIntegralAssembler( 00219 MixedDimensionMesh<ELEMENT_DIM,SPACE_DIM>* pMesh, unsigned numQuadPoints) 00220 : AbstractFeAssemblerCommon<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>(), 00221 mpMesh(pMesh) 00222 { 00223 assert(pMesh); 00224 assert(numQuadPoints > 0); 00225 assert(CAN_ASSEMBLE_VECTOR || CAN_ASSEMBLE_MATRIX); 00226 00227 mpCableQuadRule = new GaussianQuadratureRule<CABLE_ELEMENT_DIM>(numQuadPoints); 00228 00229 // Not supporting this yet - if a nonlinear assembler on cable elements is required, uncomment code 00230 // in AssembleOnCableElement below (search for NONLINEAR) 00231 assert(INTERPOLATION_LEVEL!=NONLINEAR); 00232 } 00233 00234 00235 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL> 00236 void AbstractFeCableIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::DoAssemble() 00237 { 00238 HeartEventHandler::EventType assemble_event; 00239 if (this->mAssembleMatrix) 00240 { 00241 assemble_event = HeartEventHandler::ASSEMBLE_SYSTEM; 00242 } 00243 else 00244 { 00245 assemble_event = HeartEventHandler::ASSEMBLE_RHS; 00246 } 00247 00248 if (this->mAssembleMatrix && this->mMatrixToAssemble==NULL) 00249 { 00250 EXCEPTION("Matrix to be assembled has not been set"); 00251 } 00252 if (this->mAssembleVector && this->mVectorToAssemble==NULL) 00253 { 00254 EXCEPTION("Vector to be assembled has not been set"); 00255 } 00256 00257 // This has to be below PrepareForAssembleSystem as in that method the ODEs are solved in cardiac problems 00258 HeartEventHandler::BeginEvent(assemble_event); 00259 00260 // Zero the matrix/vector if it is to be assembled 00261 if (this->mAssembleVector && this->mZeroVectorBeforeAssembly) 00262 { 00263 PetscVecTools::Zero(this->mVectorToAssemble); 00264 } 00265 if (this->mAssembleMatrix && this->mZeroMatrixBeforeAssembly) 00266 { 00267 PetscMatTools::Zero(this->mMatrixToAssemble); 00268 } 00269 00270 const size_t STENCIL_SIZE=PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES; 00271 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem; 00272 c_vector<double, STENCIL_SIZE> b_elem; 00273 00274 // Loop over elements 00275 if (this->mAssembleMatrix || this->mAssembleVector) 00276 { 00277 for (typename MixedDimensionMesh<CABLE_ELEMENT_DIM, SPACE_DIM>::CableElementIterator iter = mpMesh->GetCableElementIteratorBegin(); 00278 iter != mpMesh->GetCableElementIteratorEnd(); 00279 ++iter) 00280 { 00281 Element<CABLE_ELEMENT_DIM, SPACE_DIM>& r_element = *(*iter); 00282 00283 // Test for ownership first, since it's pointless to test the criterion on something which we might know nothing about. 00284 if ( r_element.GetOwnership() == true && ElementAssemblyCriterion(r_element)==true ) 00285 { 00286 AssembleOnCableElement(r_element, a_elem, b_elem); 00287 00288 unsigned p_indices[STENCIL_SIZE]; 00289 r_element.GetStiffnessMatrixGlobalIndices(PROBLEM_DIM, p_indices); 00290 00291 if (this->mAssembleMatrix) 00292 { 00293 PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->mMatrixToAssemble, p_indices, a_elem); 00294 } 00295 00296 if (this->mAssembleVector) 00297 { 00298 PetscVecTools::AddMultipleValues<STENCIL_SIZE>(this->mVectorToAssemble, p_indices, b_elem); 00299 } 00300 } 00301 } 00302 } 00303 00304 HeartEventHandler::EndEvent(assemble_event); 00305 } 00306 00308 // Implementation - AssembleOnCableElement and smaller 00310 00311 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL> 00312 void AbstractFeCableIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::ComputeTransformedBasisFunctionDerivatives( 00313 const ChastePoint<CABLE_ELEMENT_DIM>& rPoint, 00314 const c_matrix<double, CABLE_ELEMENT_DIM, SPACE_DIM>& rInverseJacobian, 00315 c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rReturnValue) 00316 { 00317 static c_matrix<double, CABLE_ELEMENT_DIM, NUM_CABLE_ELEMENT_NODES> grad_phi; 00318 00319 LinearBasisFunction<CABLE_ELEMENT_DIM>::ComputeBasisFunctionDerivatives(rPoint, grad_phi); 00320 rReturnValue = prod(trans(rInverseJacobian), grad_phi); 00321 } 00322 00323 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL> 00324 void AbstractFeCableIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::AssembleOnCableElement( 00325 Element<CABLE_ELEMENT_DIM,SPACE_DIM>& rElement, 00326 c_matrix<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES >& rAElem, 00327 c_vector<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES>& rBElem) 00328 { 00334 c_matrix<double, SPACE_DIM, CABLE_ELEMENT_DIM> jacobian; 00335 c_matrix<double, CABLE_ELEMENT_DIM, SPACE_DIM> inverse_jacobian; 00336 double jacobian_determinant; 00337 00340 // mpMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian); 00342 rElement.CalculateInverseJacobian(jacobian, jacobian_determinant, inverse_jacobian); 00343 00344 if (this->mAssembleMatrix) 00345 { 00346 rAElem.clear(); 00347 } 00348 00349 if (this->mAssembleVector) 00350 { 00351 rBElem.clear(); 00352 } 00353 00354 const unsigned num_nodes = rElement.GetNumNodes(); 00355 00356 // Allocate memory for the basis functions values and derivative values 00357 c_vector<double, NUM_CABLE_ELEMENT_NODES> phi; 00358 c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES> grad_phi; 00359 00360 // Loop over Gauss points 00361 for (unsigned quad_index=0; quad_index < mpCableQuadRule->GetNumQuadPoints(); quad_index++) 00362 { 00363 const ChastePoint<CABLE_ELEMENT_DIM>& quad_point = mpCableQuadRule->rGetQuadPoint(quad_index); 00364 00365 CableBasisFunction::ComputeBasisFunctions(quad_point, phi); 00366 00367 if ( this->mAssembleMatrix || INTERPOLATION_LEVEL==NONLINEAR ) 00368 { 00369 ComputeTransformedBasisFunctionDerivatives(quad_point, inverse_jacobian, grad_phi); 00370 } 00371 00372 // Location of the Gauss point in the original element will be stored in x 00373 // Where applicable, u will be set to the value of the current solution at x 00374 ChastePoint<SPACE_DIM> x(0,0,0); 00375 00376 c_vector<double,PROBLEM_DIM> u = zero_vector<double>(PROBLEM_DIM); 00377 c_matrix<double,PROBLEM_DIM,SPACE_DIM> grad_u = zero_matrix<double>(PROBLEM_DIM,SPACE_DIM); 00378 00379 // Allow the concrete version of the assembler to interpolate any desired quantities 00380 this->ResetInterpolatedQuantities(); 00381 00382 // Interpolation 00383 for (unsigned i=0; i<num_nodes; i++) 00384 { 00385 const Node<SPACE_DIM>* p_node = rElement.GetNode(i); 00386 00387 if (INTERPOLATION_LEVEL != CARDIAC) // don't interpolate X if cardiac problem 00388 { 00389 const c_vector<double, SPACE_DIM>& r_node_loc = p_node->rGetLocation(); 00390 // interpolate x 00391 x.rGetLocation() += phi(i)*r_node_loc; 00392 } 00393 00394 // Interpolate u and grad u if a current solution or guess exists 00395 unsigned node_global_index = rElement.GetNodeGlobalIndex(i); 00396 if (this->mCurrentSolutionOrGuessReplicated.GetSize() > 0) 00397 { 00398 for (unsigned index_of_unknown=0; index_of_unknown<(INTERPOLATION_LEVEL!=CARDIAC ? PROBLEM_DIM : 1); index_of_unknown++) 00399 { 00400 /* 00401 * If we have a solution (e.g. this is a dynamic problem) then 00402 * interpolate the value at this quadrature point. 00403 * 00404 * NOTE: the following assumes that if, say, there are two unknowns 00405 * u and v, they are stored in the current solution vector as 00406 * [U1 V1 U2 V2 ... U_n V_n]. 00407 */ 00408 double u_at_node = this->GetCurrentSolutionOrGuessValue(node_global_index, index_of_unknown); 00409 u(index_of_unknown) += phi(i)*u_at_node; 00410 00412 // if (INTERPOLATION_LEVEL==NONLINEAR) // don't need to construct grad_phi or grad_u in other cases 00413 // { 00414 // for (unsigned j=0; j<SPACE_DIM; j++) 00415 // { 00416 // grad_u(index_of_unknown,j) += grad_phi(j,i)*u_at_node; 00417 // } 00418 // } 00419 } 00420 } 00421 00422 // Allow the concrete version of the assembler to interpolate any desired quantities 00423 IncrementInterpolatedQuantities(phi(i), p_node); 00424 } 00425 00426 double wJ = jacobian_determinant * mpCableQuadRule->GetWeight(quad_index); 00427 00428 // Create rAElem and rBElem 00429 if (this->mAssembleMatrix) 00430 { 00431 noalias(rAElem) += ComputeCableMatrixTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ; 00432 } 00433 00434 if (this->mAssembleVector) 00435 { 00436 noalias(rBElem) += ComputeCableVectorTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ; 00437 } 00438 } 00439 } 00440 00441 #endif /*ABSTRACTFECABLEINTEGRALASSEMBLER_HPP_*/