MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Member Functions | Protected Attributes | Friends | List of all members
mfem::TMOP_Integrator Class Reference

A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor. More...

#include <tmop.hpp>

Inheritance diagram for mfem::TMOP_Integrator:
[legend]
Collaboration diagram for mfem::TMOP_Integrator:
[legend]

Public Member Functions

 TMOP_Integrator (TMOP_QualityMetric *m, TargetConstructor *tc)
 
 ~TMOP_Integrator ()
 
void ReleasePADeviceMemory ()
 
void SetIntegrationRules (IntegrationRules &irules, int order)
 Prescribe a set of integration rules; relevant for mixed meshes. More...
 
void SetCoefficient (Coefficient &w1)
 Sets a scaling Coefficient for the quality metric term of the integrator. More...
 
void EnableLimiting (const GridFunction &n0, const GridFunction &dist, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
 Limiting of the mesh displacements (general version). More...
 
void EnableLimiting (const GridFunction &n0, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
 Adds a limiting term to the integrator with limiting distance function (dist in the general version of the method) equal to 1. More...
 
void EnableAdaptiveLimiting (const GridFunction &z0, Coefficient &coeff, AdaptivityEvaluator &ae)
 Restriction of the node positions to certain regions. More...
 
void EnableAdaptiveLimiting (const ParGridFunction &z0, Coefficient &coeff, AdaptivityEvaluator &ae)
 Parallel support for adaptive limiting. More...
 
void SetLimitingNodes (const GridFunction &n0)
 Update the original/reference nodes used for limiting. More...
 
virtual double GetElementEnergy (const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
 Computes the integral of W(Jacobian(Trt)) over a target zone. More...
 
virtual void AssembleElementVector (const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
 Perform the local action of the NonlinearFormIntegrator. More...
 
virtual void AssembleElementGrad (const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
 Assemble the local gradient matrix. More...
 
virtual void AssemblePA (const FiniteElementSpace &)
 Method defining partial assembly. More...
 
virtual void AssembleGradPA (const Vector &, const FiniteElementSpace &)
 Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at the state x. More...
 
virtual double GetLocalStateEnergyPA (const Vector &) const
 Compute the local (to the MPI rank) energy with partial assembly. More...
 
virtual void AddMultPA (const Vector &, Vector &) const
 Method for partially assembled action. More...
 
virtual void AddMultGradPA (const Vector &, Vector &) const
 Method for partially assembled gradient action. More...
 
virtual void AssembleGradDiagonalPA (Vector &) const
 Method for computing the diagonal of the gradient with partial assembly. More...
 
DiscreteAdaptTCGetDiscreteAdaptTC () const
 
void EnableNormalization (const GridFunction &x)
 Computes the normalization factors of the metric and limiting integrals using the mesh position given by x. More...
 
void ParEnableNormalization (const ParGridFunction &x)
 
void EnableFiniteDifferences (const GridFunction &x)
 Enables FD-based approximation and computes dx. More...
 
void EnableFiniteDifferences (const ParGridFunction &x)
 
void SetFDhScale (double dxscale_)
 
bool GetFDFlag () const
 
double GetFDh () const
 
void SetExactActionFlag (bool flag_)
 Flag to control if exact action of Integration is effected. More...
 
- Public Member Functions inherited from mfem::NonlinearFormIntegrator
virtual void SetIntRule (const IntegrationRule *ir)
 Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == NULL). More...
 
void SetIntegrationRule (const IntegrationRule &ir)
 Prescribe a fixed IntegrationRule to use. More...
 
void SetPAMemoryType (MemoryType mt)
 
const IntegrationRuleGetIntegrationRule () const
 Get the integration rule of the integrator (possibly NULL). More...
 
virtual void AssembleFaceVector (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect)
 Perform the local action of the NonlinearFormIntegrator resulting from a face integral term. More...
 
virtual void AssembleFaceGrad (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, DenseMatrix &elmat)
 Assemble the local action of the gradient of the NonlinearFormIntegrator resulting from a face integral term. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
virtual bool SupportsCeed () const
 Indicates whether this integrator can use a Ceed backend. More...
 
virtual void AssembleMF (const FiniteElementSpace &fes)
 Method defining fully unassembled operator. More...
 
virtual void AddMultMF (const Vector &x, Vector &y) const
 
ceed::Operator & GetCeedOp ()
 
virtual ~NonlinearFormIntegrator ()
 

Protected Member Functions

void ComputeNormalizationEnergies (const GridFunction &x, double &metric_energy, double &lim_energy)
 
void AssembleElementVectorExact (const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
 
void AssembleElementGradExact (const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
 
void AssembleElementVectorFD (const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
 
void AssembleElementGradFD (const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
 
void AssembleElemVecAdaptLim (const FiniteElement &el, const Vector &weights, IsoparametricTransformation &Tpr, const IntegrationRule &ir, DenseMatrix &m)
 
void AssembleElemGradAdaptLim (const FiniteElement &el, const Vector &weights, IsoparametricTransformation &Tpr, const IntegrationRule &ir, DenseMatrix &m)
 
double GetFDDerivative (const FiniteElement &el, ElementTransformation &T, Vector &elfun, const int nodenum, const int idir, const double baseenergy, bool update_stored)
 
void ComputeFDh (const Vector &x, const FiniteElementSpace &fes)
 Determines the perturbation, h, for FD-based approximation. More...
 
void ComputeFDh (const Vector &x, const ParFiniteElementSpace &pfes)
 
void ComputeMinJac (const Vector &x, const FiniteElementSpace &fes)
 
void UpdateAfterMeshChange (const Vector &new_x)
 
void DisableLimiting ()
 
const IntegrationRuleEnergyIntegrationRule (const FiniteElement &el) const
 
const IntegrationRuleActionIntegrationRule (const FiniteElement &el) const
 
const IntegrationRuleGradientIntegrationRule (const FiniteElement &el) const
 
void AssembleGradPA_2D (const Vector &) const
 
void AssembleGradPA_3D (const Vector &) const
 
void AssembleGradPA_C0_2D (const Vector &) const
 
void AssembleGradPA_C0_3D (const Vector &) const
 
double GetLocalStateEnergyPA_2D (const Vector &) const
 
double GetLocalStateEnergyPA_C0_2D (const Vector &) const
 
double GetLocalStateEnergyPA_3D (const Vector &) const
 
double GetLocalStateEnergyPA_C0_3D (const Vector &) const
 
void AddMultPA_2D (const Vector &, Vector &) const
 
void AddMultPA_3D (const Vector &, Vector &) const
 
void AddMultPA_C0_2D (const Vector &, Vector &) const
 
void AddMultPA_C0_3D (const Vector &, Vector &) const
 
void AddMultGradPA_2D (const Vector &, Vector &) const
 
void AddMultGradPA_3D (const Vector &, Vector &) const
 
void AddMultGradPA_C0_2D (const Vector &, Vector &) const
 
void AddMultGradPA_C0_3D (const Vector &, Vector &) const
 
void AssembleDiagonalPA_2D (Vector &) const
 
void AssembleDiagonalPA_3D (Vector &) const
 
void AssembleDiagonalPA_C0_2D (Vector &) const
 
void AssembleDiagonalPA_C0_3D (Vector &) const
 
void AssemblePA_Limiting ()
 
void ComputeAllElementTargets (const Vector &xe=Vector()) const
 
- Protected Member Functions inherited from mfem::NonlinearFormIntegrator
 NonlinearFormIntegrator (const IntegrationRule *ir=NULL)
 

Protected Attributes

TMOP_QualityMetricmetric
 
const TargetConstructortargetC
 
IntegrationRulesIntegRules
 
int integ_order
 
Coefficientcoeff1
 
double metric_normal
 
const GridFunctionnodes0
 
Coefficientcoeff0
 
const GridFunctionlim_dist
 
TMOP_LimiterFunctionlim_func
 
double lim_normal
 
const GridFunctionzeta_0
 
GridFunctionzeta
 
Coefficientcoeff_zeta
 
AdaptivityEvaluatoradapt_eval
 
DiscreteAdaptTCdiscr_tc
 
bool fdflag
 
double dx
 
double dxscale
 
bool fd_call_flag
 
bool exact_action
 
Array< Vector * > ElemDer
 
Array< Vector * > ElemPertEnergy
 
DenseMatrix DSh
 
DenseMatrix DS
 
DenseMatrix Jrt
 
DenseMatrix Jpr
 
DenseMatrix Jpt
 
DenseMatrix P
 
DenseMatrix PMatI
 
DenseMatrix PMatO
 
struct {
   bool   enabled
 
   int   dim
 
   int   ne
 
   int   nq
 
   DenseTensor   Jtr
 
   bool   Jtr_needs_update
 
   bool   Jtr_debug_grad
 
   Vector   E
 
   Vector   O
 
   Vector   X0
 
   Vector   H
 
   Vector   C0
 
   Vector   LD
 
   Vector   H0
 
   const DofToQuad *   maps
 
   const DofToQuad *   maps_lim = nullptr
 
   const GeometricFactors *   geom
 
   const FiniteElementSpace *   fes
 
   const IntegrationRule *   ir
 
PA
 
- Protected Attributes inherited from mfem::NonlinearFormIntegrator
const IntegrationRuleIntRule
 
ceed::Operator * ceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 

Friends

class TMOPNewtonSolver
 
class TMOPComboIntegrator
 

Detailed Description

A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.

Represents \( \int W(Jpt) dx \) over a target zone, where W is the metric's strain energy density function, and Jpt is the Jacobian of the target->physical coordinates transformation. The virtual target zone is defined by the TargetConstructor.

Definition at line 1198 of file tmop.hpp.

Constructor & Destructor Documentation

mfem::TMOP_Integrator::TMOP_Integrator ( TMOP_QualityMetric m,
TargetConstructor tc 
)
inline
Parameters
[in]mTMOP_QualityMetric that will be integrated (not owned).
[in]tcTarget-matrix construction algorithm to use (not owned).

Definition at line 1400 of file tmop.hpp.

mfem::TMOP_Integrator::~TMOP_Integrator ( )

Definition at line 2203 of file tmop.cpp.

Member Function Documentation

const IntegrationRule& mfem::TMOP_Integrator::ActionIntegrationRule ( const FiniteElement el) const
inlineprotected

Definition at line 1357 of file tmop.hpp.

virtual void mfem::TMOP_Integrator::AddMultGradPA ( const Vector x,
Vector y 
) const
virtual

Method for partially assembled gradient action.

All arguments are E-vectors. This method can be called only after the method AssembleGradPA() has been called.

Parameters
[in]xThe gradient Operator is applied to the Vector x.
[in,out]yThe result Vector: \( y += G x \).

Reimplemented from mfem::NonlinearFormIntegrator.

void mfem::TMOP_Integrator::AddMultGradPA_2D ( const Vector ,
Vector  
) const
protected
void mfem::TMOP_Integrator::AddMultGradPA_3D ( const Vector ,
Vector  
) const
protected
void mfem::TMOP_Integrator::AddMultGradPA_C0_2D ( const Vector ,
Vector  
) const
protected
void mfem::TMOP_Integrator::AddMultGradPA_C0_3D ( const Vector ,
Vector  
) const
protected
virtual void mfem::TMOP_Integrator::AddMultPA ( const Vector x,
Vector y 
) const
virtual

Method for partially assembled action.

Perform the action of integrator on the input x and add the result to the output y. Both x and y are E-vectors, i.e. they represent the element-wise discontinuous version of the FE space.

This method can be called only after the method AssemblePA() has been called.

Reimplemented from mfem::NonlinearFormIntegrator.

void mfem::TMOP_Integrator::AddMultPA_2D ( const Vector ,
Vector  
) const
protected
void mfem::TMOP_Integrator::AddMultPA_3D ( const Vector ,
Vector  
) const
protected
void mfem::TMOP_Integrator::AddMultPA_C0_2D ( const Vector ,
Vector  
) const
protected
void mfem::TMOP_Integrator::AddMultPA_C0_3D ( const Vector ,
Vector  
) const
protected
void mfem::TMOP_Integrator::AssembleDiagonalPA_2D ( Vector ) const
protected
void mfem::TMOP_Integrator::AssembleDiagonalPA_3D ( Vector ) const
protected
void mfem::TMOP_Integrator::AssembleDiagonalPA_C0_2D ( Vector ) const
protected
void mfem::TMOP_Integrator::AssembleDiagonalPA_C0_3D ( Vector ) const
protected
void mfem::TMOP_Integrator::AssembleElementGrad ( const FiniteElement el,
ElementTransformation Tr,
const Vector elfun,
DenseMatrix elmat 
)
virtual

Assemble the local gradient matrix.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 2395 of file tmop.cpp.

void mfem::TMOP_Integrator::AssembleElementGradExact ( const FiniteElement el,
ElementTransformation T,
const Vector elfun,
DenseMatrix elmat 
)
protected

Definition at line 2548 of file tmop.cpp.

void mfem::TMOP_Integrator::AssembleElementGradFD ( const FiniteElement el,
ElementTransformation T,
const Vector elfun,
DenseMatrix elmat 
)
protected

Definition at line 2843 of file tmop.cpp.

void mfem::TMOP_Integrator::AssembleElementVector ( const FiniteElement el,
ElementTransformation Tr,
const Vector elfun,
Vector elvect 
)
virtual

Perform the local action of the NonlinearFormIntegrator.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 2381 of file tmop.cpp.

void mfem::TMOP_Integrator::AssembleElementVectorExact ( const FiniteElement el,
ElementTransformation T,
const Vector elfun,
Vector elvect 
)
protected

Definition at line 2410 of file tmop.cpp.

void mfem::TMOP_Integrator::AssembleElementVectorFD ( const FiniteElement el,
ElementTransformation T,
const Vector elfun,
Vector elvect 
)
protected

Definition at line 2779 of file tmop.cpp.

void mfem::TMOP_Integrator::AssembleElemGradAdaptLim ( const FiniteElement el,
const Vector weights,
IsoparametricTransformation Tpr,
const IntegrationRule ir,
DenseMatrix m 
)
protected

Definition at line 2693 of file tmop.cpp.

void mfem::TMOP_Integrator::AssembleElemVecAdaptLim ( const FiniteElement el,
const Vector weights,
IsoparametricTransformation Tpr,
const IntegrationRule ir,
DenseMatrix m 
)
protected

Definition at line 2654 of file tmop.cpp.

virtual void mfem::TMOP_Integrator::AssembleGradDiagonalPA ( Vector diag) const
virtual

Method for computing the diagonal of the gradient with partial assembly.

The result Vector diag is an E-Vector. This method can be called only after the method AssembleGradPA() has been called.

Parameters
[in,out]diagThe result Vector: \( diag += diag(G) \).

Reimplemented from mfem::NonlinearFormIntegrator.

virtual void mfem::TMOP_Integrator::AssembleGradPA ( const Vector x,
const FiniteElementSpace fes 
)
virtual

Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at the state x.

The result of the partial assembly is stored internally so that it can be used later in the methods AddMultGradPA() and AssembleGradDiagonalPA(). The state Vector x is an E-vector.

Reimplemented from mfem::NonlinearFormIntegrator.

void mfem::TMOP_Integrator::AssembleGradPA_2D ( const Vector ) const
protected
void mfem::TMOP_Integrator::AssembleGradPA_3D ( const Vector ) const
protected
void mfem::TMOP_Integrator::AssembleGradPA_C0_2D ( const Vector ) const
protected
void mfem::TMOP_Integrator::AssembleGradPA_C0_3D ( const Vector ) const
protected
virtual void mfem::TMOP_Integrator::AssemblePA ( const FiniteElementSpace fes)
virtual

Method defining partial assembly.

The result of the partial assembly is stored internally so that it can be used later in the methods AddMultPA().

Reimplemented from mfem::NonlinearFormIntegrator.

void mfem::TMOP_Integrator::AssemblePA_Limiting ( )
protected
void mfem::TMOP_Integrator::ComputeAllElementTargets ( const Vector xe = Vector()) const
protected
void mfem::TMOP_Integrator::ComputeFDh ( const Vector x,
const FiniteElementSpace fes 
)
protected

Determines the perturbation, h, for FD-based approximation.

Definition at line 3052 of file tmop.cpp.

void mfem::TMOP_Integrator::ComputeFDh ( const Vector x,
const ParFiniteElementSpace pfes 
)
protected

Definition at line 3059 of file tmop.cpp.

void mfem::TMOP_Integrator::ComputeMinJac ( const Vector x,
const FiniteElementSpace fes 
)
protected

Definition at line 3009 of file tmop.cpp.

void mfem::TMOP_Integrator::ComputeNormalizationEnergies ( const GridFunction x,
double &  metric_energy,
double &  lim_energy 
)
protected

Definition at line 2957 of file tmop.cpp.

void mfem::TMOP_Integrator::DisableLimiting ( )
inlineprotected

Definition at line 1342 of file tmop.hpp.

void mfem::TMOP_Integrator::EnableAdaptiveLimiting ( const GridFunction z0,
Coefficient coeff,
AdaptivityEvaluator ae 
)

Restriction of the node positions to certain regions.

Adds the term \( \int c (z(x) - z_0(x_0))^2 \), where z0(x0) is a given function on the starting mesh, and z(x) is its image on the new mesh. Minimizing this, means that a node at x0 is allowed to move to a position x(x0) only if z(x) ~ z0(x0). Such term can be used for tangential mesh relaxation.

Parameters
[in]z0Function z0 that controls the adaptive limiting.
[in]coeffCoefficient c for the above integral.
[in]aeAdaptivityEvaluator to compute z(x) from z0(x0).

Definition at line 2239 of file tmop.cpp.

void mfem::TMOP_Integrator::EnableAdaptiveLimiting ( const ParGridFunction z0,
Coefficient coeff,
AdaptivityEvaluator ae 
)

Parallel support for adaptive limiting.

Definition at line 2256 of file tmop.cpp.

void mfem::TMOP_Integrator::EnableFiniteDifferences ( const GridFunction x)

Enables FD-based approximation and computes dx.

Definition at line 3070 of file tmop.cpp.

void mfem::TMOP_Integrator::EnableFiniteDifferences ( const ParGridFunction x)

Definition at line 3084 of file tmop.cpp.

void mfem::TMOP_Integrator::EnableLimiting ( const GridFunction n0,
const GridFunction dist,
Coefficient w0,
TMOP_LimiterFunction lfunc = NULL 
)

Limiting of the mesh displacements (general version).

Adds the term \( \int w_0 f(x, x_0, d) dx \), where f is a measure of the displacement between x and x_0, given the max allowed displacement d.

Parameters
[in]n0Original mesh node coordinates (x0 above).
[in]distAllowed displacement in physical space (d above).
[in]w0Coefficient scaling the limiting integral.
[in]lfuncTMOP_LimiterFunction defining the function f. If NULL, a TMOP_QuadraticLimiter will be used. The TMOP_Integrator assumes ownership of this pointer.

Definition at line 2214 of file tmop.cpp.

void mfem::TMOP_Integrator::EnableLimiting ( const GridFunction n0,
Coefficient w0,
TMOP_LimiterFunction lfunc = NULL 
)

Adds a limiting term to the integrator with limiting distance function (dist in the general version of the method) equal to 1.

Definition at line 2228 of file tmop.cpp.

void mfem::TMOP_Integrator::EnableNormalization ( const GridFunction x)

Computes the normalization factors of the metric and limiting integrals using the mesh position given by x.

Definition at line 2938 of file tmop.cpp.

const IntegrationRule& mfem::TMOP_Integrator::EnergyIntegrationRule ( const FiniteElement el) const
inlineprotected

Definition at line 1348 of file tmop.hpp.

DiscreteAdaptTC* mfem::TMOP_Integrator::GetDiscreteAdaptTC ( ) const
inline

Definition at line 1503 of file tmop.hpp.

double mfem::TMOP_Integrator::GetElementEnergy ( const FiniteElement el,
ElementTransformation T,
const Vector elfun 
)
virtual

Computes the integral of W(Jacobian(Trt)) over a target zone.

Parameters
[in]elType of FiniteElement.
[in]TMesh element transformation.
[in]elfunPhysical coordinates of the zone.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 2273 of file tmop.cpp.

double mfem::TMOP_Integrator::GetFDDerivative ( const FiniteElement el,
ElementTransformation T,
Vector elfun,
const int  nodenum,
const int  idir,
const double  baseenergy,
bool  update_stored 
)
protected

Definition at line 2757 of file tmop.cpp.

bool mfem::TMOP_Integrator::GetFDFlag ( ) const
inline

Definition at line 1519 of file tmop.hpp.

double mfem::TMOP_Integrator::GetFDh ( ) const
inline

Definition at line 1520 of file tmop.hpp.

virtual double mfem::TMOP_Integrator::GetLocalStateEnergyPA ( const Vector x) const
virtual

Compute the local (to the MPI rank) energy with partial assembly.

Here the state x is an E-vector. This method can be called only after the method AssemblePA() has been called.

Reimplemented from mfem::NonlinearFormIntegrator.

double mfem::TMOP_Integrator::GetLocalStateEnergyPA_2D ( const Vector ) const
protected
double mfem::TMOP_Integrator::GetLocalStateEnergyPA_3D ( const Vector ) const
protected
double mfem::TMOP_Integrator::GetLocalStateEnergyPA_C0_2D ( const Vector ) const
protected
double mfem::TMOP_Integrator::GetLocalStateEnergyPA_C0_3D ( const Vector ) const
protected
const IntegrationRule& mfem::TMOP_Integrator::GradientIntegrationRule ( const FiniteElement el) const
inlineprotected

Definition at line 1362 of file tmop.hpp.

void mfem::TMOP_Integrator::ParEnableNormalization ( const ParGridFunction x)

Definition at line 2946 of file tmop.cpp.

void mfem::TMOP_Integrator::ReleasePADeviceMemory ( )

Release the device memory of large PA allocations. This will copy device memory back to the host before releasing.

Definition at line 2193 of file tmop.cpp.

void mfem::TMOP_Integrator::SetCoefficient ( Coefficient w1)
inline

Sets a scaling Coefficient for the quality metric term of the integrator.

With this addition, the integrator becomes \( \int w1 W(Jpt) dx \).

Note that the Coefficient is evaluated in the physical configuration and not in the target configuration which may be undefined.

Definition at line 1430 of file tmop.hpp.

void mfem::TMOP_Integrator::SetExactActionFlag ( bool  flag_)
inline

Flag to control if exact action of Integration is effected.

Definition at line 1523 of file tmop.hpp.

void mfem::TMOP_Integrator::SetFDhScale ( double  dxscale_)
inline

Definition at line 1518 of file tmop.hpp.

void mfem::TMOP_Integrator::SetIntegrationRules ( IntegrationRules irules,
int  order 
)
inline

Prescribe a set of integration rules; relevant for mixed meshes.

This function has priority over SetIntRule(), if both are called.

Definition at line 1418 of file tmop.hpp.

void mfem::TMOP_Integrator::SetLimitingNodes ( const GridFunction n0)
inline

Update the original/reference nodes used for limiting.

Definition at line 1471 of file tmop.hpp.

void mfem::TMOP_Integrator::UpdateAfterMeshChange ( const Vector new_x)
protected

Definition at line 3042 of file tmop.cpp.

Friends And Related Function Documentation

friend class TMOPComboIntegrator
friend

Definition at line 1202 of file tmop.hpp.

friend class TMOPNewtonSolver
friend

Definition at line 1201 of file tmop.hpp.

Member Data Documentation

AdaptivityEvaluator* mfem::TMOP_Integrator::adapt_eval
protected

Definition at line 1232 of file tmop.hpp.

Vector mfem::TMOP_Integrator::C0
mutable

Definition at line 1293 of file tmop.hpp.

Coefficient* mfem::TMOP_Integrator::coeff0
protected

Definition at line 1220 of file tmop.hpp.

Coefficient* mfem::TMOP_Integrator::coeff1
protected

Definition at line 1212 of file tmop.hpp.

Coefficient* mfem::TMOP_Integrator::coeff_zeta
protected

Definition at line 1231 of file tmop.hpp.

int mfem::TMOP_Integrator::dim

Definition at line 1289 of file tmop.hpp.

DiscreteAdaptTC* mfem::TMOP_Integrator::discr_tc
protected

Definition at line 1234 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::DS
protected

Definition at line 1260 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::DSh
protected

Definition at line 1260 of file tmop.hpp.

double mfem::TMOP_Integrator::dx
protected

Definition at line 1238 of file tmop.hpp.

double mfem::TMOP_Integrator::dxscale
protected

Definition at line 1239 of file tmop.hpp.

Vector mfem::TMOP_Integrator::E
mutable

Definition at line 1293 of file tmop.hpp.

Array<Vector *> mfem::TMOP_Integrator::ElemDer
protected

Definition at line 1247 of file tmop.hpp.

Array<Vector *> mfem::TMOP_Integrator::ElemPertEnergy
protected

Definition at line 1248 of file tmop.hpp.

bool mfem::TMOP_Integrator::enabled

Definition at line 1288 of file tmop.hpp.

bool mfem::TMOP_Integrator::exact_action
protected

Definition at line 1245 of file tmop.hpp.

bool mfem::TMOP_Integrator::fd_call_flag
protected

Definition at line 1242 of file tmop.hpp.

bool mfem::TMOP_Integrator::fdflag
protected

Definition at line 1237 of file tmop.hpp.

const FiniteElementSpace* mfem::TMOP_Integrator::fes

Definition at line 1297 of file tmop.hpp.

const GeometricFactors* mfem::TMOP_Integrator::geom

Definition at line 1296 of file tmop.hpp.

Vector mfem::TMOP_Integrator::H
mutable

Definition at line 1293 of file tmop.hpp.

Vector mfem::TMOP_Integrator::H0
mutable

Definition at line 1293 of file tmop.hpp.

int mfem::TMOP_Integrator::integ_order
protected

Definition at line 1209 of file tmop.hpp.

IntegrationRules* mfem::TMOP_Integrator::IntegRules
protected

Definition at line 1208 of file tmop.hpp.

const IntegrationRule* mfem::TMOP_Integrator::ir

Definition at line 1298 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::Jpr
protected

Definition at line 1260 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::Jpt
protected

Definition at line 1260 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::Jrt
protected

Definition at line 1260 of file tmop.hpp.

DenseTensor mfem::TMOP_Integrator::Jtr
mutable

Definition at line 1290 of file tmop.hpp.

bool mfem::TMOP_Integrator::Jtr_debug_grad
mutable

Definition at line 1292 of file tmop.hpp.

bool mfem::TMOP_Integrator::Jtr_needs_update
mutable

Definition at line 1291 of file tmop.hpp.

Vector mfem::TMOP_Integrator::LD
mutable

Definition at line 1293 of file tmop.hpp.

const GridFunction* mfem::TMOP_Integrator::lim_dist
protected

Definition at line 1222 of file tmop.hpp.

TMOP_LimiterFunction* mfem::TMOP_Integrator::lim_func
protected

Definition at line 1224 of file tmop.hpp.

double mfem::TMOP_Integrator::lim_normal
protected

Definition at line 1226 of file tmop.hpp.

const DofToQuad* mfem::TMOP_Integrator::maps

Definition at line 1294 of file tmop.hpp.

const DofToQuad* mfem::TMOP_Integrator::maps_lim = nullptr

Definition at line 1295 of file tmop.hpp.

TMOP_QualityMetric* mfem::TMOP_Integrator::metric
protected

Definition at line 1204 of file tmop.hpp.

double mfem::TMOP_Integrator::metric_normal
protected

Definition at line 1214 of file tmop.hpp.

int mfem::TMOP_Integrator::ne

Definition at line 1289 of file tmop.hpp.

const GridFunction* mfem::TMOP_Integrator::nodes0
protected

Definition at line 1219 of file tmop.hpp.

int mfem::TMOP_Integrator::nq

Definition at line 1289 of file tmop.hpp.

Vector mfem::TMOP_Integrator::O
mutable

Definition at line 1293 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::P
protected

Definition at line 1260 of file tmop.hpp.

struct { ... } mfem::TMOP_Integrator::PA
DenseMatrix mfem::TMOP_Integrator::PMatI
protected

Definition at line 1260 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::PMatO
protected

Definition at line 1260 of file tmop.hpp.

const TargetConstructor* mfem::TMOP_Integrator::targetC
protected

Definition at line 1205 of file tmop.hpp.

Vector mfem::TMOP_Integrator::X0
mutable

Definition at line 1293 of file tmop.hpp.

GridFunction* mfem::TMOP_Integrator::zeta
protected

Definition at line 1230 of file tmop.hpp.

const GridFunction* mfem::TMOP_Integrator::zeta_0
protected

Definition at line 1229 of file tmop.hpp.


The documentation for this class was generated from the following files: