Abstract class for local mesh quality metrics in the target-matrix optimization paradigm (TMOP) by P. Knupp et al.
More...
#include <tmop.hpp>
Inherits mfem::HyperelasticModel.
Inherited by mfem::TMOP_Metric_001, mfem::TMOP_Metric_002, mfem::TMOP_Metric_007, mfem::TMOP_Metric_009, mfem::TMOP_Metric_022, mfem::TMOP_Metric_050, mfem::TMOP_Metric_055, mfem::TMOP_Metric_056, mfem::TMOP_Metric_058, mfem::TMOP_Metric_077, mfem::TMOP_Metric_211, mfem::TMOP_Metric_252, mfem::TMOP_Metric_301, mfem::TMOP_Metric_302, mfem::TMOP_Metric_303, mfem::TMOP_Metric_315, mfem::TMOP_Metric_316, mfem::TMOP_Metric_321, and mfem::TMOP_Metric_352.
Abstract class for local mesh quality metrics in the target-matrix optimization paradigm (TMOP) by P. Knupp et al.
Definition at line 24 of file tmop.hpp.
mfem::TMOP_QualityMetric::TMOP_QualityMetric |
( |
| ) |
|
|
inline |
virtual mfem::TMOP_QualityMetric::~TMOP_QualityMetric |
( |
| ) |
|
|
inlinevirtual |
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the local gradient matrix 'A'.
- Parameters
-
[in] | Jpt | Represents the target->physical transformation Jacobian matrix. |
[in] | DS | Gradient of the basis matrix (dof x dim). |
[in] | weight | Quadrature weight coefficient for the point. |
[in,out] | A | Local gradient matrix where the contribution from this point will be added. |
Computes weight * d(dW_dxi)_d(xj) at the current point, for all i and j, where x1 ... xn are the FE dofs. This function is usually defined using the matrix invariants and their derivatives.
Implements mfem::HyperelasticModel.
Implemented in mfem::TMOP_Metric_352, mfem::TMOP_Metric_321, mfem::TMOP_Metric_316, mfem::TMOP_Metric_315, mfem::TMOP_Metric_303, mfem::TMOP_Metric_302, mfem::TMOP_Metric_301, mfem::TMOP_Metric_252, mfem::TMOP_Metric_211, mfem::TMOP_Metric_077, mfem::TMOP_Metric_058, mfem::TMOP_Metric_056, mfem::TMOP_Metric_055, mfem::TMOP_Metric_050, mfem::TMOP_Metric_022, mfem::TMOP_Metric_009, mfem::TMOP_Metric_007, mfem::TMOP_Metric_002, and mfem::TMOP_Metric_001.
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
- Parameters
-
[in] | Jpt | Represents the target->physical transformation Jacobian matrix. |
[out] | P | The evaluated 1st Piola-Kirchhoff stress tensor. |
Implements mfem::HyperelasticModel.
Implemented in mfem::TMOP_Metric_352, mfem::TMOP_Metric_321, mfem::TMOP_Metric_316, mfem::TMOP_Metric_315, mfem::TMOP_Metric_303, mfem::TMOP_Metric_302, mfem::TMOP_Metric_301, mfem::TMOP_Metric_252, mfem::TMOP_Metric_211, mfem::TMOP_Metric_077, mfem::TMOP_Metric_058, mfem::TMOP_Metric_056, mfem::TMOP_Metric_055, mfem::TMOP_Metric_050, mfem::TMOP_Metric_022, mfem::TMOP_Metric_009, mfem::TMOP_Metric_007, mfem::TMOP_Metric_002, and mfem::TMOP_Metric_001.
virtual double mfem::TMOP_QualityMetric::EvalW |
( |
const DenseMatrix & |
Jpt | ) |
const |
|
pure virtual |
Evaluate the strain energy density function, W = W(Jpt).
- Parameters
-
[in] | Jpt | Represents the target->physical transformation Jacobian matrix. |
Implements mfem::HyperelasticModel.
Implemented in mfem::TMOP_Metric_352, mfem::TMOP_Metric_321, mfem::TMOP_Metric_316, mfem::TMOP_Metric_315, mfem::TMOP_Metric_303, mfem::TMOP_Metric_302, mfem::TMOP_Metric_301, mfem::TMOP_Metric_252, mfem::TMOP_Metric_211, mfem::TMOP_Metric_077, mfem::TMOP_Metric_058, mfem::TMOP_Metric_056, mfem::TMOP_Metric_055, mfem::TMOP_Metric_050, mfem::TMOP_Metric_022, mfem::TMOP_Metric_009, mfem::TMOP_Metric_007, mfem::TMOP_Metric_002, and mfem::TMOP_Metric_001.
void mfem::TMOP_QualityMetric::SetTargetJacobian |
( |
const DenseMatrix & |
_Jtr | ) |
|
|
inline |
Specify the reference-element -> target-element Jacobian matrix for the point of interest.
The specified Jacobian matrix, Jtr, can be used by metrics that cannot be written just as a function of the target->physical Jacobian matrix, Jpt.
Definition at line 44 of file tmop.hpp.
Jacobian of the reference-element to target-element transformation.
Definition at line 27 of file tmop.hpp.
The documentation for this class was generated from the following file: