Chaste  Release::3.4
AbstractFeVolumeIntegralAssembler.hpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34 */
35 
36 #ifndef ABSTRACTFEVOLUMEINTEGRALASSEMBLER_HPP_
37 #define ABSTRACTFEVOLUMEINTEGRALASSEMBLER_HPP_
38 
39 #include "AbstractFeAssemblerCommon.hpp"
40 #include "GaussianQuadratureRule.hpp"
41 #include "BoundaryConditionsContainer.hpp"
42 #include "PetscVecTools.hpp"
43 #include "PetscMatTools.hpp"
44 
77 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
79  public AbstractFeAssemblerCommon<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,CAN_ASSEMBLE_VECTOR,CAN_ASSEMBLE_MATRIX,INTERPOLATION_LEVEL>
80 {
81 protected:
84 
87 
90 
111  const c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
112  c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rReturnValue);
113 
119  void DoAssemble();
120 
121 protected:
122 
142  virtual c_matrix<double,PROBLEM_DIM*(ELEMENT_DIM+1),PROBLEM_DIM*(ELEMENT_DIM+1)> ComputeMatrixTerm(
143  c_vector<double, ELEMENT_DIM+1>& rPhi,
144  c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rGradPhi,
146  c_vector<double,PROBLEM_DIM>& rU,
147  c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
149  {
150  // If this line is reached this means this method probably hasn't been over-ridden correctly in
151  // the concrete class
153  return zero_matrix<double>(PROBLEM_DIM*(ELEMENT_DIM+1),PROBLEM_DIM*(ELEMENT_DIM+1));
154  }
155 
175  virtual c_vector<double,PROBLEM_DIM*(ELEMENT_DIM+1)> ComputeVectorTerm(
176  c_vector<double, ELEMENT_DIM+1>& rPhi,
177  c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rGradPhi,
179  c_vector<double,PROBLEM_DIM>& rU,
180  c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
182  {
183  // If this line is reached this means this method probably hasn't been over-ridden correctly in
184  // the concrete class
186  return zero_vector<double>(PROBLEM_DIM*(ELEMENT_DIM+1));
187  }
188 
189 
204  virtual void AssembleOnElement(Element<ELEMENT_DIM,SPACE_DIM>& rElement,
205  c_matrix<double, PROBLEM_DIM*(ELEMENT_DIM+1), PROBLEM_DIM*(ELEMENT_DIM+1) >& rAElem,
206  c_vector<double, PROBLEM_DIM*(ELEMENT_DIM+1)>& rBElem);
207 
216  {
217  return true;
218  }
219 
220 
221 public:
222 
229 
234  {
235  delete mpQuadRule;
236  }
237 };
238 
239 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
242  : AbstractFeAssemblerCommon<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>(),
243  mpMesh(pMesh)
244 {
245  assert(pMesh);
246  // Default to 2nd order quadrature. Our default basis functions are piecewise linear
247  // which means that we are integrating functions which in the worst case (mass matrix)
248  // are quadratic.
250 }
251 
252 
253 
254 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
256 {
257  assert(this->mAssembleMatrix || this->mAssembleVector);
258 
259  HeartEventHandler::EventType assemble_event;
260  if (this->mAssembleMatrix)
261  {
262  assemble_event = HeartEventHandler::ASSEMBLE_SYSTEM;
263  }
264  else
265  {
266  assemble_event = HeartEventHandler::ASSEMBLE_RHS;
267  }
268 
269  if (this->mAssembleMatrix && this->mMatrixToAssemble==NULL)
270  {
271  EXCEPTION("Matrix to be assembled has not been set");
272  }
273  if (this->mAssembleVector && this->mVectorToAssemble==NULL)
274  {
275  EXCEPTION("Vector to be assembled has not been set");
276  }
277 
278  HeartEventHandler::BeginEvent(assemble_event);
279 
280  // Zero the matrix/vector if it is to be assembled
281  if (this->mAssembleVector && this->mZeroVectorBeforeAssembly)
282  {
283  PetscVecTools::Zero(this->mVectorToAssemble);
284  }
285  if (this->mAssembleMatrix && this->mZeroMatrixBeforeAssembly)
286  {
287  PetscMatTools::Zero(this->mMatrixToAssemble);
288  }
289 
290  const size_t STENCIL_SIZE=PROBLEM_DIM*(ELEMENT_DIM+1);
291  c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
292  c_vector<double, STENCIL_SIZE> b_elem;
293 
294  // Loop over elements
295  for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = mpMesh->GetElementIteratorBegin();
296  iter != mpMesh->GetElementIteratorEnd();
297  ++iter)
298  {
299  Element<ELEMENT_DIM, SPACE_DIM>& r_element = *iter;
300 
301  // Test for ownership first, since it's pointless to test the criterion on something which we might know nothing about.
302  if ( r_element.GetOwnership() == true && ElementAssemblyCriterion(r_element)==true )
303  {
304  AssembleOnElement(r_element, a_elem, b_elem);
305 
306  unsigned p_indices[STENCIL_SIZE];
307  r_element.GetStiffnessMatrixGlobalIndices(PROBLEM_DIM, p_indices);
308 
309  if (this->mAssembleMatrix)
310  {
311  PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->mMatrixToAssemble, p_indices, a_elem);
312  }
313 
314  if (this->mAssembleVector)
315  {
316  PetscVecTools::AddMultipleValues<STENCIL_SIZE>(this->mVectorToAssemble, p_indices, b_elem);
317  }
318  }
319  }
320 
321  HeartEventHandler::EndEvent(assemble_event);
322 }
323 
324 
326 // Implementation - AssembleOnElement and smaller
328 
329 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
331  const ChastePoint<ELEMENT_DIM>& rPoint,
332  const c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
333  c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rReturnValue)
334 {
335  assert(ELEMENT_DIM < 4 && ELEMENT_DIM > 0);
336  static c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> grad_phi;
337 
339  rReturnValue = prod(trans(rInverseJacobian), grad_phi);
340 }
341 
342 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
345  c_matrix<double, PROBLEM_DIM*(ELEMENT_DIM+1), PROBLEM_DIM*(ELEMENT_DIM+1) >& rAElem,
346  c_vector<double, PROBLEM_DIM*(ELEMENT_DIM+1)>& rBElem)
347 {
353  c_matrix<double, SPACE_DIM, ELEMENT_DIM> jacobian;
354  c_matrix<double, ELEMENT_DIM, SPACE_DIM> inverse_jacobian;
355  double jacobian_determinant;
356 
357  mpMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
358 
359  if (this->mAssembleMatrix)
360  {
361  rAElem.clear();
362  }
363 
364  if (this->mAssembleVector)
365  {
366  rBElem.clear();
367  }
368 
369  const unsigned num_nodes = rElement.GetNumNodes();
370 
371  // Allocate memory for the basis functions values and derivative values
372  c_vector<double, ELEMENT_DIM+1> phi;
373  c_matrix<double, SPACE_DIM, ELEMENT_DIM+1> grad_phi;
374 
375  // Loop over Gauss points
376  for (unsigned quad_index=0; quad_index < mpQuadRule->GetNumQuadPoints(); quad_index++)
377  {
378  const ChastePoint<ELEMENT_DIM>& quad_point = mpQuadRule->rGetQuadPoint(quad_index);
379 
380  BasisFunction::ComputeBasisFunctions(quad_point, phi);
381 
382  if ( this->mAssembleMatrix || INTERPOLATION_LEVEL==NONLINEAR )
383  {
384  ComputeTransformedBasisFunctionDerivatives(quad_point, inverse_jacobian, grad_phi);
385  }
386 
387  // Location of the Gauss point in the original element will be stored in x
388  // Where applicable, u will be set to the value of the current solution at x
389  ChastePoint<SPACE_DIM> x(0,0,0);
390 
391  c_vector<double,PROBLEM_DIM> u = zero_vector<double>(PROBLEM_DIM);
392  c_matrix<double,PROBLEM_DIM,SPACE_DIM> grad_u = zero_matrix<double>(PROBLEM_DIM,SPACE_DIM);
393 
394  // Allow the concrete version of the assembler to interpolate any desired quantities
395  this->ResetInterpolatedQuantities();
396 
397  // Interpolation
398  for (unsigned i=0; i<num_nodes; i++)
399  {
400  const Node<SPACE_DIM>* p_node = rElement.GetNode(i);
401 
402  if (INTERPOLATION_LEVEL != CARDIAC) // don't even interpolate X if cardiac problem
403  {
404  const c_vector<double, SPACE_DIM>& r_node_loc = p_node->rGetLocation();
405  // interpolate x
406  x.rGetLocation() += phi(i)*r_node_loc;
407  }
408 
409  // Interpolate u and grad u if a current solution or guess exists
410  unsigned node_global_index = rElement.GetNodeGlobalIndex(i);
411  if (this->mCurrentSolutionOrGuessReplicated.GetSize() > 0)
412  {
413  for (unsigned index_of_unknown=0; index_of_unknown<(INTERPOLATION_LEVEL!=CARDIAC ? PROBLEM_DIM : 1); index_of_unknown++)
414  {
415  /*
416  * If we have a solution (e.g. this is a dynamic problem) then
417  * interpolate the value at this quadrature point.
418  *
419  * NOTE: the following assumes that if, say, there are two unknowns
420  * u and v, they are stored in the current solution vector as
421  * [U1 V1 U2 V2 ... U_n V_n].
422  */
423  double u_at_node = this->GetCurrentSolutionOrGuessValue(node_global_index, index_of_unknown);
424  u(index_of_unknown) += phi(i)*u_at_node;
425 
426  if (INTERPOLATION_LEVEL==NONLINEAR) // don't need to construct grad_phi or grad_u in other cases
427  {
428  for (unsigned j=0; j<SPACE_DIM; j++)
429  {
430  grad_u(index_of_unknown,j) += grad_phi(j,i)*u_at_node;
431  }
432  }
433  }
434  }
435 
436  // Allow the concrete version of the assembler to interpolate any desired quantities
437  this->IncrementInterpolatedQuantities(phi(i), p_node);
438  if ( this->mAssembleMatrix || INTERPOLATION_LEVEL==NONLINEAR )
439  {
440  this->IncrementInterpolatedGradientQuantities(grad_phi, i, p_node);
441  }
442  }
443 
444  double wJ = jacobian_determinant * mpQuadRule->GetWeight(quad_index);
445 
446  // Create rAElem and rBElem
447  if (this->mAssembleMatrix)
448  {
449  noalias(rAElem) += ComputeMatrixTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
450  }
451 
452  if (this->mAssembleVector)
453  {
454  noalias(rBElem) += ComputeVectorTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
455  }
456  }
457 }
458 
459 
460 #endif /*ABSTRACTFEVOLUMEINTEGRALASSEMBLER_HPP_*/
Definition: Node.hpp:58
c_vector< double, DIM > & rGetLocation()
Definition: ChastePoint.cpp:76
unsigned GetNodeGlobalIndex(unsigned localIndex) const
#define EXCEPTION(message)
Definition: Exception.hpp:143
LinearBasisFunction< ELEMENT_DIM > BasisFunction
static void ComputeBasisFunctionDerivatives(const ChastePoint< ELEMENT_DIM > &rPoint, c_matrix< double, ELEMENT_DIM, ELEMENT_DIM+1 > &rReturnValue)
#define NEVER_REACHED
Definition: Exception.hpp:206
virtual c_vector< double, PROBLEM_DIM *(ELEMENT_DIM+1)> ComputeVectorTerm(c_vector< double, ELEMENT_DIM+1 > &rPhi, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement)
Node< SPACE_DIM > * GetNode(unsigned localIndex) const
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
bool GetOwnership() const
static void Zero(Mat matrix)
GaussianQuadratureRule< ELEMENT_DIM > * mpQuadRule
unsigned GetNumNodes() const
void GetStiffnessMatrixGlobalIndices(unsigned problemDim, unsigned *pIndices) const
virtual void AssembleOnElement(Element< ELEMENT_DIM, SPACE_DIM > &rElement, c_matrix< double, PROBLEM_DIM *(ELEMENT_DIM+1), PROBLEM_DIM *(ELEMENT_DIM+1) > &rAElem, c_vector< double, PROBLEM_DIM *(ELEMENT_DIM+1)> &rBElem)
virtual bool ElementAssemblyCriterion(Element< ELEMENT_DIM, SPACE_DIM > &rElement)
static void Zero(Vec vector)
void ComputeTransformedBasisFunctionDerivatives(const ChastePoint< ELEMENT_DIM > &rPoint, const c_matrix< double, ELEMENT_DIM, SPACE_DIM > &rInverseJacobian, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rReturnValue)
const c_vector< double, SPACE_DIM > & rGetLocation() const
Definition: Node.cpp:140
virtual c_matrix< double, PROBLEM_DIM *(ELEMENT_DIM+1), PROBLEM_DIM *(ELEMENT_DIM+1)> ComputeMatrixTerm(c_vector< double, ELEMENT_DIM+1 > &rPhi, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement)
unsigned GetIndex() const
AbstractFeVolumeIntegralAssembler(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)