MFEM  v4.2.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 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...
 
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
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 &irule)
 Prescribe a fixed IntegrationRule to use. 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 &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
virtual void AddMultPA (const Vector &x, Vector &y) const
 Method for partially assembled action. More...
 
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
 
- 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
 
- Protected Attributes inherited from mfem::NonlinearFormIntegrator
const IntegrationRuleIntRule
 

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

mfem::TMOP_Integrator::~TMOP_Integrator ( )

Definition at line 1885 of file tmop.cpp.

Member Function Documentation

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

Definition at line 1004 of file tmop.hpp.

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

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

Definition at line 2229 of file tmop.cpp.

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

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

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

Definition at line 2091 of file tmop.cpp.

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

Definition at line 2460 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 2374 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 2335 of file tmop.cpp.

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

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

Definition at line 2729 of file tmop.cpp.

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

Definition at line 2736 of file tmop.cpp.

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

Definition at line 2690 of file tmop.cpp.

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

Definition at line 2638 of file tmop.cpp.

void mfem::TMOP_Integrator::DisableLimiting ( )
inlineprotected

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

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

Parallel support for adaptive limiting.

Definition at line 1938 of file tmop.cpp.

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

Enables FD-based approximation and computes dx.

Definition at line 2747 of file tmop.cpp.

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

Definition at line 2761 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 1896 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 1903 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 2619 of file tmop.cpp.

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

Definition at line 995 of file tmop.hpp.

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

Definition at line 1103 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 1955 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 2438 of file tmop.cpp.

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

Definition at line 1119 of file tmop.hpp.

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

Definition at line 1120 of file tmop.hpp.

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

Definition at line 1009 of file tmop.hpp.

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

Definition at line 2627 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 1044 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 1123 of file tmop.hpp.

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

Definition at line 1118 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 1032 of file tmop.hpp.

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

Update the original/reference nodes used for limiting.

Definition at line 1085 of file tmop.hpp.

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

Definition at line 2723 of file tmop.cpp.

Friends And Related Function Documentation

friend class TMOPComboIntegrator
friend

Definition at line 888 of file tmop.hpp.

friend class TMOPNewtonSolver
friend

Definition at line 887 of file tmop.hpp.

Member Data Documentation

AdaptivityEvaluator* mfem::TMOP_Integrator::adapt_eval
protected

Definition at line 918 of file tmop.hpp.

Coefficient* mfem::TMOP_Integrator::coeff0
protected

Definition at line 906 of file tmop.hpp.

Coefficient* mfem::TMOP_Integrator::coeff1
protected

Definition at line 898 of file tmop.hpp.

Coefficient* mfem::TMOP_Integrator::coeff_zeta
protected

Definition at line 917 of file tmop.hpp.

DiscreteAdaptTC* mfem::TMOP_Integrator::discr_tc
protected

Definition at line 920 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::DS
protected

Definition at line 946 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::DSh
protected

Definition at line 946 of file tmop.hpp.

double mfem::TMOP_Integrator::dx
protected

Definition at line 924 of file tmop.hpp.

double mfem::TMOP_Integrator::dxscale
protected

Definition at line 925 of file tmop.hpp.

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

Definition at line 933 of file tmop.hpp.

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

Definition at line 934 of file tmop.hpp.

bool mfem::TMOP_Integrator::exact_action
protected

Definition at line 931 of file tmop.hpp.

bool mfem::TMOP_Integrator::fd_call_flag
protected

Definition at line 928 of file tmop.hpp.

bool mfem::TMOP_Integrator::fdflag
protected

Definition at line 923 of file tmop.hpp.

int mfem::TMOP_Integrator::integ_order
protected

Definition at line 895 of file tmop.hpp.

IntegrationRules* mfem::TMOP_Integrator::IntegRules
protected

Definition at line 894 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::Jpr
protected

Definition at line 946 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::Jpt
protected

Definition at line 946 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::Jrt
protected

Definition at line 946 of file tmop.hpp.

const GridFunction* mfem::TMOP_Integrator::lim_dist
protected

Definition at line 908 of file tmop.hpp.

TMOP_LimiterFunction* mfem::TMOP_Integrator::lim_func
protected

Definition at line 910 of file tmop.hpp.

double mfem::TMOP_Integrator::lim_normal
protected

Definition at line 912 of file tmop.hpp.

TMOP_QualityMetric* mfem::TMOP_Integrator::metric
protected

Definition at line 890 of file tmop.hpp.

double mfem::TMOP_Integrator::metric_normal
protected

Definition at line 900 of file tmop.hpp.

const GridFunction* mfem::TMOP_Integrator::nodes0
protected

Definition at line 905 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::P
protected

Definition at line 946 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::PMatI
protected

Definition at line 946 of file tmop.hpp.

DenseMatrix mfem::TMOP_Integrator::PMatO
protected

Definition at line 946 of file tmop.hpp.

const TargetConstructor* mfem::TMOP_Integrator::targetC
protected

Definition at line 891 of file tmop.hpp.

GridFunction* mfem::TMOP_Integrator::zeta
protected

Definition at line 916 of file tmop.hpp.

const GridFunction* mfem::TMOP_Integrator::zeta_0
protected

Definition at line 915 of file tmop.hpp.


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