MFEM
v4.4.0
Finite element discretization library
|
#include <tmop.hpp>
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... | ||||
FiniteElementSpace * | GetTSpecFESpace () | |||
Get the FESpace associated with tspec. More... | ||||
GridFunction * | GetTSpecData () | |||
Get the entire tspec. More... | ||||
void | UpdateAfterMeshTopologyChange () | |||
Update all discrete fields based on tspec and update for AMR. More... | ||||
ParFiniteElementSpace * | GetTSpecParFESpace () | |||
void | ParUpdateAfterMeshTopologyChange () | |||
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 Vector & | GetTspecPert1H () | |||
const Vector & | GetTspecPert2H () | |||
const Vector & | GetTspecPertMixH () | |||
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.
| ||||
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 GridFunction * | GetNodes () 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 |
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 |
FiniteElementSpace * | tspec_fesv |
FiniteElementSpace * | coarse_tspec_fesv |
GridFunction * | tspec_gf |
ParFiniteElementSpace * | ptspec_fesv |
ParGridFunction * | tspec_pgf |
int | amr_el |
double | lim_min_size |
bool | good_tspec |
bool | good_tspec_grad |
bool | good_tspec_hess |
AdaptivityEvaluator * | adapt_eval |
Protected Attributes inherited from mfem::TargetConstructor | |
const GridFunction * | nodes |
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... | |
|
inline |
|
virtual |
Computes reference-to-target transformation Jacobians for all quadrature points in all elements.
[in] | fes | The nodal FE space |
[in] | ir | The quadrature rule to use for all elements |
[in] | xe | E-vector with the current physical coordinates/positions; this parameter is used only when needed by the target constructor, see UsesPhysicalCoordinates() |
[out] | Jtr | The computed ref->target Jacobian matrices. |
Reimplemented from mfem::TargetConstructor.
|
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.
|
virtual |
Reimplemented from mfem::TargetConstructor.
|
protected |
|
protected |
void mfem::DiscreteAdaptTC::GetDiscreteTargetSpec | ( | GridFunction & | tspec_, |
int | idx | ||
) |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
void mfem::DiscreteAdaptTC::ParUpdateAfterMeshTopologyChange | ( | ) |
|
inline |
|
inline |
|
inline |
void mfem::DiscreteAdaptTC::RestoreTargetSpecificationAtNode | ( | ElementTransformation & | T, |
int | nodenum | ||
) |
|
inline |
|
protected |
|
inline |
|
virtual |
|
virtual |
|
virtual |
|
virtual |
|
virtual |
|
inline |
|
virtual |
|
virtual |
|
virtual |
|
virtual |
|
virtual |
|
protected |
|
protected |
void mfem::DiscreteAdaptTC::SetTspecDataForDerefinement | ( | FiniteElementSpace * | fes | ) |
void mfem::DiscreteAdaptTC::SetTspecFromIntRule | ( | int | e_id, |
const IntegrationRule & | intrule | ||
) |
void mfem::DiscreteAdaptTC::UpdateAfterMeshTopologyChange | ( | ) |
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.
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.
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.
void mfem::DiscreteAdaptTC::UpdateTargetSpecificationAtNode | ( | const FiniteElement & | el, |
ElementTransformation & | T, | ||
int | nodenum, | ||
int | idir, | ||
const Vector & | IntData | ||
) |
|
protected |
|
protected |
|
mutableprotected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |