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 | 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 UpdateTargetSpecification (const Vector &new_x, bool use_flag=false)
 
void UpdateTargetSpecification (Vector &new_x, Vector &IntData)
 
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 use_flag=false)
 
void UpdateHessianTargetSpecification (const Vector &x, double dx, bool use_flag=false)
 
void SetAdaptivityEvaluator (AdaptivityEvaluator *ae)
 
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 ComputeElementTargetsGradient (const IntegrationRule &ir, const Vector &elfun, IsoparametricTransformation &Tpr, DenseTensor &dJtr) const
 
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 ()
 
void SetNodes (const GridFunction &n)
 Set 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...
 
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 ()
 
void SetTspecAtIndex (int idx, const ParGridFunction &tspec_)
 
void FinalizeParDiscreteTargetSpec (const ParGridFunction &tspec_)
 
- Protected Member Functions inherited from mfem::TargetConstructor
bool Parallel () const
 
bool Parallel () const
 
void ComputeAvgVolume () 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
 
DenseTensor Jtrcomp
 
const FiniteElementSpacetspec_fes
 
const FiniteElementSpacetspec_fesv
 
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
 
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 736 of file tmop.hpp.

Constructor & Destructor Documentation

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

Definition at line 777 of file tmop.hpp.

virtual mfem::DiscreteAdaptTC::~DiscreteAdaptTC ( )
inlinevirtual

Definition at line 786 of file tmop.hpp.

Member Function Documentation

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

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

Reimplemented from mfem::TargetConstructor.

Definition at line 1422 of file tmop.cpp.

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

Definition at line 1018 of file tmop.cpp.

void mfem::DiscreteAdaptTC::FinalizeSerialDiscreteTargetSpec ( )
protected

Definition at line 1165 of file tmop.cpp.

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

Definition at line 855 of file tmop.hpp.

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

Definition at line 856 of file tmop.hpp.

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

Definition at line 857 of file tmop.hpp.

void mfem::DiscreteAdaptTC::ResetUpdateFlags ( )
inline

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

Definition at line 819 of file tmop.hpp.

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

Definition at line 1223 of file tmop.cpp.

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

Definition at line 849 of file tmop.hpp.

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

Definition at line 1090 of file tmop.cpp.

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

Definition at line 1065 of file tmop.cpp.

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

Definition at line 1074 of file tmop.cpp.

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

Definition at line 1049 of file tmop.cpp.

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

Definition at line 1057 of file tmop.cpp.

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

Definition at line 1083 of file tmop.cpp.

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

Definition at line 1147 of file tmop.cpp.

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

Definition at line 1156 of file tmop.cpp.

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

Definition at line 1130 of file tmop.cpp.

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

Definition at line 1139 of file tmop.cpp.

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

Definition at line 1181 of file tmop.cpp.

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

Definition at line 1118 of file tmop.cpp.

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

Definition at line 1037 of file tmop.cpp.

void mfem::DiscreteAdaptTC::UpdateGradientTargetSpecification ( const Vector x,
double  dx,
bool  use_flag = false 
)

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

Definition at line 1769 of file tmop.cpp.

void mfem::DiscreteAdaptTC::UpdateHessianTargetSpecification ( const Vector x,
double  dx,
bool  use_flag = false 
)

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

Definition at line 1795 of file tmop.cpp.

void mfem::DiscreteAdaptTC::UpdateTargetSpecification ( const Vector new_x,
bool  use_flag = false 
)

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

Definition at line 1188 of file tmop.cpp.

void mfem::DiscreteAdaptTC::UpdateTargetSpecification ( Vector new_x,
Vector IntData 
)

Definition at line 1200 of file tmop.cpp.

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

Definition at line 1206 of file tmop.cpp.

Member Data Documentation

AdaptivityEvaluator* mfem::DiscreteAdaptTC::adapt_eval
protected

Definition at line 766 of file tmop.hpp.

int mfem::DiscreteAdaptTC::aspectratioidx
protected

Definition at line 741 of file tmop.hpp.

bool mfem::DiscreteAdaptTC::good_tspec
protected

Definition at line 762 of file tmop.hpp.

bool mfem::DiscreteAdaptTC::good_tspec_grad
protected

Definition at line 762 of file tmop.hpp.

bool mfem::DiscreteAdaptTC::good_tspec_hess
protected

Definition at line 762 of file tmop.hpp.

DenseTensor mfem::DiscreteAdaptTC::Jtrcomp
mutableprotected

Definition at line 753 of file tmop.hpp.

int mfem::DiscreteAdaptTC::ncomp
protected

Definition at line 741 of file tmop.hpp.

int mfem::DiscreteAdaptTC::orientationidx
protected

Definition at line 741 of file tmop.hpp.

int mfem::DiscreteAdaptTC::sizeidx
protected

Definition at line 741 of file tmop.hpp.

int mfem::DiscreteAdaptTC::skewidx
protected

Definition at line 741 of file tmop.hpp.

Vector mfem::DiscreteAdaptTC::tspec
protected

Definition at line 742 of file tmop.hpp.

const FiniteElementSpace* mfem::DiscreteAdaptTC::tspec_fes
protected

Definition at line 757 of file tmop.hpp.

const FiniteElementSpace* mfem::DiscreteAdaptTC::tspec_fesv
protected

Definition at line 758 of file tmop.hpp.

Vector mfem::DiscreteAdaptTC::tspec_pert1h
protected

Definition at line 744 of file tmop.hpp.

Vector mfem::DiscreteAdaptTC::tspec_pert2h
protected

Definition at line 745 of file tmop.hpp.

Vector mfem::DiscreteAdaptTC::tspec_pertmix
protected

Definition at line 746 of file tmop.hpp.

Vector mfem::DiscreteAdaptTC::tspec_sav
protected

Definition at line 743 of file tmop.hpp.


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