MFEM  v4.6.0
Finite element discretization library
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
mfem::DiscreteAdaptTC Class Reference

#include <tmop.hpp>

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

Public Member Functions

 DiscreteAdaptTC (TargetType ttype)
 
virtual ~DiscreteAdaptTC ()
 
void ResetUpdateFlags ()
 Used in combination with the Update methods to avoid extra computations. More...
 
void GetDiscreteTargetSpec (GridFunction &tspec_, int idx)
 Get one of the discrete fields from tspec. More...
 
FiniteElementSpaceGetTSpecFESpace ()
 Get the FESpace associated with tspec. More...
 
GridFunctionGetTSpecData ()
 Get the entire tspec. More...
 
void UpdateAfterMeshTopologyChange ()
 Update all discrete fields based on tspec and update for AMR. More...
 
ParFiniteElementSpaceGetTSpecParFESpace ()
 
void ParUpdateAfterMeshTopologyChange ()
 
void UpdateTargetSpecification (const Vector &new_x, bool reuse_flag=false, int new_x_ordering=Ordering::byNODES)
 
void UpdateTargetSpecification (Vector &new_x, Vector &IntData, int new_x_ordering=Ordering::byNODES)
 
void UpdateTargetSpecificationAtNode (const FiniteElement &el, ElementTransformation &T, int nodenum, int idir, const Vector &IntData)
 
void RestoreTargetSpecificationAtNode (ElementTransformation &T, int nodenum)
 
void UpdateGradientTargetSpecification (const Vector &x, double dx, bool reuse_flag=false, int x_ordering=Ordering::byNODES)
 
void UpdateHessianTargetSpecification (const Vector &x, double dx, bool reuse_flag=false, int x_ordering=Ordering::byNODES)
 
void SetAdaptivityEvaluator (AdaptivityEvaluator *ae)
 
const AdaptivityEvaluatorGetAdaptivityEvaluator () const
 
const VectorGetTspecPert1H ()
 
const VectorGetTspecPert2H ()
 
const VectorGetTspecPertMixH ()
 
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. Note that this function assumes that UpdateTargetSpecification() has been called with the position vector corresponding to 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
 
void SetTspecFromIntRule (int e_id, const IntegrationRule &intrule)
 
void SetMinSizeForTargets (double min_size_)
 
void SetTspecDataForDerefinement (FiniteElementSpace *fes)
 Computes target specification data with respect to the coarse FE space. More...
 
void ResetRefinementTspecData ()
 
void ResetDerefinementTspecData ()
 
void SetRefinementSubElement (int amr_el_)
 
Target specification methods.

The following methods are used to specify geometric parameters of the targets when these parameters are given by discrete FE functions. Note that every GridFunction given to the Set methods must use a H1_FECollection of the same order. The number of components must correspond to the type of geometric parameter and dimension.

Parameters
[in]tspec_Input values of a geometric parameter. Note that the methods in this class support only functions that use H1_FECollection collection of the same order.
virtual void SetSerialDiscreteTargetSpec (const GridFunction &tspec_)
 
virtual void SetSerialDiscreteTargetSize (const GridFunction &tspec_)
 
virtual void SetSerialDiscreteTargetSkew (const GridFunction &tspec_)
 
virtual void SetSerialDiscreteTargetAspectRatio (const GridFunction &tspec_)
 
virtual void SetSerialDiscreteTargetOrientation (const GridFunction &tspec_)
 
virtual void SetParDiscreteTargetSpec (const ParGridFunction &tspec_)
 
virtual void SetParDiscreteTargetSize (const ParGridFunction &tspec_)
 
virtual void SetParDiscreteTargetSkew (const ParGridFunction &tspec_)
 
virtual void SetParDiscreteTargetAspectRatio (const ParGridFunction &tspec_)
 
virtual void SetParDiscreteTargetOrientation (const ParGridFunction &tspec_)
 
- Public Member Functions inherited from mfem::TargetConstructor
 TargetConstructor (TargetType ttype)
 Constructor for use in serial. More...
 
 TargetConstructor (TargetType ttype, MPI_Comm mpicomm)
 Constructor for use in parallel. More...
 
virtual ~TargetConstructor ()
 
bool Parallel () const
 
MPI_Comm GetComm () const
 
bool Parallel () const
 
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...
 

Protected Member Functions

void SetDiscreteTargetBase (const GridFunction &tspec_)
 
void SetTspecAtIndex (int idx, const GridFunction &tspec_)
 
void FinalizeSerialDiscreteTargetSpec (const GridFunction &tspec_)
 
void SetTspecAtIndex (int idx, const ParGridFunction &tspec_)
 
void FinalizeParDiscreteTargetSpec (const ParGridFunction &tspec_)
 
- Protected Member Functions inherited from mfem::TargetConstructor
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
 
template<>
bool ComputeAllElementTargets (const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &, DenseTensor &Jtr) const
 
template<>
bool ComputeAllElementTargets (const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &, DenseTensor &Jtr) const
 

Protected Attributes

int ncomp
 
int sizeidx
 
int skewidx
 
int aspectratioidx
 
int orientationidx
 
Vector tspec
 
Vector tspec_sav
 
Vector tspec_pert1h
 
Vector tspec_pert2h
 
Vector tspec_pertmix
 
DenseMatrix tspec_refine
 
Vector tspec_derefine
 
DenseTensor Jtrcomp
 
FiniteElementSpacetspec_fesv
 
FiniteElementSpacecoarse_tspec_fesv
 
GridFunctiontspec_gf
 
ParFiniteElementSpaceptspec_fesv
 
ParGridFunctiontspec_pgf
 
int amr_el
 
double lim_min_size
 
bool good_tspec
 
bool good_tspec_grad
 
bool good_tspec_hess
 
AdaptivityEvaluatoradapt_eval
 
- Protected Attributes inherited from mfem::TargetConstructor
const GridFunctionnodes
 
double avg_volume
 
double volume_scale
 
const TargetType target_type
 
bool uses_phys_coords
 
MPI_Comm comm
 

Additional Inherited Members

- Public Types inherited from mfem::TargetConstructor
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...
 

Detailed Description

Definition at line 1510 of file tmop.hpp.

Constructor & Destructor Documentation

◆ DiscreteAdaptTC()

mfem::DiscreteAdaptTC::DiscreteAdaptTC ( TargetType  ttype)
inline

Definition at line 1568 of file tmop.hpp.

◆ ~DiscreteAdaptTC()

mfem::DiscreteAdaptTC::~DiscreteAdaptTC ( )
virtual

Definition at line 2771 of file tmop.cpp.

Member Function Documentation

◆ ComputeAllElementTargets()

void mfem::DiscreteAdaptTC::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 from mfem::TargetConstructor.

Definition at line 130 of file tmop_pa_da3.cpp.

◆ ComputeElementTargets()

void mfem::DiscreteAdaptTC::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. Note that this function assumes that UpdateTargetSpecification() has been called with the position vector corresponding to elfun.

Reimplemented from mfem::TargetConstructor.

Definition at line 2107 of file tmop.cpp.

◆ ComputeElementTargetsGradient()

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

Reimplemented from mfem::TargetConstructor.

Definition at line 2322 of file tmop.cpp.

◆ FinalizeParDiscreteTargetSpec()

void mfem::DiscreteAdaptTC::FinalizeParDiscreteTargetSpec ( const ParGridFunction tspec_)
protected

Definition at line 1791 of file tmop.cpp.

◆ FinalizeSerialDiscreteTargetSpec()

void mfem::DiscreteAdaptTC::FinalizeSerialDiscreteTargetSpec ( const GridFunction tspec_)
protected

Definition at line 1970 of file tmop.cpp.

◆ GetAdaptivityEvaluator()

const AdaptivityEvaluator* mfem::DiscreteAdaptTC::GetAdaptivityEvaluator ( ) const
inline

Definition at line 1664 of file tmop.hpp.

◆ GetDiscreteTargetSpec()

void mfem::DiscreteAdaptTC::GetDiscreteTargetSpec ( GridFunction tspec_,
int  idx 
)

Get one of the discrete fields from tspec.

Definition at line 1991 of file tmop.cpp.

◆ GetTSpecData()

GridFunction* mfem::DiscreteAdaptTC::GetTSpecData ( )
inline

Get the entire tspec.

Definition at line 1618 of file tmop.hpp.

◆ GetTSpecFESpace()

FiniteElementSpace* mfem::DiscreteAdaptTC::GetTSpecFESpace ( )
inline

Get the FESpace associated with tspec.

Definition at line 1616 of file tmop.hpp.

◆ GetTSpecParFESpace()

ParFiniteElementSpace* mfem::DiscreteAdaptTC::GetTSpecParFESpace ( )
inline

Definition at line 1623 of file tmop.hpp.

◆ GetTspecPert1H()

const Vector& mfem::DiscreteAdaptTC::GetTspecPert1H ( )
inline

Definition at line 1669 of file tmop.hpp.

◆ GetTspecPert2H()

const Vector& mfem::DiscreteAdaptTC::GetTspecPert2H ( )
inline

Definition at line 1670 of file tmop.hpp.

◆ GetTspecPertMixH()

const Vector& mfem::DiscreteAdaptTC::GetTspecPertMixH ( )
inline

Definition at line 1671 of file tmop.hpp.

◆ ParUpdateAfterMeshTopologyChange()

void mfem::DiscreteAdaptTC::ParUpdateAfterMeshTopologyChange ( )

Definition at line 1816 of file tmop.cpp.

◆ ResetDerefinementTspecData()

void mfem::DiscreteAdaptTC::ResetDerefinementTspecData ( )
inline

Definition at line 1713 of file tmop.hpp.

◆ ResetRefinementTspecData()

void mfem::DiscreteAdaptTC::ResetRefinementTspecData ( )
inline

Definition at line 1706 of file tmop.hpp.

◆ ResetUpdateFlags()

void mfem::DiscreteAdaptTC::ResetUpdateFlags ( )
inline

Used in combination with the Update methods to avoid extra computations.

Definition at line 1610 of file tmop.hpp.

◆ RestoreTargetSpecificationAtNode()

void mfem::DiscreteAdaptTC::RestoreTargetSpecificationAtNode ( ElementTransformation T,
int  nodenum 
)

Definition at line 2059 of file tmop.cpp.

◆ SetAdaptivityEvaluator()

void mfem::DiscreteAdaptTC::SetAdaptivityEvaluator ( AdaptivityEvaluator ae)
inline

Definition at line 1658 of file tmop.hpp.

◆ SetDiscreteTargetBase()

void mfem::DiscreteAdaptTC::SetDiscreteTargetBase ( const GridFunction tspec_)
protected

Definition at line 1893 of file tmop.cpp.

◆ SetMinSizeForTargets()

void mfem::DiscreteAdaptTC::SetMinSizeForTargets ( double  min_size_)
inline

Definition at line 1700 of file tmop.hpp.

◆ SetParDiscreteTargetAspectRatio()

void mfem::DiscreteAdaptTC::SetParDiscreteTargetAspectRatio ( const ParGridFunction tspec_)
virtual

Definition at line 1867 of file tmop.cpp.

◆ SetParDiscreteTargetOrientation()

void mfem::DiscreteAdaptTC::SetParDiscreteTargetOrientation ( const ParGridFunction tspec_)
virtual

Definition at line 1877 of file tmop.cpp.

◆ SetParDiscreteTargetSize()

void mfem::DiscreteAdaptTC::SetParDiscreteTargetSize ( const ParGridFunction tspec_)
virtual

Definition at line 1847 of file tmop.cpp.

◆ SetParDiscreteTargetSkew()

void mfem::DiscreteAdaptTC::SetParDiscreteTargetSkew ( const ParGridFunction tspec_)
virtual

Definition at line 1857 of file tmop.cpp.

◆ SetParDiscreteTargetSpec()

void mfem::DiscreteAdaptTC::SetParDiscreteTargetSpec ( const ParGridFunction tspec_)
virtual

Definition at line 1887 of file tmop.cpp.

◆ SetRefinementSubElement()

