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

Base class representing target-matrix construction algorithms for mesh optimization via the target-matrix optimization paradigm (TMOP). More...

#include <tmop.hpp>

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

Public Types

enum  TargetType {
  IDEAL_SHAPE_UNIT_SIZE, IDEAL_SHAPE_EQUAL_SIZE, IDEAL_SHAPE_GIVEN_SIZE, GIVEN_SHAPE_AND_SIZE,
  GIVEN_FULL
}
 Target-matrix construction algorithms supported by this class. More...
 

Public Member Functions

 TargetConstructor (TargetType ttype)
 Constructor for use in serial. More...
 
 TargetConstructor (TargetType ttype, MPI_Comm mpicomm)
 Constructor for use in parallel. More...
 
virtual ~TargetConstructor ()
 
void SetNodes (const GridFunction &n)
 Set the nodes to be used in the target-matrix construction. More...
 
const GridFunctionGetNodes () const
 Get the nodes to be used in the target-matrix construction. More...
 
void SetVolumeScale (double vol_scale)
 Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1. More...
 
TargetType GetTargetType () const
 
bool UsesPhysicalCoordinates () const
 Return true if the methods ComputeElementTargets(), ComputeAllElementTargets(), and ComputeElementTargetsGradient() use the physical node coordinates provided by the parameters 'elfun', or 'xe'. More...
 
virtual bool ContainsVolumeInfo () const
 Checks if the target matrices contain non-trivial size specification. More...
 
virtual void ComputeElementTargets (int e_id, const FiniteElement &fe, const IntegrationRule &ir, const Vector &elfun, DenseTensor &Jtr) const
 Given an element and quadrature rule, computes ref->target transformation Jacobians for each quadrature point in the element. The physical positions of the element's nodes are given by elfun. More...
 
