Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
AbstractContinuumMechanicsAssembler.hpp
1/*
2
3Copyright (c) 2005-2024, University of Oxford.
4All rights reserved.
5
6University of Oxford means the Chancellor, Masters and Scholars of the
7University of Oxford, having an administrative office at Wellington
8Square, Oxford OX1 2JD, UK.
9
10This file is part of Chaste.
11
12Redistribution and use in source and binary forms, with or without
13modification, 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
23THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33
34*/
35
36#ifndef ABSTRACTCONTINUUMMECHANICSASSEMBLER_HPP_
37#define ABSTRACTCONTINUUMMECHANICSASSEMBLER_HPP_
38
39#include "AbstractFeAssemblerInterface.hpp"
40#include "AbstractTetrahedralMesh.hpp"
41#include "QuadraticMesh.hpp"
42#include "DistributedQuadraticMesh.hpp"
43#include "LinearBasisFunction.hpp"
44#include "QuadraticBasisFunction.hpp"
45#include "ReplicatableVector.hpp"
46#include "DistributedVector.hpp"
47#include "PetscTools.hpp"
48#include "PetscVecTools.hpp"
49#include "PetscMatTools.hpp"
50#include "GaussianQuadratureRule.hpp"
51
52
86template<unsigned DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX>
87class AbstractContinuumMechanicsAssembler : public AbstractFeAssemblerInterface<CAN_ASSEMBLE_VECTOR,CAN_ASSEMBLE_MATRIX>
88{
89protected:
93 static const bool BLOCK_SYMMETRIC_MATRIX = true; //generalise to non-block symmetric matrices later (when needed maybe)
94
96 static const unsigned NUM_VERTICES_PER_ELEMENT = DIM+1;
97
99 static const unsigned NUM_NODES_PER_ELEMENT = (DIM+1)*(DIM+2)/2; // assuming quadratic
100
105
108
111
114
122
123
143 virtual c_matrix<double,SPATIAL_BLOCK_SIZE_ELEMENTAL,SPATIAL_BLOCK_SIZE_ELEMENTAL> ComputeSpatialSpatialMatrixTerm(
144 c_vector<double, NUM_NODES_PER_ELEMENT>& rQuadPhi,
145 c_matrix<double, DIM, NUM_NODES_PER_ELEMENT>& rGradQuadPhi,
146 c_vector<double,DIM>& rX,
147 Element<DIM,DIM>* pElement)
148 {
150 }
151
174 virtual c_matrix<double,SPATIAL_BLOCK_SIZE_ELEMENTAL,PRESSURE_BLOCK_SIZE_ELEMENTAL> ComputeSpatialPressureMatrixTerm(
175 c_vector<double, NUM_NODES_PER_ELEMENT>& rQuadPhi,
176 c_matrix<double, DIM, NUM_NODES_PER_ELEMENT>& rGradQuadPhi,
177 c_vector<double, NUM_VERTICES_PER_ELEMENT>& rLinearPhi,
178 c_matrix<double, DIM, NUM_VERTICES_PER_ELEMENT>& rGradLinearPhi,
179 c_vector<double,DIM>& rX,
180 Element<DIM,DIM>* pElement)
181 {
183 }
184
185
205 virtual c_matrix<double,PRESSURE_BLOCK_SIZE_ELEMENTAL,PRESSURE_BLOCK_SIZE_ELEMENTAL> ComputePressurePressureMatrixTerm(
206 c_vector<double, NUM_VERTICES_PER_ELEMENT>& rLinearPhi,
207 c_matrix<double, DIM, NUM_VERTICES_PER_ELEMENT>& rGradLinearPhi,
208 c_vector<double,DIM>& rX,
209 Element<DIM,DIM>* pElement)
210 {
212 }
213
214
237 virtual c_vector<double,SPATIAL_BLOCK_SIZE_ELEMENTAL> ComputeSpatialVectorTerm(
238 c_vector<double, NUM_NODES_PER_ELEMENT>& rQuadPhi,
239 c_matrix<double, DIM, NUM_NODES_PER_ELEMENT>& rGradQuadPhi,
240 c_vector<double,DIM>& rX,
241 Element<DIM,DIM>* pElement) = 0;
242
243
266 virtual c_vector<double,PRESSURE_BLOCK_SIZE_ELEMENTAL> ComputePressureVectorTerm(
267 c_vector<double, NUM_VERTICES_PER_ELEMENT>& rLinearPhi,
268 c_matrix<double, DIM, NUM_VERTICES_PER_ELEMENT>& rGradLinearPhi,
269 c_vector<double,DIM>& rX,
270 Element<DIM,DIM>* pElement)
271 {
272 return zero_vector<double>(PRESSURE_BLOCK_SIZE_ELEMENTAL);
273 }
274
287 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElem,
288 c_vector<double, STENCIL_SIZE>& rBElem);
289
290public:
295 : AbstractFeAssemblerInterface<CAN_ASSEMBLE_VECTOR,CAN_ASSEMBLE_MATRIX>(),
296 mpMesh(pMesh)
297 {
298 assert(pMesh);
299
300 // Check that the mesh is quadratic
301 QuadraticMesh<DIM>* p_quad_mesh = dynamic_cast<QuadraticMesh<DIM>* >(pMesh);
302 DistributedQuadraticMesh<DIM>* p_distributed_quad_mesh = dynamic_cast<DistributedQuadraticMesh<DIM>* >(pMesh);
303
304 if ((p_quad_mesh == NULL) && (p_distributed_quad_mesh == NULL))
305 {
306 EXCEPTION("Continuum mechanics assemblers require a quadratic mesh");
307 }
308
309 // In general the Jacobian for a mechanics problem is non-polynomial.
310 // We therefore use the highest order integration rule available
312 }
313
314// void SetCurrentSolution(Vec currentSolution);
315
320 {
321 delete mpQuadRule;
322 }
323};
324
325
327//template<unsigned DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX>
328//void AbstractContinuumMechanicsAssembler<DIM,CAN_ASSEMBLE_VECTOR,CAN_ASSEMBLE_MATRIX>::SetCurrentSolution(Vec currentSolution)
329//{
330// assert(currentSolution != NULL);
331//
332// // Replicate the current solution and store so can be used in AssembleOnElement
333// HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
334// mCurrentSolutionOrGuessReplicated.ReplicatePetscVector(currentSolution);
335// HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
336//
337// // The AssembleOnElement type methods will determine if a current solution or
338// // current guess exists by looking at the size of the replicated vector, so
339// // check the size is zero if there isn't a current solution.
340// assert(mCurrentSolutionOrGuessReplicated.GetSize() > 0);
341//}
342
343template<unsigned DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX>
345{
346 assert(this->mAssembleMatrix || this->mAssembleVector);
347 if (this->mAssembleMatrix)
348 {
349 if (this->mMatrixToAssemble == NULL)
350 {
351 EXCEPTION("Matrix to be assembled has not been set");
352 }
353 if (PetscMatTools::GetSize(this->mMatrixToAssemble) != (DIM+1)*mpMesh->GetNumNodes())
354 {
355 EXCEPTION("Matrix provided to be assembled has size " << PetscMatTools::GetSize(this->mMatrixToAssemble) << ", not expected size of " << (DIM+1)*mpMesh->GetNumNodes() << " ((dim+1)*num_nodes)");
356 }
357 }
358
359 if (this->mAssembleVector)
360 {
361 if (this->mVectorToAssemble == NULL)
362 {
363 EXCEPTION("Vector to be assembled has not been set");
364 }
365 if (PetscVecTools::GetSize(this->mVectorToAssemble) != (DIM+1)*mpMesh->GetNumNodes())
366 {
367 EXCEPTION("Vector provided to be assembled has size " << PetscVecTools::GetSize(this->mVectorToAssemble) << ", not expected size of " << (DIM+1)*mpMesh->GetNumNodes() << " ((dim+1)*num_nodes)");
368 }
369 }
370
371 // Zero the matrix/vector if it is to be assembled
372 if (this->mAssembleVector && this->mZeroVectorBeforeAssembly)
373 {
374 // Note PetscVecTools::Finalise(this->mVectorToAssemble); on an unused matrix
375 // would "compress" data and make any pre-allocated entries redundant.
376 PetscVecTools::Zero(this->mVectorToAssemble);
377 }
378 if (this->mAssembleMatrix && this->mZeroMatrixBeforeAssembly)
379 {
380 // Note PetscMatTools::Finalise(this->mMatrixToAssemble); on an unused matrix
381 // would "compress" data and make any pre-allocated entries redundant.
382 PetscMatTools::Zero(this->mMatrixToAssemble);
383 }
384
385 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem = zero_matrix<double>(STENCIL_SIZE,STENCIL_SIZE);
386 c_vector<double, STENCIL_SIZE> b_elem = zero_vector<double>(STENCIL_SIZE);
387
388
389 // Loop over elements
390 for (typename AbstractTetrahedralMesh<DIM, DIM>::ElementIterator iter = mpMesh->GetElementIteratorBegin();
391 iter != mpMesh->GetElementIteratorEnd();
392 ++iter)
393 {
394 Element<DIM, DIM>& r_element = *iter;
395
396 // Test for ownership first, since it's pointless to test the criterion on something which we might know nothing about.
397 if (r_element.GetOwnership() == true /*&& ElementAssemblyCriterion(r_element)==true*/)
398 {
399 // LCOV_EXCL_START
400 // note: if assemble matrix only
401 if (CommandLineArguments::Instance()->OptionExists("-mech_very_verbose") && this->mAssembleMatrix)
402 {
403 std::cout << "\r[" << PetscTools::GetMyRank() << "]: Element " << r_element.GetIndex() << " of " << mpMesh->GetNumElements() << std::flush;
404 }
405 // LCOV_EXCL_STOP
406
407 AssembleOnElement(r_element, a_elem, b_elem);
408
409 // Note that a different ordering is used for the elemental matrix compared to the global matrix.
410 // See comments about ordering above.
411 unsigned p_indices[STENCIL_SIZE];
412 // Work out the mapping for spatial terms
413 for (unsigned i=0; i<NUM_NODES_PER_ELEMENT; i++)
414 {
415 for (unsigned j=0; j<DIM; j++)
416 {
417 // DIM+1 on the right-hand side here is the problem dimension
418 p_indices[DIM*i+j] = (DIM+1)*r_element.GetNodeGlobalIndex(i) + j;
419 }
420 }
421 // Work out the mapping for pressure terms
422 for (unsigned i=0; i<NUM_VERTICES_PER_ELEMENT; i++)
423 {
424 p_indices[DIM*NUM_NODES_PER_ELEMENT + i] = (DIM+1)*r_element.GetNodeGlobalIndex(i)+DIM;
425 }
426
427
428 if (this->mAssembleMatrix)
429 {
430 PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->mMatrixToAssemble, p_indices, a_elem);
431 }
432
433 if (this->mAssembleVector)
434 {
435 PetscVecTools::AddMultipleValues<STENCIL_SIZE>(this->mVectorToAssemble, p_indices, b_elem);
436 }
437 }
438 }
439}
440
441template<unsigned DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX>
443 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElem,
444 c_vector<double, STENCIL_SIZE>& rBElem)
445{
446 static c_matrix<double,DIM,DIM> jacobian;
447 static c_matrix<double,DIM,DIM> inverse_jacobian;
448 double jacobian_determinant;
449
450 mpMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
451
452 if (this->mAssembleMatrix)
453 {
454 rAElem.clear();
455 }
456
457 if (this->mAssembleVector)
458 {
459 rBElem.clear();
460 }
461
462 // Allocate memory for the basis functions values and derivative values
463 static c_vector<double, NUM_VERTICES_PER_ELEMENT> linear_phi;
464 static c_vector<double, NUM_NODES_PER_ELEMENT> quad_phi;
465 static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi;
466 static c_matrix<double, DIM, NUM_VERTICES_PER_ELEMENT> grad_linear_phi;
467
468 c_vector<double,DIM> body_force;
469
470 // Loop over Gauss points
471 for (unsigned quadrature_index=0; quadrature_index < mpQuadRule->GetNumQuadPoints(); quadrature_index++)
472 {
473 double wJ = jacobian_determinant * mpQuadRule->GetWeight(quadrature_index);
474 const ChastePoint<DIM>& quadrature_point = mpQuadRule->rGetQuadPoint(quadrature_index);
475
476 // Set up basis function info
477 LinearBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, linear_phi);
478 QuadraticBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, quad_phi);
479 QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_quad_phi);
480 LinearBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_linear_phi);
481
482 // interpolate X (ie physical location of this quad point).
483 c_vector<double,DIM> X = zero_vector<double>(DIM);
484 for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
485 {
486 for (unsigned j=0; j<DIM; j++)
487 {
488 X(j) += linear_phi(vertex_index)*rElement.GetNode(vertex_index)->rGetLocation()(j);
489 }
490 }
491
492 if (this->mAssembleVector)
493 {
494 c_vector<double,SPATIAL_BLOCK_SIZE_ELEMENTAL> b_spatial
495 = ComputeSpatialVectorTerm(quad_phi, grad_quad_phi, X, &rElement);
496 c_vector<double,PRESSURE_BLOCK_SIZE_ELEMENTAL> b_pressure = ComputePressureVectorTerm(linear_phi, grad_linear_phi, X, &rElement);
497
498 for (unsigned i=0; i<SPATIAL_BLOCK_SIZE_ELEMENTAL; i++)
499 {
500 rBElem(i) += b_spatial(i)*wJ;
501 }
502
503 for (unsigned i=0; i<PRESSURE_BLOCK_SIZE_ELEMENTAL; i++)
504 {
505 rBElem(SPATIAL_BLOCK_SIZE_ELEMENTAL + i) += b_pressure(i)*wJ;
506 }
507 }
508
509 if (this->mAssembleMatrix)
510 {
511 c_matrix<double,SPATIAL_BLOCK_SIZE_ELEMENTAL,SPATIAL_BLOCK_SIZE_ELEMENTAL> a_spatial_spatial
512 = ComputeSpatialSpatialMatrixTerm(quad_phi, grad_quad_phi, X, &rElement);
513
514 c_matrix<double,SPATIAL_BLOCK_SIZE_ELEMENTAL,PRESSURE_BLOCK_SIZE_ELEMENTAL> a_spatial_pressure
515 = ComputeSpatialPressureMatrixTerm(quad_phi, grad_quad_phi, linear_phi, grad_linear_phi, X, &rElement);
516
517 c_matrix<double,PRESSURE_BLOCK_SIZE_ELEMENTAL,SPATIAL_BLOCK_SIZE_ELEMENTAL> a_pressure_spatial;
518 if (!BLOCK_SYMMETRIC_MATRIX)
519 {
520 NEVER_REACHED; // to-come: non-mixed problems
521 //a_pressure_spatial = ComputeSpatialPressureMatrixTerm(quad_phi, grad_quad_phi, lin_phi, grad_lin_phi, x, &rElement);
522 }
523
524 c_matrix<double,PRESSURE_BLOCK_SIZE_ELEMENTAL,PRESSURE_BLOCK_SIZE_ELEMENTAL> a_pressure_pressure
525 = ComputePressurePressureMatrixTerm(linear_phi, grad_linear_phi, X, &rElement);
526
527 for (unsigned i=0; i<SPATIAL_BLOCK_SIZE_ELEMENTAL; i++)
528 {
529 for (unsigned j=0; j<SPATIAL_BLOCK_SIZE_ELEMENTAL; j++)
530 {
531 rAElem(i,j) += a_spatial_spatial(i,j)*wJ;
532 }
533
534 for (unsigned j=0; j<PRESSURE_BLOCK_SIZE_ELEMENTAL; j++)
535 {
536 rAElem(i, SPATIAL_BLOCK_SIZE_ELEMENTAL + j) += a_spatial_pressure(i,j)*wJ;
537 }
538 }
539
540 for (unsigned i=0; i<PRESSURE_BLOCK_SIZE_ELEMENTAL; i++)
541 {
542 if (BLOCK_SYMMETRIC_MATRIX)
543 {
544 for (unsigned j=0; j<SPATIAL_BLOCK_SIZE_ELEMENTAL; j++)
545 {
546 rAElem(SPATIAL_BLOCK_SIZE_ELEMENTAL + i, j) += a_spatial_pressure(j,i)*wJ;
547 }
548 }
549 else
550 {
551 NEVER_REACHED; // to-come: non-mixed problems
552 }
553
554 for (unsigned j=0; j<PRESSURE_BLOCK_SIZE_ELEMENTAL; j++)
555 {
556 rAElem(SPATIAL_BLOCK_SIZE_ELEMENTAL + i, SPATIAL_BLOCK_SIZE_ELEMENTAL + j) += a_pressure_pressure(i,j)*wJ;
557 }
558 }
559 }
560 }
561}
562
563#endif // ABSTRACTCONTINUUMMECHANICSASSEMBLER_HPP_
#define EXCEPTION(message)
#define NEVER_REACHED
virtual c_matrix< double, SPATIAL_BLOCK_SIZE_ELEMENTAL, PRESSURE_BLOCK_SIZE_ELEMENTAL > ComputeSpatialPressureMatrixTerm(c_vector< double, NUM_NODES_PER_ELEMENT > &rQuadPhi, c_matrix< double, DIM, NUM_NODES_PER_ELEMENT > &rGradQuadPhi, c_vector< double, NUM_VERTICES_PER_ELEMENT > &rLinearPhi, c_matrix< double, DIM, NUM_VERTICES_PER_ELEMENT > &rGradLinearPhi, c_vector< double, DIM > &rX, Element< DIM, DIM > *pElement)
virtual c_vector< double, SPATIAL_BLOCK_SIZE_ELEMENTAL > ComputeSpatialVectorTerm(c_vector< double, NUM_NODES_PER_ELEMENT > &rQuadPhi, c_matrix< double, DIM, NUM_NODES_PER_ELEMENT > &rGradQuadPhi, c_vector< double, DIM > &rX, Element< DIM, DIM > *pElement)=0
virtual c_vector< double, PRESSURE_BLOCK_SIZE_ELEMENTAL > ComputePressureVectorTerm(c_vector< double, NUM_VERTICES_PER_ELEMENT > &rLinearPhi, c_matrix< double, DIM, NUM_VERTICES_PER_ELEMENT > &rGradLinearPhi, c_vector< double, DIM > &rX, Element< DIM, DIM > *pElement)
void AssembleOnElement(Element< DIM, DIM > &rElement, c_matrix< double, STENCIL_SIZE, STENCIL_SIZE > &rAElem, c_vector< double, STENCIL_SIZE > &rBElem)
virtual c_matrix< double, SPATIAL_BLOCK_SIZE_ELEMENTAL, SPATIAL_BLOCK_SIZE_ELEMENTAL > ComputeSpatialSpatialMatrixTerm(c_vector< double, NUM_NODES_PER_ELEMENT > &rQuadPhi, c_matrix< double, DIM, NUM_NODES_PER_ELEMENT > &rGradQuadPhi, c_vector< double, DIM > &rX, Element< DIM, DIM > *pElement)
AbstractContinuumMechanicsAssembler(AbstractTetrahedralMesh< DIM, DIM > *pMesh)
virtual c_matrix< double, PRESSURE_BLOCK_SIZE_ELEMENTAL, PRESSURE_BLOCK_SIZE_ELEMENTAL > ComputePressurePressureMatrixTerm(c_vector< double, NUM_VERTICES_PER_ELEMENT > &rLinearPhi, c_matrix< double, DIM, NUM_VERTICES_PER_ELEMENT > &rGradLinearPhi, c_vector< double, DIM > &rX, Element< DIM, DIM > *pElement)
Node< SPACE_DIM > * GetNode(unsigned localIndex) const
unsigned GetNodeGlobalIndex(unsigned localIndex) const
bool GetOwnership() const
unsigned GetIndex() const
bool OptionExists(const std::string &rOption)
static CommandLineArguments * Instance()
static void ComputeBasisFunctions(const ChastePoint< ELEMENT_DIM > &rPoint, c_vector< double, ELEMENT_DIM+1 > &rReturnValue)
static void ComputeTransformedBasisFunctionDerivatives(const ChastePoint< ELEMENT_DIM > &rPoint, const c_matrix< double, ELEMENT_DIM, ELEMENT_DIM > &rInverseJacobian, c_matrix< double, ELEMENT_DIM, ELEMENT_DIM+1 > &rReturnValue)
static void Zero(Mat matrix)
static unsigned GetSize(Mat matrix)
static unsigned GetMyRank()
static unsigned GetSize(Vec vector)
static void Zero(Vec vector)
static void ComputeBasisFunctions(const ChastePoint< ELEMENT_DIM > &rPoint, c_vector< double,(ELEMENT_DIM+1) *(ELEMENT_DIM+2)/2 > &rReturnValue)
static void ComputeTransformedBasisFunctionDerivatives(const ChastePoint< ELEMENT_DIM > &rPoint, const c_matrix< double, ELEMENT_DIM, ELEMENT_DIM > &rInverseJacobian, c_matrix< double, ELEMENT_DIM,(ELEMENT_DIM+1) *(ELEMENT_DIM+2)/2 > &rReturnValue)