MFEM  v4.5.2
Finite element discretization library
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_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. 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 EnableSurfaceFitting (const GridFunction &s0, const Array< bool > &smarker, Coefficient &coeff, AdaptivityEvaluator &ae)
 Fitting of certain DOFs to the zero level set of a function. More...
 
void EnableSurfaceFitting (const ParGridFunction &s0, const Array< bool > &smarker, Coefficient &coeff, AdaptivityEvaluator &ae)
 Parallel support for surface fitting. More...
 
void GetSurfaceFittingErrors (double &err_avg, double &err_max)
 
bool IsSurfaceFittingEnabled ()
 
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 double GetRefinementElementEnergy (const FiniteElement &el, ElementTransformation &T, const Vector &elfun, const IntegrationRule &irule)
 Computes the mean of the energies of the given element's children. More...
 
virtual double GetDerefinementElementEnergy (const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
 
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...
 
TMOP_QualityMetricGetAMRQualityMetric ()
 
void UpdateAfterMeshTopologyChange ()
 
void ParUpdateAfterMeshTopologyChange ()
 
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...
 
void UpdateSurfaceFittingWeight (double factor)
 Update the surface fitting weight as surf_fit_coeff *= factor;. More...
 
double GetSurfaceFittingWeight ()
 Get the surface fitting weight. More...
 
void ComputeUntangleMetricQuantiles (const Vector &x, const FiniteElementSpace &fes)
 
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
- 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::OperatorGetCeedOp ()
 
virtual ~NonlinearFormIntegrator ()
 

Protected Member Functions

void ComputeNormalizationEnergies (const GridFunction &x, double &metric_energy, double &lim_energy, double &surf_fit_gf_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, 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)
 
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 ComputeMinJac (const Vector &x, const FiniteElementSpace &fes)
 
void UpdateAfterMeshPositionChange (const Vector &new_x, int new_x_ordering=Ordering::byNODES)
 
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
 
double ComputeMinDetT (const Vector &x, const FiniteElementSpace &fes)
 
double ComputeUntanglerMaxMuBarrier (const Vector &x, const FiniteElementSpace &fes)
 
- Protected Member Functions inherited from mfem::NonlinearFormIntegrator
 NonlinearFormIntegrator (const IntegrationRule *ir=NULL)
 

Protected Attributes

TMOP_QualityMetrich_metric
 
TMOP_QualityMetricmetric
 
const TargetConstructortargetC
 
IntegrationRulesIntegRules
 
int integ_order
 
Coefficientmetric_coeff
 
double metric_normal
 
const GridFunctionlim_nodes0
 
Coefficientlim_coeff
 
const GridFunctionlim_dist
 
TMOP_LimiterFunctionlim_func
 
double lim_normal
 
const GridFunctionadapt_lim_gf0
 
const ParGridFunctionadapt_lim_pgf0
 
GridFunctionadapt_lim_gf
 
Coefficientadapt_lim_coeff
 
AdaptivityEvaluatoradapt_lim_eval
 
GridFunctionsurf_fit_gf
 
const Array< bool > * surf_fit_marker
 
Coefficientsurf_fit_coeff
 
AdaptivityEvaluatorsurf_fit_eval
 
double surf_fit_normal
 
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::OperatorceedOp
 
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 1638 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 1870 of file tmop.hpp.

◆ TMOP_Integrator() [2/2]

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

Definition at line 1885 of file tmop.hpp.

◆ ~TMOP_Integrator()

mfem::TMOP_Integrator::~TMOP_Integrator ( )

Definition at line 2774 of file tmop.cpp.

Member Function Documentation

◆ ActionIntegrationRule()

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

Definition at line 1819 of file tmop.hpp.

◆ AddMultGradPA()

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.

Definition at line 272 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
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.

Definition at line 250 of file tmop_pa.cpp.

◆ AddMultPA_2D()

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

Definition at line 167 of file tmop_pa_p2.cpp.

◆ AddMultPA_3D()

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

Definition at line 192 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 147 of file tmop_pa_h2d.cpp.

◆ AssembleDiagonalPA_3D()

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

Definition at line 136 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 Tr,
const Vector elfun,
DenseMatrix elmat 
)
virtual

Assemble the local gradient matrix.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 3229 of file tmop.cpp.

◆ AssembleElementGradExact()

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

Definition at line 3383 of file tmop.cpp.

◆ AssembleElementGradFD()

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

Definition at line 3806 of file tmop.cpp.

◆ AssembleElementVector()

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

◆ AssembleElementVectorExact()

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

Definition at line 3244 of file tmop.cpp.

◆ AssembleElementVectorFD()

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

Definition at line 3740 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 3527 of file tmop.cpp.