virtual void ComputeAllElementTargets (const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
 Computes reference-to-target transformation Jacobians for all quadrature points in all elements. More...
 
virtual void ComputeElementTargetsGradient (const IntegrationRule &ir, const Vector &elfun, IsoparametricTransformation &Tpr, DenseTensor &dJtr) const
 

Protected Member Functions

bool Parallel () const
 
bool Parallel () const
 
void ComputeAvgVolume () const
 
template<int DIM>
bool ComputeAllElementTargets (const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
 
void ComputeAllElementTargets_Fallback (const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
 

Protected Attributes

const GridFunctionnodes
 
double avg_volume
 
double volume_scale
 
const TargetType target_type
 
bool uses_phys_coords
 
MPI_Comm comm
 

Detailed Description

Base class representing target-matrix construction algorithms for mesh optimization via the target-matrix optimization paradigm (TMOP).

This class is used by class TMOP_Integrator to construct the target Jacobian matrices (reference-element to target-element) at quadrature points. It supports a set of algorithms chosen by the TargetType enumeration.

New target-matrix construction algorithms can be defined by deriving new classes and overriding the methods ComputeElementTargets() and ContainsVolumeInfo().

Definition at line 868 of file tmop.hpp.

Member Enumeration Documentation

Target-matrix construction algorithms supported by this class.

Enumerator
IDEAL_SHAPE_UNIT_SIZE 

Ideal shape, unit size; the nodes are not used.

IDEAL_SHAPE_EQUAL_SIZE 

Ideal shape, equal size/volume; the given nodes define the total target volume; for each mesh element, the target volume is the average volume multiplied by the volume scale, set with SetVolumeScale().

IDEAL_SHAPE_GIVEN_SIZE 

Ideal shape, given size/volume; the given nodes define the target volume at all quadrature points.

GIVEN_SHAPE_AND_SIZE 

Given shape, given size/volume; the given nodes define the exact target Jacobian matrix at all quadrature points.

GIVEN_FULL 

Full target tensor is specified at every quadrature point.

Definition at line 872 of file tmop.hpp.

Constructor & Destructor Documentation

mfem::TargetConstructor::TargetConstructor ( TargetType  ttype)
inline

Constructor for use in serial.

Definition at line 923 of file tmop.hpp.

mfem::TargetConstructor::TargetConstructor ( TargetType  ttype,
MPI_Comm  mpicomm 
)
inline

Constructor for use in parallel.

Definition at line 933 of file tmop.hpp.

virtual mfem::TargetConstructor::~TargetConstructor ( )
inlinevirtual

Definition at line 937 of file tmop.hpp.

Member Function Documentation

template<int DIM>
bool mfem::TargetConstructor::ComputeAllElementTargets ( const FiniteElementSpace fes,
const IntegrationRule ir,
const Vector xe,
DenseTensor Jtr 
) const
protected
virtual void mfem::TargetConstructor::ComputeAllElementTargets ( const FiniteElementSpace fes,
const IntegrationRule ir,
const Vector xe,
DenseTensor Jtr 
) const
virtual

Computes reference-to-target transformation Jacobians for all quadrature points in all elements.

Parameters
[in]fesThe nodal FE space
[in]irThe quadrature rule to use for all elements
[in]xeE-vector with the current physical coordinates/positions; this parameter is used only when needed by the target constructor, see UsesPhysicalCoordinates()
[out]JtrThe computed ref->target Jacobian matrices.

Reimplemented in mfem::DiscreteAdaptTC, and mfem::AnalyticAdaptTC.

void mfem::TargetConstructor::ComputeAllElementTargets_Fallback ( const FiniteElementSpace fes,
const IntegrationRule ir,
const Vector xe,
DenseTensor Jtr 
) const
protected

Definition at line 1056 of file tmop.cpp.

void mfem::TargetConstructor::ComputeAvgVolume ( ) const
protected

Definition at line 1014 of file tmop.cpp.

void mfem::TargetConstructor::ComputeElementTargets ( int  e_id,
const FiniteElement fe,
const IntegrationRule ir,
const Vector elfun,
DenseTensor Jtr 
) const
virtual

Given an element and quadrature rule, computes ref->target transformation Jacobians for each quadrature point in the element. The physical positions of the element's nodes are given by elfun.

Reimplemented in mfem::DiscreteAdaptTC, and mfem::AnalyticAdaptTC.

Definition at line 1142 of file tmop.cpp.

void mfem::TargetConstructor::ComputeElementTargetsGradient ( const IntegrationRule ir,
const Vector elfun,
IsoparametricTransformation Tpr,
DenseTensor dJtr 
) const
virtual

Reimplemented in mfem::DiscreteAdaptTC, and mfem::AnalyticAdaptTC.

Definition at line 1213 of file tmop.cpp.

bool mfem::TargetConstructor::ContainsVolumeInfo ( ) const
virtual

Checks if the target matrices contain non-trivial size specification.

Definition at line 1128 of file tmop.cpp.

const GridFunction* mfem::TargetConstructor::GetNodes ( ) const
inline

Get the nodes to be used in the target-matrix construction.

Definition at line 947 of file tmop.hpp.

TargetType mfem::TargetConstructor::GetTargetType ( ) const
inline

Definition at line 952 of file tmop.hpp.

bool mfem::TargetConstructor::Parallel ( ) const
inlineprotected

Definition at line 900 of file tmop.hpp.

bool mfem::TargetConstructor::Parallel ( ) const
inlineprotected

Definition at line 902 of file tmop.hpp.

void mfem::TargetConstructor::SetNodes ( const GridFunction n)
inline

Set the nodes to be used in the target-matrix construction.

This method should be called every time the target nodes are updated externally and recomputation of the target average volume is needed. The nodes are used by all target types except IDEAL_SHAPE_UNIT_SIZE.

Definition at line 944 of file tmop.hpp.

void mfem::TargetConstructor::SetVolumeScale ( double  vol_scale)
inline

Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.

Definition at line 950 of file tmop.hpp.

bool mfem::TargetConstructor::UsesPhysicalCoordinates ( ) const
inline

Return true if the methods ComputeElementTargets(), ComputeAllElementTargets(), and ComputeElementTargetsGradient() use the physical node coordinates provided by the parameters 'elfun', or 'xe'.

Definition at line 957 of file tmop.hpp.

Member Data Documentation

double mfem::TargetConstructor::avg_volume
mutableprotected

Definition at line 893 of file tmop.hpp.

MPI_Comm mfem::TargetConstructor::comm
protected

Definition at line 899 of file tmop.hpp.

const GridFunction* mfem::TargetConstructor::nodes
protected

Definition at line 892 of file tmop.hpp.

const TargetType mfem::TargetConstructor::target_type
protected

Definition at line 895 of file tmop.hpp.

bool mfem::TargetConstructor::uses_phys_coords
protected

Definition at line 896 of file tmop.hpp.

double mfem::TargetConstructor::volume_scale
protected

Definition at line 894 of file tmop.hpp.


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