void mfem::DiscreteAdaptTC::SetRefinementSubElement ( int  amr_el_)
inline

Definition at line 1721 of file tmop.hpp.

◆ SetSerialDiscreteTargetAspectRatio()

void mfem::DiscreteAdaptTC::SetSerialDiscreteTargetAspectRatio ( const GridFunction tspec_)
virtual

Definition at line 1950 of file tmop.cpp.

◆ SetSerialDiscreteTargetOrientation()

void mfem::DiscreteAdaptTC::SetSerialDiscreteTargetOrientation ( const GridFunction tspec_)
virtual

Definition at line 1960 of file tmop.cpp.

◆ SetSerialDiscreteTargetSize()

void mfem::DiscreteAdaptTC::SetSerialDiscreteTargetSize ( const GridFunction tspec_)
virtual

Definition at line 1930 of file tmop.cpp.

◆ SetSerialDiscreteTargetSkew()

void mfem::DiscreteAdaptTC::SetSerialDiscreteTargetSkew ( const GridFunction tspec_)
virtual

Definition at line 1940 of file tmop.cpp.

◆ SetSerialDiscreteTargetSpec()

void mfem::DiscreteAdaptTC::SetSerialDiscreteTargetSpec ( const GridFunction tspec_)
virtual

Definition at line 2016 of file tmop.cpp.

◆ SetTspecAtIndex() [1/2]

void mfem::DiscreteAdaptTC::SetTspecAtIndex ( int  idx,
const GridFunction tspec_ 
)
protected

Definition at line 1917 of file tmop.cpp.

◆ SetTspecAtIndex() [2/2]

void mfem::DiscreteAdaptTC::SetTspecAtIndex ( int  idx,
const ParGridFunction tspec_ 
)
protected

Definition at line 1834 of file tmop.cpp.

◆ SetTspecDataForDerefinement()

void mfem::DiscreteAdaptTC::SetTspecDataForDerefinement ( FiniteElementSpace fes)

Computes target specification data with respect to the coarse FE space.

Definition at line 2099 of file tmop.cpp.

◆ SetTspecFromIntRule()

void mfem::DiscreteAdaptTC::SetTspecFromIntRule ( int  e_id,
const IntegrationRule intrule 
)

Definition at line 2073 of file tmop.cpp.

◆ UpdateAfterMeshTopologyChange()

void mfem::DiscreteAdaptTC::UpdateAfterMeshTopologyChange ( )

Update all discrete fields based on tspec and update for AMR.

Definition at line 2005 of file tmop.cpp.

◆ UpdateGradientTargetSpecification()

void mfem::DiscreteAdaptTC::UpdateGradientTargetSpecification ( const Vector x,
double  dx,
bool  reuse_flag = false,
int  x_ordering = Ordering::byNODES 
)

Used for finite-difference based computations. Computes the target specifications after a mesh perturbation in x or y direction. If reuse_flag is true, repeated calls won't do anything until ResetUpdateFlags() is called.

Definition at line 2669 of file tmop.cpp.

◆ UpdateHessianTargetSpecification()

void mfem::DiscreteAdaptTC::UpdateHessianTargetSpecification ( const Vector x,
double  dx,
bool  reuse_flag = false,
int  x_ordering = Ordering::byNODES 
)

Used for finite-difference based computations. Computes the target specifications after two mesh perturbations in x and/or y direction. If reuse_flag is true, repeated calls won't do anything until ResetUpdateFlags() is called.

Definition at line 2705 of file tmop.cpp.

◆ UpdateTargetSpecification() [1/2]

void mfem::DiscreteAdaptTC::UpdateTargetSpecification ( const Vector new_x,
bool  reuse_flag = false,
int  new_x_ordering = Ordering::byNODES 
)

Used to update the target specification after the mesh has changed. The new mesh positions are given by new_x. If reuse_flag is true, repeated calls won't do anything until ResetUpdateFlags() is called.

Definition at line 2022 of file tmop.cpp.

◆ UpdateTargetSpecification() [2/2]

void mfem::DiscreteAdaptTC::UpdateTargetSpecification ( Vector new_x,
Vector IntData,
int  new_x_ordering = Ordering::byNODES 
)

