15 #include "../linalg/invariants.hpp" 78 virtual int Id()
const {
return 0; }
81 class TargetConstructor;
132 MFEM_VERIFY(
tmop_q_arr.Size() == weights.
Size(),
"Incorrect #weights");
182 double detT_ep_ = 0.0001,
183 double muT_ep_ = 0.0001,
190 "Worst-case optimization has not been fully developed!");
194 MFEM_VERIFY(m_id == 4 || m_id == 14 || m_id == 66,
195 "Incorrect input barrier metric -- must be 4 / 14 / 66");
202 { MFEM_ABORT(
"Not implemented"); }
206 { MFEM_ABORT(
"Not implemented"); }
235 virtual int Id()
const {
return 1; }
246 { MFEM_ABORT(
"Not implemented"); }
250 { MFEM_ABORT(
"Not implemented"); }
261 { MFEM_ABORT(
"Not implemented"); }
265 { MFEM_ABORT(
"Not implemented"); }
276 { MFEM_ABORT(
"Not implemented"); }
280 { MFEM_ABORT(
"Not implemented"); }
291 { MFEM_ABORT(
"Not implemented"); }
295 { MFEM_ABORT(
"Not implemented"); }
317 virtual int Id()
const {
return 2; }
336 virtual int Id()
const {
return 4; }
354 virtual int Id()
const {
return 7; }
381 { MFEM_ABORT(
"Not implemented"); }
385 { MFEM_ABORT(
"Not implemented"); }
500 virtual int Id()
const {
return 66; }
525 virtual int Id()
const {
return 77; }
545 virtual int Id()
const {
return 80; }
559 { MFEM_ABORT(
"Not implemented"); }
563 { MFEM_ABORT(
"Not implemented"); }
583 virtual int Id()
const {
return 90; }
604 virtual int Id()
const {
return 94; }
616 { MFEM_ABORT(
"Not implemented"); }
620 { MFEM_ABORT(
"Not implemented"); }
699 virtual int Id()
const {
return 302; }
720 virtual int Id()
const {
return 303; }
741 virtual int Id()
const {
return 304; }
781 virtual int Id()
const {
return 313; }
799 virtual int Id()
const {
return 315; }
839 virtual int Id()
const {
return 318; }
860 virtual int Id()
const {
return 321; }
881 virtual int Id()
const {
return 322; }
902 virtual int Id()
const {
return 323; }
922 virtual int Id()
const {
return 328; }
941 virtual int Id()
const {
return 332; }
982 virtual int Id()
const {
return 334; }
1005 virtual int Id()
const {
return 338; }
1025 virtual int Id()
const {
return 347; }
1068 virtual int Id()
const {
return 360; }
1083 { MFEM_ABORT(
"Not implemented"); }
1087 { MFEM_ABORT(
"Not implemented"); }
1101 { MFEM_ABORT(
"Not implemented"); }
1105 { MFEM_ABORT(
"Not implemented"); }
1119 { MFEM_ABORT(
"Not implemented"); }
1123 { MFEM_ABORT(
"Not implemented"); }
1137 { MFEM_ABORT(
"Not implemented"); }
1141 { MFEM_ABORT(
"Not implemented"); }
1171 virtual double Eval(
const Vector &x,
const Vector &x0,
double d)
const = 0;
1193 MFEM_ASSERT(x.
Size() == x0.
Size(),
"Bad input.");
1201 MFEM_ASSERT(x.
Size() == x0.
Size(),
"Bad input.");
1204 subtract(1.0 / (dist * dist), x, x0, d1);
1210 MFEM_ASSERT(x.
Size() == x0.
Size(),
"Bad input.");
1212 d2.
Diag(1.0 / (dist * dist), x.
Size());
1224 MFEM_ASSERT(x.
Size() == x0.
Size(),
"Bad input.");
1232 MFEM_ASSERT(x.
Size() == x0.
Size(),
"Bad input.");
1235 double dist_squared = dist*dist;
1237 dist_squared, x, x0, d1);
1243 MFEM_ASSERT(x.
Size() == x0.
Size(),
"Bad input.");
1246 double dist_squared = dist*dist;
1247 double dist_squared_squared = dist_squared*dist_squared;
1252 d2(0,0) = ((400.0*tmp(0)*tmp(0)*
f)/dist_squared_squared)+(20.0*
f/dist_squared);
1253 d2(1,1) = ((400.0*tmp(1)*tmp(1)*
f)/dist_squared_squared)+(20.0*
f/dist_squared);
1254 d2(0,1) = (400.0*tmp(0)*tmp(1)*
f)/dist_squared_squared;
1259 d2(0,2) = (400.0*tmp(0)*tmp(2)*
f)/dist_squared_squared;
1260 d2(1,2) = (400.0*tmp(1)*tmp(2)*
f)/dist_squared_squared;
1263 d2(2,2) = ((400.0*tmp(2)*tmp(2)*
f)/dist_squared_squared)+(20.0*
f/dist_squared);
1271 class FiniteElementCollection;
1272 class FiniteElementSpace;
1273 class ParFiniteElementSpace;
1311 const Vector &init_field) = 0;
1386 comm = MPI_COMM_NULL;
1507 class ParGridFunction;
1638 int nodenum,
int idir,
1648 bool reuse_flag =
false,
1655 bool reuse_flag =
false,
1724 class TMOPNewtonSolver;
1858 double &metric_energy,
double &lim_energy,
1859 double &surf_fit_gf_energy);
1899 Vector &elfun,
const int nodenum,
const int idir,
1900 const double baseenergy,
bool update_stored);
1991 {
PA.enabled =
false; }
2016 "This function must be called before EnableNormalization, as " 2017 "the normalization computations must know how to integrate.");
2142 double &err_avg,
double &err_max);
2250 for (
int i = 0; i <
tmopi.Size(); i++) {
delete tmopi[i]; }
2308 const TargetConstructor &tc,
2309 const Mesh &mesh, GridFunction &metric_gf);
Abstract class for all finite elements.
virtual int Id() const
Return the metric ID.
void ComputeUntangleMetricQuantiles(const Vector &x, const FiniteElementSpace &fes)
InvariantsEvaluator2D< double > ie
GridFunction * surf_fit_hess
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 quadratu...
virtual double GetDerefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
virtual void SetParDiscreteTargetSize(const ParGridFunction &tspec_)
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 quadratu...
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Shifted barrier form of metric 56 (area, ideal barrier metric), 2D.
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
3D shifted barrier form of metric 316 (not typed).
void GetWeights(Array< double > &weights) const
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...
InvariantsEvaluator3D< double > ie
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)=0
const GridFunction * surf_fit_pos
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const =0
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
InvariantsEvaluator3D< double > ie
void AssembleElemVecSurfFit(const FiniteElement &el_x, IsoparametricTransformation &Tpr, DenseMatrix &mat)
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Class for an integration rule - an Array of IntegrationPoint.
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
VectorCoefficient * vector_tspec
3D non-barrier Shape (S) metric.
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
void EnableFiniteDifferences(const GridFunction &x)
Enables FD-based approximation and computes dx.
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
void UpdateTargetSpecification(const Vector &new_x, bool reuse_flag=false, int new_x_ordering=Ordering::byNODES)
virtual int Id() const
Return the metric ID.
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
const TargetType target_type
virtual int Id() const
Return the metric ID.
Class for grid function - Vector with associated FE space.
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
void ComputeAvgMetrics(const GridFunction &nodes, const TargetConstructor &tc, Vector &averages) const
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
3D barrier Shape+Size (VS) metric, well-posed (invex).
InvariantsEvaluator3D< double > ie
InvariantsEvaluator3D< double > ie
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
InvariantsEvaluator2D< double > ie
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
double GetSurfaceFittingWeight()
Get the surface fitting weight.
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
Base class for vector Coefficients that optionally depend on time and space.
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual ~TMOP_Metric_334()
virtual double Eval(const Vector &x, const Vector &x0, double dist) const
Returns the limiting function, f(x, x0, d).
void ComputeBalancedWeights(const GridFunction &nodes, const TargetConstructor &tc, Vector &weights) const
const GridFunction * lim_dist
virtual void SetSerialDiscreteTargetOrientation(const GridFunction &tspec_)
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
Compute the local energy.
void AddMultGradPA_2D(const Vector &, Vector &) const
bool UsesPhysicalCoordinates() const
Return true if the methods ComputeElementTargets(), ComputeAllElementTargets(), and ComputeElementTar...
struct mfem::TMOP_Integrator::@23 PA
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field, int new_nodes_ordering=Ordering::byNODES)=0
virtual ~TMOP_Metric_333()
virtual void EvalGrad(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip, int comp)=0
Evaluate the derivative of the matrix coefficient with respect to comp in the element described by T ...
virtual int Id() const
Return the metric ID.
virtual int Id() const
Return the metric ID.
InvariantsEvaluator2D< double > ie
void EnableNormalization(const GridFunction &x)
Computes the normalization factors of the metric and limiting integrals using the mesh position given...
TMOP_QualityMetric & GetAMRQualityMetric()
void AssembleGradPA_2D(const Vector &) const
void ComputeAllElementTargets(const Vector &xe=Vector()) const
const IntegrationRule & GradientIntegrationRule(const FiniteElement &el) const
void SetSize(int s)
Resize the vector to size s.
void InterpolateTMOP_QualityMetric(TMOP_QualityMetric &metric, const TargetConstructor &tc, const Mesh &mesh, GridFunction &metric_gf)
Interpolates the metric's values at the nodes of metric_gf.
3D barrier Shape (S) metric, well-posed (polyconvex & invex).
TMOP_QualityMetric * metric
void SetAdaptivityEvaluator(AdaptivityEvaluator *ae)
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual ~TMOP_Metric_332()
InvariantsEvaluator2D< double > ie
virtual ~TMOP_LimiterFunction()
Virtual destructor.
Base class for limiting functions to be used in class TMOP_Integrator.
virtual void ComputeElementTargetsGradient(const IntegrationRule &ir, const Vector &elfun, IsoparametricTransformation &Tpr, DenseTensor &dJtr) const
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual double EvalW(const DenseMatrix &Jpt) const =0
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Array< int > surf_fit_marker_dof_index
virtual void Eval_d2(const Vector &x, const Vector &x0, double dist, DenseMatrix &d2) const =0
Returns the Hessian of the limiting function f(x, x0, d) with respect to x.
virtual int Id() const
Return the metric ID.
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
3D non-barrier Skew metric.
const Array< TMOP_Integrator * > & GetTMOPIntegrators() const
void ParUpdateAfterMeshTopologyChange()
IntegrationRules * IntegRules
InvariantsEvaluator2D< double > ie
TMOP_QualityMetric * sz_metric
TMOP_QualityMetric * sh_metric
DiscreteAdaptTC * GetDiscreteAdaptTC() const
void SetFDhScale(double dxscale_)
virtual void Eval_d2(const Vector &x, const Vector &x0, double dist, DenseMatrix &d2) const
Returns the Hessian of the limiting function f(x, x0, d) with respect to x.
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const =0
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
InvariantsEvaluator2D< double > ie
int Size() const
Returns the size of the vector.
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
void SetMinSizeForTargets(double min_size_)
void AddMultPA_3D(const Vector &, Vector &) const
virtual int Id() const
Return the metric ID.
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Container class for integration rules.
ParFiniteElementSpace * pfes
Data type dense matrix using column-major storage.
3D barrier Shape (S) metric, well-posed (polyconvex & invex).
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
void SetVolumeScale(double vol_scale)
Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
void EnableAdaptiveLimiting(const GridFunction &z0, Coefficient &coeff, AdaptivityEvaluator &ae)
Restriction of the node positions to certain regions.
TMOP_QualityMetric * sz_metric
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual int Id() const
Return the metric ID.
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
InvariantsEvaluator3D< double > ie
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
void EnableNormalization(const GridFunction &x)
Normalization factor that considers all integrators in the combination.
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
TMOP_QualityMetric & tmop_metric
3D compound barrier Shape+Size (VS) metric (polyconvex, balanced).
2D barrier Shape+Size (VS) metric (polyconvex).
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
bool IsSurfaceFittingEnabled()
3D barrier Shape+Size (VS) metric, well-posed (invex).
const AdaptivityEvaluator * GetAdaptivityEvaluator() const
Abstract parallel finite element space.
virtual void AssemblePA(const FiniteElementSpace &)
Method defining partial assembly.
TMOP_Metric_080(double gamma)
InvariantsEvaluator2D< double > ie
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
void SetSerialMetaInfo(const Mesh &m, const FiniteElementSpace &f)
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
const IntegrationRule * ir
void SetExactActionFlag(bool flag_)
Flag to control if exact action of Integration is effected.
InvariantsEvaluator2D< double > ie
void AssembleElementVectorFD(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
const GridFunction * nodes
TMOP_WorstCaseUntangleOptimizer_Metric(TMOP_QualityMetric &tmop_metric_, int exponent_=1, double alpha_=1.5, double detT_ep_=0.0001, double muT_ep_=0.0001, BarrierType btype_=BarrierType::None, WorstCaseType wctype_=WorstCaseType::None)
virtual int Id() const
Return the metric ID.
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void SetSerialDiscreteTargetSkew(const GridFunction &tspec_)
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
void UpdateSurfaceFittingWeight(double factor)
Update the surface fitting weight as surf_fit_coeff *= factor;.
AdaptivityEvaluator * surf_fit_eval_bg_grad
void SetRefinementSubElement(int amr_el_)
2D non-barrier Aspect ratio metric.
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
void UpdateAfterMeshPositionChange(const Vector &x_new, const FiniteElementSpace &x_fes)
InvariantsEvaluator2D< double > ie
virtual double GetLocalStateEnergyPA(const Vector &) const
Compute the local (to the MPI rank) energy with partial assembly.
void IntegrateOverTarget(bool integ_over_target_)
2D barrier Shape+Size (VS) metric (not polyconvex).
virtual ~TMOP_Metric_090()
TMOP_QualityMetric * sh_metric
Default limiter function in TMOP_Integrator.
virtual double GetLocalStateEnergyPA(const Vector &) const
Compute the local (to the MPI rank) energy with partial assembly.
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
void AssembleElementGradExact(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
virtual void AddMultGradPA(const Vector &, Vector &) const
Method for partially assembled gradient action.
Array< Vector * > ElemDer
3D barrier Shape+Size (VS) metric, well-posed (invex).
virtual ~TMOP_AMetric_126()
const Vector & GetTspecPertMixH()
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
void EnableLimiting(const GridFunction &n0, const GridFunction &dist, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
Adds the limiting term to the first integrator. Disables it for the rest.
InvariantsEvaluator2D< double > ie
TMOP_QualityMetric * sz_metric
2D Shifted barrier form of shape metric (mu_2).
TMOP_Metric_313(double &mindet)
TMOP_Metric_352(double &t0)
Coefficient * scalar_tspec
virtual void SetSerialDiscreteTargetAspectRatio(const GridFunction &tspec_)
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
TMOPMatrixCoefficient(int dim)
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
2D barrier shape (S) metric (not polyconvex).
virtual ~TMOP_QualityMetric()
std::function< double(const Vector &)> f(double mass_coeff)
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
const IntegrationRule & EnergyIntegrationRule(const FiniteElement &el) const
void UpdateGradientTargetSpecification(const Vector &x, double dx, bool reuse_flag=false, int x_ordering=Ordering::byNODES)
virtual void AddMultPA(const Vector &, Vector &) const
Method for partially assembled action.
double ComputeUntanglerMaxMuBarrier(const Vector &x, const FiniteElementSpace &fes)
TMOP_QualityMetric * sh_metric
void AddMultPA_C0_2D(const Vector &, Vector &) const
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
const GeometricFactors * geom
3D barrier Shape (S) metric, well-posed (polyconvex & invex).
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
bool ComputeAllElementTargets(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
double GetLocalStateEnergyPA_C0_3D(const Vector &) const
2D barrier Size (V) metric (polyconvex).
Abstract class for local mesh quality metrics in the target-matrix optimization paradigm (TMOP) by P...
void AssembleGradPA_C0_2D(const Vector &) const
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
void FinalizeSerialDiscreteTargetSpec(const GridFunction &tspec_)
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Coefficient * adapt_lim_coeff
virtual ~TMOPMatrixCoefficient()
virtual ~AdaptivityEvaluator()
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
InvariantsEvaluator3D< double > ie
2D compound barrier Shape+Size (VS) metric (balanced).
virtual void AddQualityMetric(TMOP_QualityMetric *tq, double wt=1.0)
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
void SetLimitingNodes(const GridFunction &n0)
Update the original/reference nodes used for limiting.
virtual ~TMOP_Metric_080()
void SetDiscreteTargetBase(const GridFunction &tspec_)
virtual ~TMOP_Metric_066()
AdaptivityEvaluator * adapt_lim_eval
virtual void AssembleGradPA(const Vector &, const FiniteElementSpace &)
Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at th...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual double EvalWBarrier(const DenseMatrix &Jpt) const
void AddMultGradPA_3D(const Vector &, Vector &) const
DiscreteAdaptTC * discr_tc
void EnableLimiting(const GridFunction &n0, const GridFunction &dist, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
Limiting of the mesh displacements (general version).
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
void AssembleDiagonalPA_3D(Vector &) const
virtual WorstCaseType GetWorstCaseType()
InvariantsEvaluator3D< double > ie
void SetCoefficient(Coefficient &w1)
Sets a scaling Coefficient for the quality metric term of the integrator.
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
void UpdateTargetSpecificationAtNode(const FiniteElement &el, ElementTransformation &T, int nodenum, int idir, const Vector &IntData)
TMOP_QualityMetric * sh_metric
void SetTspecFromIntRule(int e_id, const IntegrationRule &intrule)
virtual int Id() const
Return the metric ID.
int Append(const T &el)
Append element 'el' to array, resize if necessary.
TargetConstructor(TargetType ttype)
Constructor for use in serial.
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
DiscreteAdaptTC(TargetType ttype)
void AddMultGradPA_C0_2D(const Vector &, Vector &) const
TMOP_QualityMetric * sh_metric
void ResetRefinementTspecData()
InvariantsEvaluator2D< double > ie
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
void ClearGeometricFactors()
void SetNodes(const GridFunction &n)
Set the nodes to be used in the target-matrix construction.
InvariantsEvaluator3D< double > ie
void AssembleElementVectorExact(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
virtual int Id() const
Return the metric ID.
TMOP_QualityMetric * sz_metric
const TargetConstructor * targetC
Array< Vector * > ElemPertEnergy
const GridFunction * lim_nodes0
const ParGridFunction * adapt_lim_pgf0
TMOP_QualityMetric * sz_metric
TMOP_QualityMetric * sh_metric
TMOP_QualityMetric * sh_metric
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
InvariantsEvaluator2D< double > ie
virtual void Eval_d2(const Vector &x, const Vector &x0, double dist, DenseMatrix &d2) const
Returns the Hessian of the limiting function f(x, x0, d) with respect to x.
virtual int Id() const
Return the metric ID.
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual void AssembleGradPA(const Vector &, const FiniteElementSpace &)
Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at th...
InvariantsEvaluator3D< double > ie
virtual void SetParDiscreteTargetAspectRatio(const ParGridFunction &tspec_)
virtual void AssemblePA(const FiniteElementSpace &)
Method defining partial assembly.
void AssembleGradPA_C0_3D(const Vector &) const
virtual ~TMOP_ExponentialLimiter()
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual ~TargetConstructor()
virtual void ComputeElementTargetsGradient(const IntegrationRule &ir, const Vector &elfun, IsoparametricTransformation &Tpr, DenseTensor &dJtr) const
void ParEnableNormalization(const ParGridFunction &x)
InvariantsEvaluator2D< double > ie
virtual int Id() const
Return the metric ID.
TMOP_QualityMetric * sz_metric
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
InvariantsEvaluator2D< double > ie
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
GridFunction * GetTSpecData()
Get the entire tspec.
FiniteElementSpace * tspec_fesv
Coefficient * metric_coeff
void AssembleElemGradSurfFit(const FiniteElement &el_x, IsoparametricTransformation &Tpr, DenseMatrix &mat)
TMOP_QualityMetric * h_metric
Array< TMOP_QualityMetric * > tmop_q_arr
virtual int Id() const
Return the metric ID.
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual BarrierType GetBarrierType()
void AssembleDiagonalPA_2D(Vector &) const
void RestoreTargetSpecificationAtNode(ElementTransformation &T, int nodenum)
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
void AssemblePA_Limiting()
TMOP_QualityMetric * sh_metric
2D barrier Shape+Orientation (OS) metric (polyconvex).
void SetTspecAtIndex(int idx, const GridFunction &tspec_)
virtual double GetRefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, const IntegrationRule &irule)
Computes the mean of the energies of the given element's children.
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Array< int > surf_fit_dof_count
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
double GetLocalStateEnergyPA_3D(const Vector &) const
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
3D compound barrier Shape+Size (VS) metric (polyconvex, balanced).
InvariantsEvaluator2D< double > ie
AdaptivityEvaluator * surf_fit_eval
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
TMOP_Metric_066(double gamma)
virtual void Eval_d1(const Vector &x, const Vector &x0, double dist, Vector &d1) const
Returns the gradient of the limiting function f(x, x0, d) with respect to x.
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual ~TMOP_QuadraticLimiter()
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
FiniteElementSpace * GetTSpecFESpace()
Get the FESpace associated with tspec.
TMOP_QualityMetric * sh_metric
TMOP_QualityMetric * sz_metric
virtual void SetTargetJacobian(const DenseMatrix &Jtr_)
Specify the reference-element -> target-element Jacobian matrix for the point of interest.
InvariantsEvaluator2D< double > ie
TMOP_Metric_211(double epsilon=1e-4)
void AssembleElemGradAdaptLim(const FiniteElement &el, IsoparametricTransformation &Tpr, const IntegrationRule &ir, const Vector &weights, DenseMatrix &m)
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
const DofToQuad * maps_lim
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
3D Size (V) untangling metric.
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual void SetTargetJacobian(const DenseMatrix &Jtr_)
Specify the reference-element -> target-element Jacobian matrix for the point of interest.
void SetLimitingNodes(const GridFunction &n0)
Update the original/reference nodes used for limiting.
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
TMOP_AMetric_126(double gamma)
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
InvariantsEvaluator3D< double > ie
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
void AddMultPA_C0_3D(const Vector &, Vector &) const
void subtract(const Vector &x, const Vector &y, Vector &z)
TMOPMatrixCoefficient * matrix_tspec
void AssembleDiagonalPA_C0_2D(Vector &) const
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
const GridFunction * GetNodes() const
Get the nodes to be used in the target-matrix construction.
InvariantsEvaluator2D< double > ie
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual bool ContainsVolumeInfo() const
Checks if the target matrices contain non-trivial size specification.
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Base class for Matrix Coefficients that optionally depend on time and space.
InvariantsEvaluator2D< double > ie
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
void ReleasePADeviceMemory(bool copy_to_host=true)
virtual int Id() const
Return the metric ID.
virtual void Eval_d1(const Vector &x, const Vector &x0, double dist, Vector &d1) const
Returns the gradient of the limiting function f(x, x0, d) with respect to x.
void Clear()
Delete the matrix data array (if owned) and reset the matrix state.
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
TMOP_Metric_252(double &t0)
Note that t0 is stored by reference.
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
void SetParMetaInfo(const ParMesh &m, const ParFiniteElementSpace &f)
Parallel version of SetSerialMetaInfo.
GridFunction * surf_fit_grad
void ResetUpdateFlags()
Used in combination with the Update methods to avoid extra computations.
virtual int Id() const
Return the metric ID.
AnalyticAdaptTC(TargetType ttype)
InvariantsEvaluator3D< double > ie
void AddTMOPIntegrator(TMOP_Integrator *ti)
Adds a new TMOP_Integrator to the combination.
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
2D non-barrier Skew metric.
virtual void AssembleGradDiagonalPA(Vector &) const
Method for computing the diagonal of the gradient with partial assembly.
virtual int Id() const
Return the metric ID.
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
void UpdateAfterMeshTopologyChange()
virtual void SetMaxMuT(double max_muT_)
virtual void SetMinDetT(double min_detT_)
InvariantsEvaluator3D< double > ie
2D barrier Shape+Size (VS) metric (not polyconvex).
void SetTransformation(ElementTransformation &)
The method HyperelasticModel::SetTransformation() is hidden for TMOP_QualityMetrics, because it is not used.
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 quadratu...
double GetLocalStateEnergyPA_C0_2D(const Vector &) const
void ComputeFDh(const Vector &x, const FiniteElementSpace &fes)
Determines the perturbation, h, for FD-based approximation.
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
void ComputeNormalizationEnergies(const GridFunction &x, double &metric_energy, double &lim_energy, double &surf_fit_gf_energy)
virtual double Eval(const Vector &x, const Vector &x0, double d) const =0
Returns the limiting function, f(x, x0, d).
Coefficient * surf_fit_coeff
Class for integration point with weight.
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual ~TMOP_Metric_347()
AdaptivityEvaluator * adapt_eval
void SetWeights(const Vector &weights)
Changes the weights of the metrics in the combination.
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
ParFiniteElementSpace * ptspec_fesv
InvariantsEvaluator3D< double > ie
const FiniteElementSpace * fes
Array< TMOP_Integrator * > tmopi
TMOP_QualityMetric * sz_metric
void ComputeMinJac(const Vector &x, const FiniteElementSpace &fes)
void AssembleDiagonalPA_C0_3D(Vector &) const
TMOP_Integrator(TMOP_QualityMetric *m, TargetConstructor *tc, TMOP_QualityMetric *hm)
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
const IntegrationRule & ActionIntegrationRule(const FiniteElement &el) const
void Diag(double c, int n)
Creates n x n diagonal matrix with diagonal elements c.
virtual void AddMultGradPA(const Vector &, Vector &) const
Method for partially assembled gradient action.
double GetFDDerivative(const FiniteElement &el, ElementTransformation &T, Vector &elfun, const int nodenum, const int idir, const double baseenergy, bool update_stored)
Exponential limiter function in TMOP_Integrator.
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
void EnableSurfaceFittingFromSource(const ParGridFunction &s_bg, ParGridFunction &s0, const Array< bool > &smarker, Coefficient &coeff, AdaptivityEvaluator &ae, const ParGridFunction &s_bg_grad, ParGridFunction &s0_grad, AdaptivityEvaluator &age, const ParGridFunction &s_bg_hess, ParGridFunction &s0_hess, AdaptivityEvaluator &ahe)
Fitting of certain DOFs in the current mesh to the zero level set of a function defined on another (f...
3D Shape (S) metric, untangling version of 303.
virtual int Id() const
Return the metric ID.
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
Computes the integral of W(Jacobian(Trt)) over a target zone.
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
void AssembleElemVecAdaptLim(const FiniteElement &el, IsoparametricTransformation &Tpr, const IntegrationRule &ir, const Vector &weights, DenseMatrix &mat)
void ParUpdateAfterMeshTopologyChange()
virtual double GetDerefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
void ComputeAvgVolume() const
2D compound barrier Shape+Size (VS) metric (balanced).
void Destroy()
Destroy a vector.
virtual void Eval_d1(const Vector &x, const Vector &x0, double dist, Vector &d1) const =0
Returns the gradient of the limiting function f(x, x0, d) with respect to x.
AdaptivityEvaluator * surf_fit_eval_bg_hess
TMOP_Metric_347(double gamma)
virtual int Id() const
Return the metric ID.
void GetSurfaceFittingErrors(const Vector &pos, double &err_avg, double &err_max)
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
3D non-barrier Aspect ratio metric.
const Vector & GetTspecPert2H()
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
const Vector & GetTspecPert1H()
void AssembleElementGradFD(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
ParFiniteElementSpace * GetTSpecParFESpace()
TMOP_QualityMetric * sz_metric
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
TMOP_QuadraticLimiter * surf_fit_limiter
InvariantsEvaluator2D< double > ie
virtual void SetParDiscreteTargetSpec(const ParGridFunction &tspec_)
InvariantsEvaluator2D< double > ie
TMOP_QualityMetric * sh_metric
Abstract class for hyperelastic models.
virtual int Id() const
Return the metric ID.
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
void UpdateHessianTargetSpecification(const Vector &x, double dx, bool reuse_flag=false, int x_ordering=Ordering::byNODES)
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
GridFunction * adapt_lim_gf
InvariantsEvaluator2D< double > ie
TMOP_Metric_333(double gamma)
TargetConstructor(TargetType ttype, MPI_Comm mpicomm)
Constructor for use in parallel.
virtual int Id() const
Return the metric ID.
GridFunction * surf_fit_gf
void AssembleGradPA_3D(const Vector &) const
double GetLocalStateEnergyPA_2D(const Vector &) const
TMOP_Integrator(TMOP_QualityMetric *m, TargetConstructor *tc)
void FinalizeParDiscreteTargetSpec(const ParGridFunction &tspec_)
ParGridFunction * tspec_pgf
void GetDiscreteTargetSpec(GridFunction &tspec_, int idx)
Get one of the discrete fields from tspec.
const Array< bool > * surf_fit_marker
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
TMOP_Metric_332(double gamma)
virtual int Id() const
Return the metric ID.
TMOP_Metric_022(double &t0)
3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
InvariantsEvaluator3D< double > ie
virtual void AssembleGradDiagonalPA(Vector &) const
Method for computing the diagonal of the gradient with partial assembly.
void SetTspecDataForDerefinement(FiniteElementSpace *fes)
Computes target specification data with respect to the coarse FE space.
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
virtual ~TMOP_Metric_094()
TargetType
Target-matrix construction algorithms supported by this class.
virtual ~TMOP_Metric_338()
void EnableSurfaceFitting(const GridFunction &s0, const Array< bool > &smarker, Coefficient &coeff, AdaptivityEvaluator &ae)
Fitting of certain DOFs to the zero level set of a function.
InvariantsEvaluator3D< double > ie
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
TMOP_Metric_311(double epsilon=1e-4)
virtual void SetSerialDiscreteTargetSpec(const GridFunction &tspec_)
InvariantsEvaluator3D< double > ie
virtual void SetParDiscreteTargetSkew(const ParGridFunction &tspec_)
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
virtual void ComputeElementTargetsGradient(const IntegrationRule &ir, const Vector &elfun, IsoparametricTransformation &Tpr, DenseTensor &dJtr) const
virtual double GetRefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, const IntegrationRule &irule)
InvariantsEvaluator2D< double > ie
virtual void SetParDiscreteTargetOrientation(const ParGridFunction &tspec_)
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...
Class for parallel grid function.
virtual int Id() const
Return the metric ID.
Base class representing target-matrix construction algorithms for mesh optimization via the target-ma...
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
virtual ~DiscreteAdaptTC()
TMOP_Metric_334(double gamma)
void AddMultPA_2D(const Vector &, Vector &) const
virtual void SetAnalyticTargetSpec(Coefficient *sspec, VectorCoefficient *vspec, TMOPMatrixCoefficient *mspec)
3D compound barrier Shape+Size (VS) metric (polyconvex).
Rank 3 tensor (array of matrices)
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Class for parallel meshes.
TMOP_QualityMetric * sh_metric
virtual void SetSerialDiscreteTargetSize(const GridFunction &tspec_)
TMOP_QualityMetric * sz_metric
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
TMOP_LimiterFunction * lim_func
InvariantsEvaluator3D< double > ie
const GridFunction * adapt_lim_gf0
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
double ComputeMinDetT(const Vector &x, const FiniteElementSpace &fes)
virtual void AddMultPA(const Vector &, Vector &) const
Method for partially assembled action.
double DistanceSquaredTo(const double *p) const
Compute the square of the Euclidean distance to another vector.
2D barrier Shape+Orientation (OS) metric (polyconvex).
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
void ParEnableNormalization(const ParGridFunction &x)
virtual int Id() const
Return the metric ID.
void ComputeAllElementTargets_Fallback(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
void UpdateAfterMeshTopologyChange()
Update all discrete fields based on tspec and update for AMR.
void ResetDerefinementTspecData()
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
3D barrier Shape (S) metric, well-posed (polyconvex & invex).
virtual double Eval(const Vector &x, const Vector &x0, double dist) const
Returns the limiting function, f(x, x0, d).
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
void AddMultGradPA_C0_3D(const Vector &, Vector &) const
virtual int Id() const
Return the metric ID.
2D non-barrier metric without a type.
TMOP_QualityMetric * sz_metric
TargetType GetTargetType() const
InvariantsEvaluator3D< double > ie
InvariantsEvaluator2D< double > ie
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual ~TMOP_Metric_328()
FiniteElementSpace * coarse_tspec_fesv
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
2D non-barrier Shape+Size+Orientation (VOS) metric (polyconvex).