◆ AssembleElemGradSurfFit()

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

Definition at line 3642 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 3491 of file tmop.cpp.

◆ AssembleElemVecSurfFit()

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

Definition at line 3588 of file tmop.cpp.

◆ AssembleGradDiagonalPA()

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.

Definition at line 225 of file tmop_pa.cpp.

◆ AssembleGradPA()

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.

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 271 of file tmop_pa_h2s.cpp.

◆ AssembleGradPA_3D()

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

Definition at line 337 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/3]

void mfem::NonlinearFormIntegrator::AssemblePA

Method defining partial assembly.

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

Definition at line 25 of file nonlininteg.cpp.

◆ AssemblePA() [2/3]

void mfem::NonlinearFormIntegrator::AssemblePA

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.

Definition at line 31 of file nonlininteg.cpp.

◆ AssemblePA() [3/3]

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.

Definition at line 179 of file tmop_pa.cpp.

◆ AssemblePA_Limiting()

void mfem::TMOP_Integrator::AssemblePA_Limiting ( )
protected

Definition at line 51 of file tmop_pa.cpp.

◆ ComputeAllElementTargets()

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

Definition at line 167 of file tmop_pa.cpp.

◆ ComputeFDh()

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

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

Definition at line 4074 of file tmop.cpp.

◆ ComputeMinDetT()

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

Definition at line 4134 of file tmop.cpp.

◆ ComputeMinJac()

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

Definition at line 4020 of file tmop.cpp.

◆ ComputeNormalizationEnergies()

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

Definition at line 3949 of file tmop.cpp.

◆ ComputeUntangleMetricQuantiles()

void mfem::TMOP_Integrator::ComputeUntangleMetricQuantiles ( const Vector x,
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 4238 of file tmop.cpp.

◆ ComputeUntanglerMaxMuBarrier()

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

Definition at line 4176 of file tmop.cpp.

◆ DisableLimiting()

void mfem::TMOP_Integrator::DisableLimiting ( )
inlineprotected

Definition at line 1804 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 2811 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 2828 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 4090 of file tmop.cpp.

◆ EnableFiniteDifferences() [2/2]

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

Definition at line 4112 of file tmop.cpp.

◆ EnableLimiting() [1/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 2786 of file tmop.cpp.

◆ EnableLimiting() [2/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 2800 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 3926 of file tmop.cpp.

◆ EnableSurfaceFitting() [1/2]

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

◆ EnableSurfaceFitting() [2/2]

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

Parallel support for surface fitting.

Definition at line 2864 of file tmop.cpp.

◆ EnergyIntegrationRule()

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

Definition at line 1810 of file tmop.hpp.

◆ GetAMRQualityMetric()

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

Definition at line 2010 of file tmop.hpp.

◆ GetDerefinementElementEnergy()

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

◆ GetDiscreteAdaptTC()

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

Definition at line 2031 of file tmop.hpp.

◆ GetElementEnergy()

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

◆ GetFDDerivative()

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

◆ GetFDFlag()

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

Definition at line 2047 of file tmop.hpp.

◆ GetFDh()

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

Definition at line 2048 of file tmop.hpp.

◆ GetLocalStateEnergyPA()

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.

Definition at line 297 of file tmop_pa.cpp.

◆ GetLocalStateEnergyPA_2D()

double mfem::TMOP_Integrator::GetLocalStateEnergyPA_2D ( const Vector X) const
protected

Definition at line 145 of file tmop_pa_w2.cpp.

◆ GetLocalStateEnergyPA_3D()

double mfem::TMOP_Integrator::GetLocalStateEnergyPA_3D ( const Vector X) const
protected

Definition at line 155 of file tmop_pa_w3.cpp.

◆ GetLocalStateEnergyPA_C0_2D()

double 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()

double mfem::TMOP_Integrator::GetLocalStateEnergyPA_C0_3D ( const Vector X) const
protected

Definition at line 143 of file tmop_pa_w3_c0.cpp.

◆ GetRefinementElementEnergy()

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

◆ GetSurfaceFittingErrors()

void mfem::TMOP_Integrator::GetSurfaceFittingErrors ( double &  err_avg,
double &  err_max 
)

Definition at line 2882 of file tmop.cpp.

◆ GetSurfaceFittingWeight()

double mfem::TMOP_Integrator::GetSurfaceFittingWeight ( )

Get the surface fitting weight.

Definition at line 3915 of file tmop.cpp.

◆ GradientIntegrationRule()

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

Definition at line 1824 of file tmop.hpp.

◆ IsSurfaceFittingEnabled()

bool mfem::TMOP_Integrator::IsSurfaceFittingEnabled ( )
inline

Definition at line 1973 of file tmop.hpp.

◆ ParEnableNormalization()

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

Definition at line 3936 of file tmop.cpp.

◆ ParUpdateAfterMeshTopologyChange()

void mfem::TMOP_Integrator::ParUpdateAfterMeshTopologyChange ( )

Definition at line 2925 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 2760 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 1908 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 2051 of file tmop.hpp.

◆ SetFDhScale()

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

Definition at line 2046 of file tmop.hpp.

◆ 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 1896 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 1976 of file tmop.hpp.

◆ UpdateAfterMeshPositionChange()

void mfem::TMOP_Integrator::UpdateAfterMeshPositionChange ( const Vector new_x,
int  new_x_ordering = Ordering::byNODES 
)
protected

Definition at line 4054 of file tmop.cpp.

◆ UpdateAfterMeshTopologyChange()

void mfem::TMOP_Integrator::UpdateAfterMeshTopologyChange ( )

Definition at line 2912 of file tmop.cpp.

◆ UpdateSurfaceFittingWeight()

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

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

Definition at line 3903 of file tmop.cpp.

Friends And Related Function Documentation

◆ TMOPComboIntegrator

friend class TMOPComboIntegrator
friend

Definition at line 1642 of file tmop.hpp.

◆ TMOPNewtonSolver

friend class TMOPNewtonSolver
friend

Definition at line 1641 of file tmop.hpp.

Member Data Documentation

◆ adapt_lim_coeff

Coefficient* mfem::TMOP_Integrator::adapt_lim_coeff
protected

Definition at line 1675 of file tmop.hpp.

◆ adapt_lim_eval

AdaptivityEvaluator* mfem::TMOP_Integrator::adapt_lim_eval
protected

Definition at line 1676 of file tmop.hpp.

◆ adapt_lim_gf

GridFunction* mfem::TMOP_Integrator::adapt_lim_gf
protected

Definition at line 1674 of file tmop.hpp.

◆ adapt_lim_gf0

const GridFunction* mfem::TMOP_Integrator::adapt_lim_gf0
protected

Definition at line 1670 of file tmop.hpp.

◆ adapt_lim_pgf0

const ParGridFunction* mfem::TMOP_Integrator::adapt_lim_pgf0
protected

Definition at line 1672 of file tmop.hpp.

◆ C0

Vector mfem::TMOP_Integrator::C0
mutable

Definition at line 1744 of file tmop.hpp.

◆ dim

int mfem::TMOP_Integrator::dim

Definition at line 1740 of file tmop.hpp.

◆ discr_tc

DiscreteAdaptTC* mfem::TMOP_Integrator::discr_tc
protected

Definition at line 1685 of file tmop.hpp.

◆ DS

DenseMatrix mfem::TMOP_Integrator::DS
protected

Definition at line 1711 of file tmop.hpp.

◆ DSh

DenseMatrix mfem::TMOP_Integrator::DSh
protected

Definition at line 1711 of file tmop.hpp.

◆ dx

double mfem::TMOP_Integrator::dx
protected

Definition at line 1689 of file tmop.hpp.

◆ dxscale

double mfem::TMOP_Integrator::dxscale
protected

Definition at line 1690 of file tmop.hpp.

◆ E

Vector mfem::TMOP_Integrator::E
mutable

Definition at line 1744 of file tmop.hpp.

◆ ElemDer

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

Definition at line 1698 of file tmop.hpp.

◆ ElemPertEnergy

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

Definition at line 1699 of file tmop.hpp.

◆ enabled

bool mfem::TMOP_Integrator::enabled

Definition at line 1739 of file tmop.hpp.

◆ exact_action

bool mfem::TMOP_Integrator::exact_action
protected

Definition at line 1696 of file tmop.hpp.

◆ fd_call_flag

bool mfem::TMOP_Integrator::fd_call_flag
protected

Definition at line 1693 of file tmop.hpp.

◆ fdflag

bool mfem::TMOP_Integrator::fdflag
protected

Definition at line 1688 of file tmop.hpp.

◆ fes

const FiniteElementSpace* mfem::TMOP_Integrator::fes

Definition at line 1748 of file tmop.hpp.

◆ geom

const GeometricFactors* mfem::TMOP_Integrator::geom

Definition at line 1747 of file tmop.hpp.

◆ H

Vector mfem::TMOP_Integrator::H
mutable

Definition at line 1744 of file tmop.hpp.

◆ H0

Vector mfem::TMOP_Integrator::H0
mutable

Definition at line 1744 of file tmop.hpp.

◆ h_metric

TMOP_QualityMetric* mfem::TMOP_Integrator::h_metric
protected

Definition at line 1644 of file tmop.hpp.

◆ integ_order

int mfem::TMOP_Integrator::integ_order
protected

Definition at line 1650 of file tmop.hpp.

◆ IntegRules

IntegrationRules* mfem::TMOP_Integrator::IntegRules
protected

Definition at line 1649 of file tmop.hpp.

◆ ir

const IntegrationRule* mfem::TMOP_Integrator::ir

Definition at line 1749 of file tmop.hpp.

◆ Jpr

DenseMatrix mfem::TMOP_Integrator::Jpr
protected

Definition at line 1711 of file tmop.hpp.

◆ Jpt

DenseMatrix mfem::TMOP_Integrator::Jpt
protected

Definition at line 1711 of file tmop.hpp.

◆ Jrt

DenseMatrix mfem::TMOP_Integrator::Jrt
protected

Definition at line 1711 of file tmop.hpp.

◆ Jtr

DenseTensor mfem::TMOP_Integrator::Jtr
mutable

Definition at line 1741 of file tmop.hpp.

◆ Jtr_debug_grad

bool mfem::TMOP_Integrator::Jtr_debug_grad
mutable

Definition at line 1743 of file tmop.hpp.

◆ Jtr_needs_update

bool mfem::TMOP_Integrator::Jtr_needs_update
mutable

Definition at line 1742 of file tmop.hpp.

◆ LD

Vector mfem::TMOP_Integrator::LD
mutable

Definition at line 1744 of file tmop.hpp.

◆ lim_coeff

Coefficient* mfem::TMOP_Integrator::lim_coeff
protected

Definition at line 1661 of file tmop.hpp.

◆ lim_dist

const GridFunction* mfem::TMOP_Integrator::lim_dist
protected

Definition at line 1663 of file tmop.hpp.

◆ lim_func

TMOP_LimiterFunction* mfem::TMOP_Integrator::lim_func
protected

Definition at line 1665 of file tmop.hpp.

◆ lim_nodes0

const GridFunction* mfem::TMOP_Integrator::lim_nodes0
protected

Definition at line 1660 of file tmop.hpp.

◆ lim_normal

double mfem::TMOP_Integrator::lim_normal
protected

Definition at line 1667 of file tmop.hpp.

◆ maps

const DofToQuad* mfem::TMOP_Integrator::maps

Definition at line 1745 of file tmop.hpp.

◆ maps_lim

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

Definition at line 1746 of file tmop.hpp.

◆ metric

TMOP_QualityMetric* mfem::TMOP_Integrator::metric
protected

Definition at line 1645 of file tmop.hpp.

◆ metric_coeff

Coefficient* mfem::TMOP_Integrator::metric_coeff
protected

Definition at line 1653 of file tmop.hpp.

◆ metric_normal

double mfem::TMOP_Integrator::metric_normal
protected

Definition at line 1655 of file tmop.hpp.

◆ ne

int mfem::TMOP_Integrator::ne

Definition at line 1740 of file tmop.hpp.

◆ nq

int mfem::TMOP_Integrator::nq

Definition at line 1740 of file tmop.hpp.

◆ O

Vector mfem::TMOP_Integrator::O
mutable

Definition at line 1744 of file tmop.hpp.

◆ P

DenseMatrix mfem::TMOP_Integrator::P
protected

Definition at line 1711 of file tmop.hpp.

◆ PA

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

◆ PMatI

DenseMatrix mfem::TMOP_Integrator::PMatI
protected

Definition at line 1711 of file tmop.hpp.

◆ PMatO

DenseMatrix mfem::TMOP_Integrator::PMatO
protected

Definition at line 1711 of file tmop.hpp.

◆ surf_fit_coeff

Coefficient* mfem::TMOP_Integrator::surf_fit_coeff
protected

Definition at line 1681 of file tmop.hpp.

◆ surf_fit_eval

AdaptivityEvaluator* mfem::TMOP_Integrator::surf_fit_eval
protected

Definition at line 1682 of file tmop.hpp.

◆ surf_fit_gf

GridFunction* mfem::TMOP_Integrator::surf_fit_gf
protected

Definition at line 1679 of file tmop.hpp.

◆ surf_fit_marker

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

Definition at line 1680 of file tmop.hpp.

◆ surf_fit_normal

double mfem::TMOP_Integrator::surf_fit_normal
protected

Definition at line 1683 of file tmop.hpp.

◆ targetC

const TargetConstructor* mfem::TMOP_Integrator::targetC
protected

Definition at line 1646 of file tmop.hpp.

◆ X0

Vector mfem::TMOP_Integrator::X0
mutable

Definition at line 1744 of file tmop.hpp.


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