Definition at line 2035 of file tmop.cpp.

◆ UpdateTargetSpecificationAtNode()

void mfem::DiscreteAdaptTC::UpdateTargetSpecificationAtNode ( const FiniteElement el,
ElementTransformation T,
int  nodenum,
int  idir,
const Vector IntData 
)

Definition at line 2042 of file tmop.cpp.

Member Data Documentation

◆ adapt_eval

AdaptivityEvaluator* mfem::DiscreteAdaptTC::adapt_eval
protected

Definition at line 1557 of file tmop.hpp.

◆ amr_el

int mfem::DiscreteAdaptTC::amr_el
protected

Definition at line 1548 of file tmop.hpp.

◆ aspectratioidx

int mfem::DiscreteAdaptTC::aspectratioidx
protected

Definition at line 1515 of file tmop.hpp.

◆ coarse_tspec_fesv

FiniteElementSpace* mfem::DiscreteAdaptTC::coarse_tspec_fesv
protected

Definition at line 1539 of file tmop.hpp.

◆ good_tspec

bool mfem::DiscreteAdaptTC::good_tspec
protected

Definition at line 1553 of file tmop.hpp.

◆ good_tspec_grad

bool mfem::DiscreteAdaptTC::good_tspec_grad
protected

Definition at line 1553 of file tmop.hpp.

◆ good_tspec_hess

bool mfem::DiscreteAdaptTC::good_tspec_hess
protected

Definition at line 1553 of file tmop.hpp.

◆ Jtrcomp

DenseTensor mfem::DiscreteAdaptTC::Jtrcomp
mutableprotected

Definition at line 1534 of file tmop.hpp.

◆ lim_min_size

double mfem::DiscreteAdaptTC::lim_min_size
protected

Definition at line 1549 of file tmop.hpp.

◆ ncomp

int mfem::DiscreteAdaptTC::ncomp
protected

Definition at line 1515 of file tmop.hpp.

◆ orientationidx

int mfem::DiscreteAdaptTC::orientationidx
protected

Definition at line 1515 of file tmop.hpp.

◆ ptspec_fesv

ParFiniteElementSpace* mfem::DiscreteAdaptTC::ptspec_fesv
protected

Definition at line 1543 of file tmop.hpp.

◆ sizeidx

int mfem::DiscreteAdaptTC::sizeidx
protected

Definition at line 1515 of file tmop.hpp.

◆ skewidx

int mfem::DiscreteAdaptTC::skewidx
protected

Definition at line 1515 of file tmop.hpp.

◆ tspec

Vector mfem::DiscreteAdaptTC::tspec
protected

Definition at line 1516 of file tmop.hpp.

◆ tspec_derefine

Vector mfem::DiscreteAdaptTC::tspec_derefine
protected

Definition at line 1530 of file tmop.hpp.

◆ tspec_fesv

FiniteElementSpace* mfem::DiscreteAdaptTC::tspec_fesv
protected

Definition at line 1538 of file tmop.hpp.

◆ tspec_gf

GridFunction* mfem::DiscreteAdaptTC::tspec_gf
protected

Definition at line 1540 of file tmop.hpp.

◆ tspec_pert1h

Vector mfem::DiscreteAdaptTC::tspec_pert1h
protected

Definition at line 1518 of file tmop.hpp.

◆ tspec_pert2h

Vector mfem::DiscreteAdaptTC::tspec_pert2h
protected

Definition at line 1519 of file tmop.hpp.

◆ tspec_pertmix

Vector mfem::DiscreteAdaptTC::tspec_pertmix
protected

Definition at line 1520 of file tmop.hpp.

◆ tspec_pgf

ParGridFunction* mfem::DiscreteAdaptTC::tspec_pgf
protected

Definition at line 1545 of file tmop.hpp.

◆ tspec_refine

DenseMatrix mfem::DiscreteAdaptTC::tspec_refine
protected

Definition at line 1527 of file tmop.hpp.

◆ tspec_sav

Vector mfem::DiscreteAdaptTC::tspec_sav
protected

Definition at line 1517 of file tmop.hpp.


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