Chaste  Release::3.4
AbstractFeCableIntegralAssembler.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 ABSTRACTFECABLEINTEGRALASSEMBLER_HPP_
37 #define ABSTRACTFECABLEINTEGRALASSEMBLER_HPP_
38 
39 #include "AbstractFeAssemblerCommon.hpp"
40 #include "GaussianQuadratureRule.hpp"
41 #include "PetscVecTools.hpp"
42 #include "PetscMatTools.hpp"
43 
50 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
51 class AbstractFeCableIntegralAssembler : public AbstractFeAssemblerCommon<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,CAN_ASSEMBLE_VECTOR,CAN_ASSEMBLE_MATRIX,INTERPOLATION_LEVEL>
52 {
53 protected:
55  static const unsigned CABLE_ELEMENT_DIM = 1;
56 
58  static const unsigned NUM_CABLE_ELEMENT_NODES = 2;
59 
62 
65 
68 
84  const c_matrix<double, CABLE_ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
85  c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rReturnValue);
86 
92  void DoAssemble();
93 
116  virtual c_matrix<double,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES> ComputeCableMatrixTerm(
117  c_vector<double, NUM_CABLE_ELEMENT_NODES>& rPhi,
118  c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rGradPhi,
120  c_vector<double,PROBLEM_DIM>& rU,
121  c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
123  {
124  // If this line is reached this means this method probably hasn't been over-ridden correctly in
125  // the concrete class
127  return zero_matrix<double>(PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES);
128  }
129 
150  virtual c_vector<double,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES> ComputeCableVectorTerm(
151  c_vector<double, NUM_CABLE_ELEMENT_NODES>& rPhi,
152  c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rGradPhi,
154  c_vector<double,PROBLEM_DIM>& rU,
155  c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
157  {
158  // If this line is reached this means this method probably hasn't been over-ridden correctly in
159  // the concrete class
161  return zero_vector<double>(PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES);
162  }
163 
179  c_matrix<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES >& rAElem,
180  c_vector<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES>& rBElem);
181 
182 
191  {
192  return true;
193  }
194 
195 public:
196 
203 
208  {
209  delete mpCableQuadRule;
210  }
211 };
212 
213 
214 
215 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
218  : AbstractFeAssemblerCommon<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>(),
219  mpMesh(pMesh)
220 {
221  assert(pMesh);
222  assert(CAN_ASSEMBLE_VECTOR || CAN_ASSEMBLE_MATRIX);
223  // Default to 2nd order quadrature. Our default basis functions are piecewise linear
224  // which means that we are integrating functions which in the worst case (mass matrix)
225  // are quadratic.
227 
228  // Not supporting this yet - if a nonlinear assembler on cable elements is required, uncomment code
229  // in AssembleOnCableElement below (search for NONLINEAR)
230  assert(INTERPOLATION_LEVEL!=NONLINEAR);
231 }
232 
233 
234 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
236 {
237  HeartEventHandler::EventType assemble_event;
238  if (this->mAssembleMatrix)
239  {
240  assemble_event = HeartEventHandler::ASSEMBLE_SYSTEM;
241  }
242  else
243  {
244  assemble_event = HeartEventHandler::ASSEMBLE_RHS;
245  }
246 
247  if (this->mAssembleMatrix && this->mMatrixToAssemble==NULL)
248  {
249  EXCEPTION("Matrix to be assembled has not been set");
250  }
251  if (this->mAssembleVector && this->mVectorToAssemble==NULL)
252  {
253  EXCEPTION("Vector to be assembled has not been set");
254  }
255 
256  // This has to be below PrepareForAssembleSystem as in that method the ODEs are solved in cardiac problems
257  HeartEventHandler::BeginEvent(assemble_event);
258 
259  // Zero the matrix/vector if it is to be assembled
260  if (this->mAssembleVector && this->mZeroVectorBeforeAssembly)
261  {
262  PetscVecTools::Zero(this->mVectorToAssemble);
263  }
264  if (this->mAssembleMatrix && this->mZeroMatrixBeforeAssembly)
265  {
266  PetscMatTools::Zero(this->mMatrixToAssemble);
267  }
268 
269  const size_t STENCIL_SIZE=PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES;
270  c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
271  c_vector<double, STENCIL_SIZE> b_elem;
272 
273  // Loop over elements
274  if (this->mAssembleMatrix || this->mAssembleVector)
275  {
276  for (typename MixedDimensionMesh<CABLE_ELEMENT_DIM, SPACE_DIM>::CableElementIterator iter = mpMesh->GetCableElementIteratorBegin();
277  iter != mpMesh->GetCableElementIteratorEnd();
278  ++iter)
279  {
280  Element<CABLE_ELEMENT_DIM, SPACE_DIM>& r_element = *(*iter);
281 
282  // Test for ownership first, since it's pointless to test the criterion on something which we might know nothing about.
283  if ( r_element.GetOwnership() == true && ElementAssemblyCriterion(r_element)==true )
284  {
285  AssembleOnCableElement(r_element, a_elem, b_elem);
286 
287  unsigned p_indices[STENCIL_SIZE];
288  r_element.GetStiffnessMatrixGlobalIndices(PROBLEM_DIM, p_indices);
289 
290  if (this->mAssembleMatrix)
291  {
292  PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->mMatrixToAssemble, p_indices, a_elem);
293  }
294 
295  if (this->mAssembleVector)
296  {
297  PetscVecTools::AddMultipleValues<STENCIL_SIZE>(this->mVectorToAssemble, p_indices, b_elem);
298  }
299  }
300  }
301  }
302 
303  HeartEventHandler::EndEvent(assemble_event);
304 }
305 
307 // Implementation - AssembleOnCableElement and smaller
309 
310 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
312  const ChastePoint<CABLE_ELEMENT_DIM>& rPoint,
313  const c_matrix<double, CABLE_ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
314  c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rReturnValue)
315 {
316  static c_matrix<double, CABLE_ELEMENT_DIM, NUM_CABLE_ELEMENT_NODES> grad_phi;
317 
319  rReturnValue = prod(trans(rInverseJacobian), grad_phi);
320 }
321 
322 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
325  c_matrix<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES >& rAElem,
326  c_vector<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES>& rBElem)
327 {
333  c_matrix<double, SPACE_DIM, CABLE_ELEMENT_DIM> jacobian;
334  c_matrix<double, CABLE_ELEMENT_DIM, SPACE_DIM> inverse_jacobian;
335  double jacobian_determinant;
336 
339  // mpMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
341  rElement.CalculateInverseJacobian(jacobian, jacobian_determinant, inverse_jacobian);
342 
343  if (this->mAssembleMatrix)
344  {
345  rAElem.clear();
346  }
347 
348  if (this->mAssembleVector)
349  {
350  rBElem.clear();
351  }
352 
353  const unsigned num_nodes = rElement.GetNumNodes();
354 
355  // Allocate memory for the basis functions values and derivative values
356  c_vector<double, NUM_CABLE_ELEMENT_NODES> phi;
357  c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES> grad_phi;
358 
359  // Loop over Gauss points
360  for (unsigned quad_index=0; quad_index < mpCableQuadRule->GetNumQuadPoints(); quad_index++)
361  {
362  const ChastePoint<CABLE_ELEMENT_DIM>& quad_point = mpCableQuadRule->rGetQuadPoint(quad_index);
363 
364  CableBasisFunction::ComputeBasisFunctions(quad_point, phi);
365 
366  if ( this->mAssembleMatrix || INTERPOLATION_LEVEL==NONLINEAR )
367  {
368  ComputeTransformedBasisFunctionDerivatives(quad_point, inverse_jacobian, grad_phi);
369  }
370 
371  // Location of the Gauss point in the original element will be stored in x
372  // Where applicable, u will be set to the value of the current solution at x
373  ChastePoint<SPACE_DIM> x(0,0,0);
374 
375  c_vector<double,PROBLEM_DIM> u = zero_vector<double>(PROBLEM_DIM);
376  c_matrix<double,PROBLEM_DIM,SPACE_DIM> grad_u = zero_matrix<double>(PROBLEM_DIM,SPACE_DIM);
377 
378  // Allow the concrete version of the assembler to interpolate any desired quantities
379  this->ResetInterpolatedQuantities();
380 
381  // Interpolation
382  for (unsigned i=0; i<num_nodes; i++)
383  {
384  const Node<SPACE_DIM>* p_node = rElement.GetNode(i);
385 
386  if (INTERPOLATION_LEVEL != CARDIAC) // don't interpolate X if cardiac problem
387  {
388  const c_vector<double, SPACE_DIM>& r_node_loc = p_node->rGetLocation();
389  // interpolate x
390  x.rGetLocation() += phi(i)*r_node_loc;
391  }
392 
393  // Interpolate u and grad u if a current solution or guess exists
394  unsigned node_global_index = rElement.GetNodeGlobalIndex(i);
395  if (this->mCurrentSolutionOrGuessReplicated.GetSize() > 0)
396  {
397  for (unsigned index_of_unknown=0; index_of_unknown<(INTERPOLATION_LEVEL!=CARDIAC ? PROBLEM_DIM : 1); index_of_unknown++)
398  {
399  /*
400  * If we have a solution (e.g. this is a dynamic problem) then
401  * interpolate the value at this quadrature point.
402  *
403  * NOTE: the following assumes that if, say, there are two unknowns
404  * u and v, they are stored in the current solution vector as
405  * [U1 V1 U2 V2 ... U_n V_n].
406  */
407  double u_at_node = this->GetCurrentSolutionOrGuessValue(node_global_index, index_of_unknown);
408  u(index_of_unknown) += phi(i)*u_at_node;
409 
411 // if (INTERPOLATION_LEVEL==NONLINEAR) // don't need to construct grad_phi or grad_u in other cases
412 // {
413 // for (unsigned j=0; j<SPACE_DIM; j++)
414 // {
415 // grad_u(index_of_unknown,j) += grad_phi(j,i)*u_at_node;
416 // }
417 // }
418  }
419  }
420 
421  // Allow the concrete version of the assembler to interpolate any desired quantities
422  this->IncrementInterpolatedQuantities(phi(i), p_node);
423  }
424 
425  double wJ = jacobian_determinant * mpCableQuadRule->GetWeight(quad_index);
426 
427  // Create rAElem and rBElem
428  if (this->mAssembleMatrix)
429  {
430  noalias(rAElem) += ComputeCableMatrixTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
431  }
432 
433  if (this->mAssembleVector)
434  {
435  noalias(rBElem) += ComputeCableVectorTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
436  }
437  }
438 }
439 
440 #endif /*ABSTRACTFECABLEINTEGRALASSEMBLER_HPP_*/
Definition: Node.hpp:58
c_vector< double, DIM > & rGetLocation()
Definition: ChastePoint.cpp:76
virtual c_matrix< double, PROBLEM_DIM *NUM_CABLE_ELEMENT_NODES, PROBLEM_DIM *NUM_CABLE_ELEMENT_NODES > ComputeCableMatrixTerm(c_vector< double, NUM_CABLE_ELEMENT_NODES > &rPhi, c_matrix< double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< CABLE_ELEMENT_DIM, SPACE_DIM > *pElement)
unsigned GetNodeGlobalIndex(unsigned localIndex) const
#define EXCEPTION(message)
Definition: Exception.hpp:143
static void ComputeBasisFunctionDerivatives(const ChastePoint< ELEMENT_DIM > &rPoint, c_matrix< double, ELEMENT_DIM, ELEMENT_DIM+1 > &rReturnValue)
#define NEVER_REACHED
Definition: Exception.hpp:206
void CalculateInverseJacobian(c_matrix< double, SPACE_DIM, ELEMENT_DIM > &rJacobian, double &rJacobianDeterminant, c_matrix< double, ELEMENT_DIM, SPACE_DIM > &rInverseJacobian)
Node< SPACE_DIM > * GetNode(unsigned localIndex) const
bool GetOwnership() const
static void Zero(Mat matrix)
virtual bool ElementAssemblyCriterion(Element< CABLE_ELEMENT_DIM, SPACE_DIM > &rElement)
unsigned GetNumNodes() const
void GetStiffnessMatrixGlobalIndices(unsigned problemDim, unsigned *pIndices) const
MixedDimensionMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
static void Zero(Vec vector)
const c_vector< double, SPACE_DIM > & rGetLocation() const
Definition: Node.cpp:140
virtual c_vector< double, PROBLEM_DIM *NUM_CABLE_ELEMENT_NODES > ComputeCableVectorTerm(c_vector< double, NUM_CABLE_ELEMENT_NODES > &rPhi, c_matrix< double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< CABLE_ELEMENT_DIM, SPACE_DIM > *pElement)
void ComputeTransformedBasisFunctionDerivatives(const ChastePoint< CABLE_ELEMENT_DIM > &rPoint, const c_matrix< double, CABLE_ELEMENT_DIM, SPACE_DIM > &rInverseJacobian, c_matrix< double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES > &rReturnValue)
AbstractFeCableIntegralAssembler(MixedDimensionMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
virtual void AssembleOnCableElement(Element< CABLE_ELEMENT_DIM, SPACE_DIM > &rElement, c_matrix< double, PROBLEM_DIM *NUM_CABLE_ELEMENT_NODES, PROBLEM_DIM *NUM_CABLE_ELEMENT_NODES > &rAElem, c_vector< double, PROBLEM_DIM *NUM_CABLE_ELEMENT_NODES > &rBElem)
std::vector< Element< 1, SPACE_DIM > * >::const_iterator CableElementIterator