MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
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_QualityMetric *hm)
 
 TMOP_Integrator (TMOP_QualityMetric *m, TargetConstructor *tc)
 
 ~TMOP_Integrator ()
 
void ReleasePADeviceMemory (bool copy_to_host=true)
 
void SetIntegrationRules (IntegrationRules &irules, int order)
 Prescribe a set of integration rules; relevant for mixed meshes.
 
void SetInitialMeshPos (const GridFunction *x0)
 
void IntegrateOverTarget (bool integ_over_target_)
 
void SetCoefficient (Coefficient &w1)
 Sets a scaling Coefficient for the quality metric term of the integrator.
 
void EnableLimiting (const GridFunction &n0, const GridFunction &dist, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
 Limiting of the mesh displacements (general version).
 
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.
 
void EnableAdaptiveLimiting (const GridFunction &z0, Coefficient &coeff, AdaptivityEvaluator &ae)
 Restriction of the node positions to certain regions.
 
void EnableAdaptiveLimiting (const ParGridFunction &z0, Coefficient &coeff, AdaptivityEvaluator &ae)
 Parallel support for adaptive limiting.
 
void EnableSurfaceFitting (const GridFunction &s0, const Array< bool > &smarker, Coefficient &coeff, AdaptivityEvaluator &ae)
 Fitting of certain DOFs to the zero level set of a function.
 
void EnableSurfaceFitting (const ParGridFunction &s0, const Array< bool > &smarker, Coefficient &coeff, AdaptivityEvaluator &ae, AdaptivityEvaluator *aegrad=NULL, AdaptivityEvaluator *aehess=NULL)
 
void EnableSurfaceFittingFromSource (const ParGridFunction &s_bg, ParGridFunction &s0, const Array< bool > &smarker, Coefficient &coeff, AdaptivityEvaluator &ae, const ParGridFunction &s_bg_grad, ParGridFunction &s0_grad, AdaptivityEvaluator &age, const ParGridFunction &s_bg_hess, ParGridFunction &s0_hess, AdaptivityEvaluator &ahe)
 Fitting of certain DOFs in the current mesh to the zero level set of a function defined on another (finer) source mesh.
 
void EnableSurfaceFitting (const GridFunction &pos, const Array< bool > &smarker, Coefficient &coeff)
 Fitting of certain DOFs to given positions in physical space.
 
void GetSurfaceFittingErrors (const Vector &d_loc, real_t &err_avg, real_t &err_max)
 
bool IsSurfaceFittingEnabled ()
 
void SetLimitingNodes (const GridFunction &n0)
 Update the original/reference nodes used for limiting.
 
real_t GetElementEnergy (const FiniteElement &el, ElementTransformation &T, const Vector &d_el) override
 Computes the integral of W(Jacobian(Trt)) over a target zone.
 
virtual real_t GetRefinementElementEnergy (const FiniteElement &el, ElementTransformation &T, const Vector &elfun, const IntegrationRule &irule)
 Computes the mean of the energies of the given element's children.
 
virtual real_t GetDerefinementElementEnergy (const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
 
void AssembleElementVector (const FiniteElement &el, ElementTransformation &T, const Vector &d_el, Vector &elvect) override
 First defivative of GetElementEnergy() w.r.t. each local H1 DOF.
 
void AssembleElementGrad (const FiniteElement &el, ElementTransformation &T, const Vector &d_el, DenseMatrix &elmat) override
 Second derivative of GetElementEnergy() w.r.t. each local H1 DOF.
 
TMOP_QualityMetricGetAMRQualityMetric ()
 
void UpdateAfterMeshTopologyChange ()
 
void ParUpdateAfterMeshTopologyChange ()
 
void AssemblePA (const FiniteElementSpace &) override
 Method defining partial assembly.
 
void AssembleGradPA (const Vector &, const FiniteElementSpace &) override
 Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at the state x.
 
real_t GetLocalStateEnergyPA (const Vector &) const override
 Compute the local (to the MPI rank) energy with partial assembly.
 
void AddMultPA (const Vector &, Vector &) const override
 Method for partially assembled action.
 
void AddMultGradPA (const Vector &, Vector &) const override
 Method for partially assembled gradient action.
 
void AssembleGradDiagonalPA (Vector &) const override
 Method for computing the diagonal of the gradient with partial assembly.
 
DiscreteAdaptTCGetDiscreteAdaptTC () const
 
void EnableNormalization (const GridFunction &x)
 Computes the normalization factors of the metric and limiting integrals using the mesh position given by x.
 
void ParEnableNormalization (const ParGridFunction &x)
 
void EnableFiniteDifferences (const GridFunction &x)
 Enables FD-based approximation and computes dx.
 
void EnableFiniteDifferences (const ParGridFunction &x)
 
void SetFDhScale (real_t scale)
 
bool GetFDFlag () const
 
real_t GetFDh () const
 
void SetExactActionFlag (bool flag_)
 Flag to control if exact action of Integration is effected.
 
void UpdateSurfaceFittingWeight (real_t factor)
 Update the surface fitting weight as surf_fit_coeff *= factor;.
 
real_t GetSurfaceFittingWeight ()
 Get the surface fitting weight.
 
void ComputeUntangleMetricQuantiles (const Vector &d, const FiniteElementSpace &fes)
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
- Public Member Functions inherited from mfem::NonlinearFormIntegrator
void SetIntegrationMode (Mode m)
 
bool Patchwise () const
 
void SetPAMemoryType (MemoryType mt)
 
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.
 
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.
 
virtual bool SupportsCeed () const
 Indicates whether this integrator can use a Ceed backend.
 
virtual void AssembleMF (const FiniteElementSpace &fes)
 Method defining fully unassembled operator.
 
virtual void AddMultMF (const Vector &x, Vector &y) const
 
ceed::OperatorGetCeedOp ()
 
virtual ~NonlinearFormIntegrator ()
 
- Public Member Functions inherited from mfem::Integrator
 Integrator (const IntegrationRule *ir=NULL)
 Create a new Integrator, optionally providing a prescribed quadrature rule to use in assembly.
 
virtual void SetIntRule (const IntegrationRule *ir)
 Prescribe a fixed IntegrationRule to use, or set to null to let the integrator choose an appropriate rule.
 
void SetIntegrationRule (const IntegrationRule &ir)
 Prescribe a fixed IntegrationRule to use. Sets the NURBS patch integration rule to null.
 
void SetNURBSPatchIntRule (NURBSMeshRules *pr)
 Sets an integration rule for use on NURBS patches.
 
bool HasNURBSPatchIntRule () const
 Check if a NURBS patch integration rule has been set.
 
const IntegrationRuleGetIntRule () const
 Directly return the IntRule pointer (possibly null) without checking for NURBS patch rules or falling back on a default.
 
const IntegrationRuleGetIntegrationRule () const
 Equivalent to GetIntRule, but retained for backward compatibility with applications.
 

Protected Member Functions

void ComputeNormalizationEnergies (const GridFunction &x, real_t &metric_energy, real_t &lim_energy, real_t &surf_fit_gf_energy)
 
void AssembleElementVectorExact (const FiniteElement &el, ElementTransformation &T, const Vector &d_el, Vector &elvect)
 
void AssembleElementGradExact (const FiniteElement &el, ElementTransformation &T, const Vector &d_el, DenseMatrix &elmat)
 
void AssembleElementVectorFD (const FiniteElement &el, ElementTransformation &T, const Vector &d_el, Vector &elvect)
 
void AssembleElementGradFD (const FiniteElement &el, ElementTransformation &T, const Vector &d_el, DenseMatrix &elmat)
 
void AssembleElemVecAdaptLim (const FiniteElement &el, IsoparametricTransformation &Tpr, const IntegrationRule &ir, const Vector &weights, DenseMatrix &mat)
 
void AssembleElemGradAdaptLim (const FiniteElement &el, IsoparametricTransformation &Tpr, const IntegrationRule &ir, const Vector &weights, DenseMatrix &m)
 
void AssembleElemVecSurfFit (const FiniteElement &el_x, IsoparametricTransformation &Tpr, DenseMatrix &mat)
 
void AssembleElemGradSurfFit (const FiniteElement &el_x, IsoparametricTransformation &Tpr, DenseMatrix &mat)
 
real_t GetFDDerivative (const FiniteElement &el, ElementTransformation &T, Vector &d_el, const int nodenum, const int idir, const real_t baseenergy, bool update_stored)
 
void ComputeFDh (const Vector &d, const FiniteElementSpace &fes)
 Determines the perturbation, h, for FD-based approximation.
 
void ComputeMinJac (const Vector &x, const FiniteElementSpace &fes)
 
void UpdateAfterMeshPositionChange (const Vector &d, const FiniteElementSpace &d_fes)
 
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
 
real_t GetLocalStateEnergyPA_2D (const Vector &) const
 
real_t GetLocalStateEnergyPA_C0_2D (const Vector &) const
 
real_t GetLocalStateEnergyPA_3D (const Vector &) const
 
real_t 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
 
void UpdateCoefficientsPA (const Vector &d_loc)
 
real_t ComputeMinDetT (const Vector &x, const FiniteElementSpace &fes)
 
real_t ComputeUntanglerMaxMuBarrier (const Vector &x, const FiniteElementSpace &fes)
 
void RemapSurfaceFittingLevelSetAtNodes (const Vector &new_x, int new_x_ordering)
 
- Protected Member Functions inherited from mfem::NonlinearFormIntegrator
 NonlinearFormIntegrator (const IntegrationRule *ir=NULL)
 
- Protected Member Functions inherited from mfem::Integrator
const IntegrationRuleGetIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const
 Returns an integration rule based on the arguments and internal state of the Integrator object.
 
const IntegrationRuleGetIntegrationRule (const FiniteElement &el, const ElementTransformation &trans) const
 Returns an integration rule based on the arguments and internal state. (Version for identical trial_fe and test_fe)
 
virtual const IntegrationRuleGetDefaultIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const
 Subclasses should override to choose a default integration rule.
 

Protected Attributes

const GridFunctionx_0
 
bool periodic = false
 
TMOP_QualityMetrich_metric
 
TMOP_QualityMetricmetric
 
const TargetConstructortargetC
 
IntegrationRulesIntegRules
 
int integ_order
 
bool integ_over_target = true
 
Coefficientmetric_coeff
 
real_t metric_normal
 
const GridFunctionlim_nodes0
 
Coefficientlim_coeff
 
const GridFunctionlim_dist
 
TMOP_LimiterFunctionlim_func
 
real_t lim_normal
 
const GridFunctionadapt_lim_gf0
 
const ParGridFunctionadapt_lim_pgf0
 
GridFunctionadapt_lim_gf
 
Coefficientadapt_lim_coeff
 
AdaptivityEvaluatoradapt_lim_eval
 
const Array< bool > * surf_fit_marker
 
Coefficientsurf_fit_coeff
 
GridFunctionsurf_fit_gf
 
AdaptivityEvaluatorsurf_fit_eval
 
TMOP_QuadraticLimitersurf_fit_limiter
 
const GridFunctionsurf_fit_pos
 
real_t surf_fit_normal
 
GridFunctionsurf_fit_grad
 
GridFunctionsurf_fit_hess
 
AdaptivityEvaluatorsurf_fit_eval_grad
 
AdaptivityEvaluatorsurf_fit_eval_hess
 
Array< int > surf_fit_dof_count
 
Array< int > surf_fit_marker_dof_index
 
DiscreteAdaptTCdiscr_tc
 
bool fdflag
 
real_t fd_h
 
real_t fd_h_scale
 
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   XL 
 
   Vector   H 
 
   Vector   C0 
 
   Vector   LD 
 
   Vector   H0 
 
   Vector   MC 
 
   const DofToQuad *   maps 
 
   const DofToQuad *   maps_lim = nullptr 
 
   const GeometricFactors *   geom 
 
   const FiniteElementSpace *   fes 
 
   const IntegrationRule *   ir 
 
PA 
 
- Protected Attributes inherited from mfem::NonlinearFormIntegrator
Mode integrationMode = Mode::ELEMENTWISE
 
ceed::OperatorceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 
- Protected Attributes inherited from mfem::Integrator
const IntegrationRuleIntRule
 
NURBSMeshRulespatchRules = nullptr
 

Friends

class TMOPNewtonSolver
 
class TMOPComboIntegrator
 

Additional Inherited Members

- Public Types inherited from mfem::NonlinearFormIntegrator
enum  Mode { ELEMENTWISE = 0 , PATCHWISE = 1 , PATCHWISE_REDUCED = 2 }
 

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 1884 of file tmop.hpp.

Constructor & Destructor Documentation

◆ TMOP_Integrator() [1/2]

mfem::TMOP_Integrator::TMOP_Integrator ( TMOP_QualityMetric * m,
TargetConstructor * tc,
TMOP_QualityMetric * hm )
inline
Parameters
[in]mTMOP_QualityMetric for r-adaptivity (not owned).
[in]tcTarget-matrix construction algorithm to use (not owned).
[in]hmTMOP_QualityMetric for h-adaptivity (not owned).

Definition at line 2152 of file tmop.hpp.

◆ TMOP_Integrator() [2/2]

mfem::TMOP_Integrator::TMOP_Integrator ( TMOP_QualityMetric * m,
TargetConstructor * tc )
inline

Definition at line 2169 of file tmop.hpp.

◆ ~TMOP_Integrator()

mfem::TMOP_Integrator::~TMOP_Integrator ( )

Definition at line 3513 of file tmop.cpp.

Member Function Documentation

◆ ActionIntegrationRule()

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

Definition at line 2094 of file tmop.hpp.

◆ AddMultGradPA()

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

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.

Definition at line 379 of file tmop_pa.cpp.

◆ AddMultGradPA_2D()

void mfem::TMOP_Integrator::AddMultGradPA_2D ( const Vector & R,
Vector & C ) const
protected

Definition at line 114 of file tmop_pa_h2m.cpp.

◆ AddMultGradPA_3D()

void mfem::TMOP_Integrator::AddMultGradPA_3D ( const Vector & R,
Vector & C ) const
protected

Definition at line 117 of file tmop_pa_h3m.cpp.

◆ AddMultGradPA_C0_2D()

void mfem::TMOP_Integrator::AddMultGradPA_C0_2D ( const Vector & R,
Vector & C ) const
protected

Definition at line 94 of file tmop_pa_h2m_c0.cpp.

◆ AddMultGradPA_C0_3D()

void mfem::TMOP_Integrator::AddMultGradPA_C0_3D ( const Vector & R,
Vector & C ) const
protected

Definition at line 98 of file tmop_pa_h3m_c0.cpp.

◆ AddMultPA()

void mfem::TMOP_Integrator::AddMultPA ( const Vector & x,
Vector & y ) const
overridevirtual

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.

Definition at line 348 of file tmop_pa.cpp.

◆ AddMultPA_2D()

void mfem::TMOP_Integrator::AddMultPA_2D ( const Vector & X,
Vector & Y ) const
protected

Definition at line 203 of file tmop_pa_p2.cpp.

◆ AddMultPA_3D()

void mfem::TMOP_Integrator::AddMultPA_3D ( const Vector & X,
Vector & Y ) const
protected

Definition at line 239 of file tmop_pa_p3.cpp.

◆ AddMultPA_C0_2D()

void mfem::TMOP_Integrator::AddMultPA_C0_2D ( const Vector & X,
Vector & Y ) const
protected

Definition at line 152 of file tmop_pa_p2_c0.cpp.

◆ AddMultPA_C0_3D()

void mfem::TMOP_Integrator::AddMultPA_C0_3D ( const Vector & X,
Vector & Y ) const
protected

Definition at line 160 of file tmop_pa_p3_c0.cpp.

◆ AssembleDiagonalPA_2D()

void mfem::TMOP_Integrator::AssembleDiagonalPA_2D ( Vector & D) const
protected

Definition at line 182 of file tmop_pa_h2d.cpp.

◆ AssembleDiagonalPA_3D()

void mfem::TMOP_Integrator::AssembleDiagonalPA_3D ( Vector & D) const
protected

Definition at line 252 of file tmop_pa_h3d.cpp.

◆ AssembleDiagonalPA_C0_2D()

void mfem::TMOP_Integrator::AssembleDiagonalPA_C0_2D ( Vector & D) const
protected

Definition at line 84 of file tmop_pa_h2d_c0.cpp.

◆ AssembleDiagonalPA_C0_3D()

void mfem::TMOP_Integrator::AssembleDiagonalPA_C0_3D ( Vector & D) const
protected

Definition at line 111 of file tmop_pa_h3d_c0.cpp.

◆ AssembleElementGrad()

void mfem::TMOP_Integrator::AssembleElementGrad ( const FiniteElement & el,
ElementTransformation & T,
const Vector & d_el,
DenseMatrix & elmat )
overridevirtual

Second derivative of GetElementEnergy() w.r.t. each local H1 DOF.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 4259 of file tmop.cpp.

◆ AssembleElementGradExact()

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

Definition at line 4450 of file tmop.cpp.

◆ AssembleElementGradFD()

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

Definition at line 4972 of file tmop.cpp.

◆ AssembleElementVector()

void mfem::TMOP_Integrator::AssembleElementVector ( const FiniteElement & el,
ElementTransformation & T,
const Vector & d_el,
Vector & elvect )
overridevirtual

First defivative of GetElementEnergy() w.r.t. each local H1 DOF.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 4245 of file tmop.cpp.

◆ AssembleElementVectorExact()

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

Definition at line 4274 of file tmop.cpp.

◆ AssembleElementVectorFD()

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

Definition at line 4894 of file tmop.cpp.

◆ AssembleElemGradAdaptLim()

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

Definition at line 4615 of file tmop.cpp.

◆ AssembleElemGradSurfFit()

void mfem::TMOP_Integrator::AssembleElemGradSurfFit ( const FiniteElement & el_x,
IsoparametricTransformation & Tpr,
DenseMatrix & mat )
protected

Definition at line 4757 of file tmop.cpp.

◆ AssembleElemVecAdaptLim()

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

Definition at line 4579 of file tmop.cpp.

◆ AssembleElemVecSurfFit()

void mfem::TMOP_Integrator::AssembleElemVecSurfFit ( const FiniteElement & el_x,
IsoparametricTransformation & Tpr,
DenseMatrix & mat )
protected

Definition at line 4676 of file tmop.cpp.

◆ AssembleGradDiagonalPA()

void mfem::TMOP_Integrator::AssembleGradDiagonalPA ( Vector & diag) const
overridevirtual

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.

Definition at line 323 of file tmop_pa.cpp.

◆ AssembleGradPA()

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

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.

Definition at line 23 of file tmop_pa.cpp.

◆ AssembleGradPA_2D()

void mfem::TMOP_Integrator::AssembleGradPA_2D ( const Vector & X) const
protected

Definition at line 349 of file tmop_pa_h2s.cpp.

◆ AssembleGradPA_3D()

void mfem::TMOP_Integrator::AssembleGradPA_3D ( const Vector & X) const
protected

Definition at line 441 of file tmop_pa_h3s.cpp.

◆ AssembleGradPA_C0_2D()

void mfem::TMOP_Integrator::AssembleGradPA_C0_2D ( const Vector & X) const
protected

Definition at line 150 of file tmop_pa_h2s_c0.cpp.

◆ AssembleGradPA_C0_3D()

void mfem::TMOP_Integrator::AssembleGradPA_C0_3D ( const Vector & X) const
protected

Definition at line 169 of file tmop_pa_h3s_c0.cpp.

◆ AssemblePA() [1/2]

void mfem::TMOP_Integrator::AssemblePA ( const FiniteElementSpace & fes)
overridevirtual

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.

Definition at line 232 of file tmop_pa.cpp.

◆ AssemblePA() [2/2]

void mfem::NonlinearFormIntegrator::AssemblePA ( const FiniteElementSpace & trial_fes,
const FiniteElementSpace & test_fes )
virtual

The result of the partial assembly is stored internally so that it can be used later in the methods AddMultPA(). Used with BilinearFormIntegrators that have different spaces.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 99 of file nonlininteg.cpp.

◆ AssemblePA_Limiting()

void mfem::TMOP_Integrator::AssemblePA_Limiting ( )
protected

Definition at line 60 of file tmop_pa.cpp.

◆ ComputeAllElementTargets()

void mfem::TMOP_Integrator::ComputeAllElementTargets ( const Vector & xe = Vector()) const
protected

Definition at line 173 of file tmop_pa.cpp.

◆ ComputeFDh()

void mfem::TMOP_Integrator::ComputeFDh ( const Vector & d,
const FiniteElementSpace & fes )
protected

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

Definition at line 5412 of file tmop.cpp.

◆ ComputeMinDetT()

real_t mfem::TMOP_Integrator::ComputeMinDetT ( const Vector & x,
const FiniteElementSpace & fes )
protected

Definition at line 5488 of file tmop.cpp.

◆ ComputeMinJac()

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

Definition at line 5201 of file tmop.cpp.

◆ ComputeNormalizationEnergies()

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

Definition at line 5128 of file tmop.cpp.

◆ ComputeUntangleMetricQuantiles()

void mfem::TMOP_Integrator::ComputeUntangleMetricQuantiles ( const Vector & d,
const FiniteElementSpace & fes )

Computes quantiles needed for UntangleMetrics. Note that in parallel, the ParFiniteElementSpace must be passed as argument for consistency across MPI ranks.

Definition at line 5592 of file tmop.cpp.

◆ ComputeUntanglerMaxMuBarrier()

real_t mfem::TMOP_Integrator::ComputeUntanglerMaxMuBarrier ( const Vector & x,
const FiniteElementSpace & fes )
protected

Definition at line 5530 of file tmop.cpp.

◆ DisableLimiting()

void mfem::TMOP_Integrator::DisableLimiting ( )
inlineprotected

Definition at line 2079 of file tmop.hpp.

◆ EnableAdaptiveLimiting() [1/2]

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 term 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 3553 of file tmop.cpp.

◆ EnableAdaptiveLimiting() [2/2]

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

Parallel support for adaptive limiting.

Definition at line 3570 of file tmop.cpp.

◆ EnableFiniteDifferences() [1/2]

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

Enables FD-based approximation and computes dx.

Definition at line 5434 of file tmop.cpp.

◆ EnableFiniteDifferences() [2/2]

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

Definition at line 5461 of file tmop.cpp.

◆ EnableLimiting() [1/2]

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 3542 of file tmop.cpp.

◆ EnableLimiting() [2/2]

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 3528 of file tmop.cpp.

◆ EnableNormalization()

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 5104 of file tmop.cpp.

◆ EnableSurfaceFitting() [1/3]

void mfem::TMOP_Integrator::EnableSurfaceFitting ( const GridFunction & pos,
const Array< bool > & smarker,
Coefficient & coeff )

Fitting of certain DOFs to given positions in physical space.

Having a set S of marked nodes (or DOFs) and their target positions in physical space x_t, we move these nodes to the target positions during the optimization process. This function adds to the TMOP functional the term \( \sum_{i \in S} c \frac{1}{2} (x_i - x_{t,i})^2 \), where \(c\) corresponds to coeff below and is evaluated at the DOF locations.

Parameters
[in]posThe desired positions for the mesh nodes.
[in]smarkerIndicates which DOFs will be aligned.
[in]coeffCoefficient c for the above integral.

Definition at line 3622 of file tmop.cpp.

◆ EnableSurfaceFitting() [2/3]

void mfem::TMOP_Integrator::EnableSurfaceFitting ( const GridFunction & s0,
const Array< bool > & smarker,
Coefficient & coeff,
AdaptivityEvaluator & ae )

Fitting of certain DOFs to the zero level set of a function.

Having a level set function s0(x0) on the starting mesh, and a set of marked nodes (or DOFs), we move these nodes to the zero level set of s0. If s(x) is the image of s0(x0) on the current mesh, this function adds to the TMOP functional the term \( \int c \bar{s}(x))^2 \), where \(\bar{s}(x)\) is the restriction of s(x) on the aligned DOFs. Minimizing this term means that a marked node at x0 is allowed to move to a position x(x0) only if s(x) ~ 0. Such term can be used for surface fitting and tangential relaxation.

Parameters
[in]s0The level set function on the initial mesh.
[in]smarkerIndicates which DOFs will be aligned.
[in]coeffCoefficient c for the above integral.
[in]aeAdaptivityEvaluator to compute s(x) from s0(x0).

Definition at line 3588 of file tmop.cpp.

◆ EnableSurfaceFitting() [3/3]

void mfem::TMOP_Integrator::EnableSurfaceFitting ( const ParGridFunction & s0,
const Array< bool > & smarker,
Coefficient & coeff,
AdaptivityEvaluator & ae,
AdaptivityEvaluator * aegrad = NULL,
AdaptivityEvaluator * aehess = NULL )

Parallel support for surface fitting to the zero level set of a function. Here, we add two optional inputs: aegrad and aehess. When provided, the first and second derivative of the input level set are computed on the initial mesh, and aegrad and aehess are used to remap grad_s(x) from grad_s0(x0) and hess_s(x) from hess_s0(x0), respectively.

Definition at line 3647 of file tmop.cpp.

◆ EnableSurfaceFittingFromSource()

void mfem::TMOP_Integrator::EnableSurfaceFittingFromSource ( const ParGridFunction & s_bg,
ParGridFunction & s0,
const Array< bool > & smarker,
Coefficient & coeff,
AdaptivityEvaluator & ae,
const ParGridFunction & s_bg_grad,
ParGridFunction & s0_grad,
AdaptivityEvaluator & age,
const ParGridFunction & s_bg_hess,
ParGridFunction & s0_hess,
AdaptivityEvaluator & ahe )

Fitting of certain DOFs in the current mesh to the zero level set of a function defined on another (finer) source mesh.

Having a level set function s_bg(x_bg) on a source/background mesh, a set of marked nodes (or DOFs) in the current mesh, we move the marked nodes to the zero level set of s_bg. This functionality is used for surface fitting and tangential relaxation.

Parameters
[in]s_bgThe level set function on the background mesh.
[in]s0The level set function (automatically) interpolated on the initial mesh.
[in]smarkerMarker for aligned DOFs in the current mesh.
[in]coeffCoefficient c for the fitting penalty term.
[in]aeInterpolates s(x) from s_bg(x_bg).
[in]s_bg_gradGradient of s_bg on the background mesh.
[in]s0_gradGradient of s0 on the initial mesh.
[in]ageInterpolates s_grad(x) from s_bg_grad(x_bg).
[in]s_bg_hessHessian of s(x) on the background mesh.
[in]s0_hessHessian of s0 on the initial mesh.
[in]aheInterpolates s_hess(x) from s_bg_hess(x_bg). See the pmesh-fitting miniapp for details on usage.

Definition at line 3754 of file tmop.cpp.

◆ EnergyIntegrationRule()

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

Definition at line 2085 of file tmop.hpp.

◆ GetAMRQualityMetric()

TMOP_QualityMetric & mfem::TMOP_Integrator::GetAMRQualityMetric ( )
inline

Definition at line 2374 of file tmop.hpp.

◆ GetDerefinementElementEnergy()

real_t mfem::TMOP_Integrator::GetDerefinementElementEnergy ( const FiniteElement & el,
ElementTransformation & T,
const Vector & elfun )
virtual

This function is similar to GetElementEnergy, but ignores components such as limiting etc. to compute the element energy.

Definition at line 4191 of file tmop.cpp.

◆ GetDiscreteAdaptTC()

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

Definition at line 2395 of file tmop.hpp.

◆ GetElementEnergy()

real_t mfem::TMOP_Integrator::GetElementEnergy ( const FiniteElement & el,
ElementTransformation & T,
const Vector & d_el )
overridevirtual

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

Parameters
[in]elType of FiniteElement (TMOP assumes H1 elements).
[in]TMesh element transformation.
[in]d_elPhysical displacement of the zone w.r.t. x_0.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 3926 of file tmop.cpp.

◆ GetFDDerivative()

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

Definition at line 4872 of file tmop.cpp.

◆ GetFDFlag()

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

Definition at line 2411 of file tmop.hpp.

◆ GetFDh()

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

Definition at line 2412 of file tmop.hpp.

◆ GetLocalStateEnergyPA()

real_t mfem::TMOP_Integrator::GetLocalStateEnergyPA ( const Vector & x) const
overridevirtual

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.

Definition at line 404 of file tmop_pa.cpp.

◆ GetLocalStateEnergyPA_2D()

real_t mfem::TMOP_Integrator::GetLocalStateEnergyPA_2D ( const Vector & X) const
protected

Definition at line 171 of file tmop_pa_w2.cpp.

◆ GetLocalStateEnergyPA_3D()

real_t mfem::TMOP_Integrator::GetLocalStateEnergyPA_3D ( const Vector & X) const
protected

Definition at line 184 of file tmop_pa_w3.cpp.

◆ GetLocalStateEnergyPA_C0_2D()

real_t mfem::TMOP_Integrator::GetLocalStateEnergyPA_C0_2D ( const Vector & X) const
protected

Definition at line 131 of file tmop_pa_w2_c0.cpp.

◆ GetLocalStateEnergyPA_C0_3D()

real_t mfem::TMOP_Integrator::GetLocalStateEnergyPA_C0_3D ( const Vector & X) const
protected

Definition at line 143 of file tmop_pa_w3_c0.cpp.

◆ GetRefinementElementEnergy()

real_t mfem::TMOP_Integrator::GetRefinementElementEnergy ( const FiniteElement & el,
ElementTransformation & T,
const Vector & elfun,
const IntegrationRule & irule )
virtual

Computes the mean of the energies of the given element's children.

In addition to the inputs for GetElementEnergy, this function requires an IntegrationRule to be specified that will give the decomposition of the given element based on the refinement type being considered.

Definition at line 4105 of file tmop.cpp.

◆ GetSurfaceFittingErrors()

void mfem::TMOP_Integrator::GetSurfaceFittingErrors ( const Vector & d_loc,
real_t & err_avg,
real_t & err_max )

Definition at line 3829 of file tmop.cpp.

◆ GetSurfaceFittingWeight()

real_t mfem::TMOP_Integrator::GetSurfaceFittingWeight ( )

Get the surface fitting weight.

Definition at line 5093 of file tmop.cpp.

◆ GradientIntegrationRule()

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

Definition at line 2099 of file tmop.hpp.

◆ IntegrateOverTarget()

void mfem::TMOP_Integrator::IntegrateOverTarget ( bool integ_over_target_)
inline

The TMOP integrals can be computed over the reference element or the target elements. This function is used to switch between the two options. By default integration is performed over the target elements.

Definition at line 2196 of file tmop.hpp.

◆ IsSurfaceFittingEnabled()

bool mfem::TMOP_Integrator::IsSurfaceFittingEnabled ( )
inline

Definition at line 2332 of file tmop.hpp.

◆ ParEnableNormalization()

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

Definition at line 5114 of file tmop.cpp.

◆ ParUpdateAfterMeshTopologyChange()

void mfem::TMOP_Integrator::ParUpdateAfterMeshTopologyChange ( )

Definition at line 3913 of file tmop.cpp.

◆ ReleasePADeviceMemory()

void mfem::TMOP_Integrator::ReleasePADeviceMemory ( bool copy_to_host = true)

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

Definition at line 3482 of file tmop.cpp.

◆ RemapSurfaceFittingLevelSetAtNodes()

void mfem::TMOP_Integrator::RemapSurfaceFittingLevelSetAtNodes ( const Vector & new_x,
int new_x_ordering )
protected

Definition at line 5235 of file tmop.cpp.

◆ SetCoefficient()

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 2211 of file tmop.hpp.

◆ SetExactActionFlag()

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

Flag to control if exact action of Integration is effected.

Definition at line 2415 of file tmop.hpp.

◆ SetFDhScale()

void mfem::TMOP_Integrator::SetFDhScale ( real_t scale)
inline

Definition at line 2410 of file tmop.hpp.

◆ SetInitialMeshPos()

void mfem::TMOP_Integrator::SetInitialMeshPos ( const GridFunction * x0)

As the integrator operates on mesh displacements, this function is needed to set the initial mesh positions. For periodic meshes, the function is called with L2 positions. Called with nullptr to unset the x_0 after the problem is solved.

Definition at line 3496 of file tmop.cpp.

◆ SetIntegrationRules()

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 2180 of file tmop.hpp.

◆ SetLimitingNodes()

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

Update the original/reference nodes used for limiting.

Definition at line 2338 of file tmop.hpp.

◆ UpdateAfterMeshPositionChange()

void mfem::TMOP_Integrator::UpdateAfterMeshPositionChange ( const Vector & d,
const FiniteElementSpace & d_fes )
protected

Definition at line 5363 of file tmop.cpp.

◆ UpdateAfterMeshTopologyChange()

void mfem::TMOP_Integrator::UpdateAfterMeshTopologyChange ( )

Definition at line 3900 of file tmop.cpp.

◆ UpdateCoefficientsPA()

void mfem::TMOP_Integrator::UpdateCoefficientsPA ( const Vector & d_loc)
protected

Definition at line 185 of file tmop_pa.cpp.

◆ UpdateSurfaceFittingWeight()

void mfem::TMOP_Integrator::UpdateSurfaceFittingWeight ( real_t factor)

Update the surface fitting weight as surf_fit_coeff *= factor;.

Definition at line 5081 of file tmop.cpp.

Friends And Related Symbol Documentation

◆ TMOPComboIntegrator

friend class TMOPComboIntegrator
friend

Definition at line 1888 of file tmop.hpp.

◆ TMOPNewtonSolver

friend class TMOPNewtonSolver
friend

Definition at line 1887 of file tmop.hpp.

Member Data Documentation

◆ adapt_lim_coeff

Coefficient* mfem::TMOP_Integrator::adapt_lim_coeff
protected

Definition at line 1930 of file tmop.hpp.

◆ adapt_lim_eval

AdaptivityEvaluator* mfem::TMOP_Integrator::adapt_lim_eval
protected

Definition at line 1931 of file tmop.hpp.

◆ adapt_lim_gf

GridFunction* mfem::TMOP_Integrator::adapt_lim_gf
protected

Definition at line 1929 of file tmop.hpp.

◆ adapt_lim_gf0

const GridFunction* mfem::TMOP_Integrator::adapt_lim_gf0
protected

Definition at line 1925 of file tmop.hpp.

◆ adapt_lim_pgf0

const ParGridFunction* mfem::TMOP_Integrator::adapt_lim_pgf0
protected

Definition at line 1927 of file tmop.hpp.

◆ C0

Vector mfem::TMOP_Integrator::C0

Definition at line 2019 of file tmop.hpp.

◆ dim

int mfem::TMOP_Integrator::dim

Definition at line 2015 of file tmop.hpp.

◆ discr_tc

DiscreteAdaptTC* mfem::TMOP_Integrator::discr_tc
protected

Definition at line 1948 of file tmop.hpp.

◆ DS

DenseMatrix mfem::TMOP_Integrator::DS
protected

Definition at line 1974 of file tmop.hpp.

◆ DSh

DenseMatrix mfem::TMOP_Integrator::DSh
protected

Definition at line 1974 of file tmop.hpp.

◆ E

Vector mfem::TMOP_Integrator::E
mutable

Definition at line 2019 of file tmop.hpp.

◆ ElemDer

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

Definition at line 1961 of file tmop.hpp.

◆ ElemPertEnergy

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

Definition at line 1962 of file tmop.hpp.

◆ enabled

bool mfem::TMOP_Integrator::enabled

Definition at line 2014 of file tmop.hpp.

◆ exact_action

bool mfem::TMOP_Integrator::exact_action
protected

Definition at line 1959 of file tmop.hpp.

◆ fd_call_flag

bool mfem::TMOP_Integrator::fd_call_flag
protected

Definition at line 1956 of file tmop.hpp.

◆ fd_h

real_t mfem::TMOP_Integrator::fd_h
protected

Definition at line 1952 of file tmop.hpp.

◆ fd_h_scale

real_t mfem::TMOP_Integrator::fd_h_scale
protected

Definition at line 1953 of file tmop.hpp.

◆ fdflag

bool mfem::TMOP_Integrator::fdflag
protected

Definition at line 1951 of file tmop.hpp.

◆ fes

const FiniteElementSpace* mfem::TMOP_Integrator::fes

Definition at line 2023 of file tmop.hpp.

◆ geom

const GeometricFactors* mfem::TMOP_Integrator::geom

Definition at line 2022 of file tmop.hpp.

◆ H

Vector mfem::TMOP_Integrator::H

Definition at line 2019 of file tmop.hpp.

◆ H0

Vector mfem::TMOP_Integrator::H0

Definition at line 2019 of file tmop.hpp.

◆ h_metric

TMOP_QualityMetric* mfem::TMOP_Integrator::h_metric
protected

Definition at line 1898 of file tmop.hpp.

◆ integ_order

int mfem::TMOP_Integrator::integ_order
protected

Definition at line 1904 of file tmop.hpp.

◆ integ_over_target

bool mfem::TMOP_Integrator::integ_over_target = true
protected

Definition at line 1905 of file tmop.hpp.

◆ IntegRules

IntegrationRules* mfem::TMOP_Integrator::IntegRules
protected

Definition at line 1903 of file tmop.hpp.

◆ ir

const IntegrationRule* mfem::TMOP_Integrator::ir

Definition at line 2024 of file tmop.hpp.

◆ Jpr

DenseMatrix mfem::TMOP_Integrator::Jpr
protected

Definition at line 1974 of file tmop.hpp.

◆ Jpt

DenseMatrix mfem::TMOP_Integrator::Jpt
protected

Definition at line 1974 of file tmop.hpp.

◆ Jrt

DenseMatrix mfem::TMOP_Integrator::Jrt
protected

Definition at line 1974 of file tmop.hpp.

◆ Jtr

DenseTensor mfem::TMOP_Integrator::Jtr
mutable

Definition at line 2016 of file tmop.hpp.

◆ Jtr_debug_grad

bool mfem::TMOP_Integrator::Jtr_debug_grad
mutable

Definition at line 2018 of file tmop.hpp.

◆ Jtr_needs_update

bool mfem::TMOP_Integrator::Jtr_needs_update
mutable

Definition at line 2017 of file tmop.hpp.

◆ LD

Vector mfem::TMOP_Integrator::LD

Definition at line 2019 of file tmop.hpp.

◆ lim_coeff

Coefficient* mfem::TMOP_Integrator::lim_coeff
protected

Definition at line 1916 of file tmop.hpp.

◆ lim_dist

const GridFunction* mfem::TMOP_Integrator::lim_dist
protected

Definition at line 1918 of file tmop.hpp.

◆ lim_func

TMOP_LimiterFunction* mfem::TMOP_Integrator::lim_func
protected

Definition at line 1920 of file tmop.hpp.

◆ lim_nodes0

const GridFunction* mfem::TMOP_Integrator::lim_nodes0
protected

Definition at line 1915 of file tmop.hpp.

◆ lim_normal

real_t mfem::TMOP_Integrator::lim_normal
protected

Definition at line 1922 of file tmop.hpp.

◆ maps

const DofToQuad* mfem::TMOP_Integrator::maps

Definition at line 2020 of file tmop.hpp.

◆ maps_lim

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

Definition at line 2021 of file tmop.hpp.

◆ MC

Vector mfem::TMOP_Integrator::MC

Definition at line 2019 of file tmop.hpp.

◆ metric

TMOP_QualityMetric* mfem::TMOP_Integrator::metric
protected

Definition at line 1899 of file tmop.hpp.

◆ metric_coeff

Coefficient* mfem::TMOP_Integrator::metric_coeff
protected

Definition at line 1908 of file tmop.hpp.

◆ metric_normal

real_t mfem::TMOP_Integrator::metric_normal
protected

Definition at line 1910 of file tmop.hpp.

◆ ne

int mfem::TMOP_Integrator::ne

Definition at line 2015 of file tmop.hpp.

◆ nq

int mfem::TMOP_Integrator::nq

Definition at line 2015 of file tmop.hpp.

◆ O

Vector mfem::TMOP_Integrator::O

Definition at line 2019 of file tmop.hpp.

◆ P

DenseMatrix mfem::TMOP_Integrator::P
protected

Definition at line 1974 of file tmop.hpp.

◆ [struct]

struct { ... } mfem::TMOP_Integrator::PA

◆ periodic

bool mfem::TMOP_Integrator::periodic = false
protected

Definition at line 1896 of file tmop.hpp.

◆ PMatI

DenseMatrix mfem::TMOP_Integrator::PMatI
protected

Definition at line 1974 of file tmop.hpp.

◆ PMatO

DenseMatrix mfem::TMOP_Integrator::PMatO
protected

Definition at line 1974 of file tmop.hpp.

◆ surf_fit_coeff

Coefficient* mfem::TMOP_Integrator::surf_fit_coeff
protected

Definition at line 1935 of file tmop.hpp.

◆ surf_fit_dof_count

Array<int> mfem::TMOP_Integrator::surf_fit_dof_count
protected

Definition at line 1945 of file tmop.hpp.

◆ surf_fit_eval

AdaptivityEvaluator* mfem::TMOP_Integrator::surf_fit_eval
protected

Definition at line 1938 of file tmop.hpp.

◆ surf_fit_eval_grad

AdaptivityEvaluator* mfem::TMOP_Integrator::surf_fit_eval_grad
protected

Definition at line 1944 of file tmop.hpp.

◆ surf_fit_eval_hess

AdaptivityEvaluator * mfem::TMOP_Integrator::surf_fit_eval_hess
protected

Definition at line 1944 of file tmop.hpp.

◆ surf_fit_gf

GridFunction* mfem::TMOP_Integrator::surf_fit_gf
protected

Definition at line 1937 of file tmop.hpp.

◆ surf_fit_grad

GridFunction* mfem::TMOP_Integrator::surf_fit_grad
protected

Definition at line 1943 of file tmop.hpp.

◆ surf_fit_hess

GridFunction * mfem::TMOP_Integrator::surf_fit_hess
protected

Definition at line 1943 of file tmop.hpp.

◆ surf_fit_limiter

TMOP_QuadraticLimiter* mfem::TMOP_Integrator::surf_fit_limiter
protected

Definition at line 1940 of file tmop.hpp.

◆ surf_fit_marker

const Array<bool>* mfem::TMOP_Integrator::surf_fit_marker
protected

Definition at line 1934 of file tmop.hpp.

◆ surf_fit_marker_dof_index

Array<int> mfem::TMOP_Integrator::surf_fit_marker_dof_index
protected

Definition at line 1946 of file tmop.hpp.

◆ surf_fit_normal

real_t mfem::TMOP_Integrator::surf_fit_normal
protected

Definition at line 1942 of file tmop.hpp.

◆ surf_fit_pos

const GridFunction* mfem::TMOP_Integrator::surf_fit_pos
protected

Definition at line 1941 of file tmop.hpp.

◆ targetC

const TargetConstructor* mfem::TMOP_Integrator::targetC
protected

Definition at line 1900 of file tmop.hpp.

◆ X0

Vector mfem::TMOP_Integrator::X0

Definition at line 2019 of file tmop.hpp.

◆ x_0

const GridFunction* mfem::TMOP_Integrator::x_0
protected

Definition at line 1895 of file tmop.hpp.

◆ XL

Vector mfem::TMOP_Integrator::XL

Definition at line 2019 of file tmop.hpp.


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