MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
tmop.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_TMOP_HPP
13 #define MFEM_TMOP_HPP
14 
15 #include "../linalg/invariants.hpp"
16 #include "nonlininteg.hpp"
17 
18 namespace mfem
19 {
20 
21 /** @brief Abstract class for local mesh quality metrics in the target-matrix
22  optimization paradigm (TMOP) by P. Knupp et al. */
24 {
25 protected:
26  const DenseMatrix *Jtr; /**< Jacobian of the reference-element to
27  target-element transformation. */
28 
29  /** @brief The method SetTransformation() is hidden for TMOP_QualityMetric%s,
30  because it is not used. */
32 
33 public:
34  TMOP_QualityMetric() : Jtr(NULL) { }
35  virtual ~TMOP_QualityMetric() { }
36 
37  /** @brief Specify the reference-element -> target-element Jacobian matrix
38  for the point of interest.
39 
40  The specified Jacobian matrix, #Jtr, can be used by metrics that cannot
41  be written just as a function of the target->physical Jacobian matrix,
42  Jpt. */
43  void SetTargetJacobian(const DenseMatrix &_Jtr) { Jtr = &_Jtr; }
44 
45  /** @brief Evaluate the strain energy density function, W = W(Jpt).
46  @param[in] Jpt Represents the target->physical transformation
47  Jacobian matrix. */
48  virtual double EvalW(const DenseMatrix &Jpt) const = 0;
49 
50  /** @brief Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
51  @param[in] Jpt Represents the target->physical transformation
52  Jacobian matrix.
53  @param[out] P The evaluated 1st Piola-Kirchhoff stress tensor. */
54  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const = 0;
55 
56  /** @brief Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor
57  and assemble its contribution to the local gradient matrix 'A'.
58  @param[in] Jpt Represents the target->physical transformation
59  Jacobian matrix.
60  @param[in] DS Gradient of the basis matrix (dof x dim).
61  @param[in] weight Quadrature weight coefficient for the point.
62  @param[in,out] A Local gradient matrix where the contribution from this
63  point will be added.
64 
65  Computes weight * d(dW_dxi)_d(xj) at the current point, for all i and j,
66  where x1 ... xn are the FE dofs. This function is usually defined using
67  the matrix invariants and their derivatives.
68  */
69  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
70  const double weight, DenseMatrix &A) const = 0;
71 };
72 
73 
74 /// Metric without a type, 2D
76 {
77 protected:
79 
80 public:
81  // W = |J|^2.
82  virtual double EvalW(const DenseMatrix &Jpt) const;
83 
84  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
85 
86  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
87  const double weight, DenseMatrix &A) const;
88 };
89 
90 /// Skew metric, 2D.
92 {
93 public:
94  // W = 0.5 (1 - cos(angle_Jpr - angle_Jtr)).
95  virtual double EvalW(const DenseMatrix &Jpt) const;
96 
97  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
98  { MFEM_ABORT("Not implemented"); }
99 
100  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
101  const double weight, DenseMatrix &A) const
102  { MFEM_ABORT("Not implemented"); }
103 };
104 
105 /// Skew metric, 3D.
107 {
108 public:
109  // W = 1/6 (3 - sum_i cos(angle_Jpr_i - angle_Jtr_i)), i = 1..3.
110  virtual double EvalW(const DenseMatrix &Jpt) const;
111 
112  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
113  { MFEM_ABORT("Not implemented"); }
114 
115  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
116  const double weight, DenseMatrix &A) const
117  { MFEM_ABORT("Not implemented"); }
118 };
119 
120 /// Aspect ratio metric, 2D.
122 {
123 public:
124  // W = 0.5 (ar_Jpr/ar_Jtr + ar_Jtr/ar_Jpr) - 1.
125  virtual double EvalW(const DenseMatrix &Jpt) const;
126 
127  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
128  { MFEM_ABORT("Not implemented"); }
129 
130  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
131  const double weight, DenseMatrix &A) const
132  { MFEM_ABORT("Not implemented"); }
133 };
134 
135 /// Aspect ratio metric, 3D.
137 {
138 public:
139  // W = 1/3 sum [0.5 (ar_Jpr_i/ar_Jtr_i + ar_Jtr_i/ar_Jpr_i) - 1], i = 1..3.
140  virtual double EvalW(const DenseMatrix &Jpt) const;
141 
142  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
143  { MFEM_ABORT("Not implemented"); }
144 
145  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
146  const double weight, DenseMatrix &A) const
147  { MFEM_ABORT("Not implemented"); }
148 };
149 
150 /// Shape, ideal barrier metric, 2D
152 {
153 protected:
155 
156 public:
157  // W = 0.5|J|^2 / det(J) - 1.
158  virtual double EvalW(const DenseMatrix &Jpt) const;
159 
160  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
161 
162  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
163  const double weight, DenseMatrix &A) const;
164 };
165 
166 /// Shape & area, ideal barrier metric, 2D
168 {
169 protected:
171 
172 public:
173  // W = |J - J^-t|^2.
174  virtual double EvalW(const DenseMatrix &Jpt) const;
175 
176  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
177 
178  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
179  const double weight, DenseMatrix &A) const;
180 };
181 
182 /// Shape & area metric, 2D
184 {
185 protected:
187 
188 public:
189  // W = det(J) * |J - J^-t|^2.
190  virtual double EvalW(const DenseMatrix &Jpt) const;
191 
192  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
193 
194  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
195  const double weight, DenseMatrix &A) const;
196 };
197 
198 /// Shifted barrier form of metric 2 (shape, ideal barrier metric), 2D
200 {
201 protected:
202  double &tau0;
204 
205 public:
206  TMOP_Metric_022(double &t0): tau0(t0) {}
207 
208  // W = 0.5(|J|^2 - 2det(J)) / (det(J) - tau0).
209  virtual double EvalW(const DenseMatrix &Jpt) const;
210 
211  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
212 
213  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
214  const double weight, DenseMatrix &A) const;
215 };
216 
217 /// Shape, ideal barrier metric, 2D
219 {
220 protected:
222 
223 public:
224  // W = 0.5|J^t J|^2 / det(J)^2 - 1.
225  virtual double EvalW(const DenseMatrix &Jpt) const;
226 
227  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
228 
229  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
230  const double weight, DenseMatrix &A) const;
231 };
232 
233 /// Area metric, 2D
235 {
236 protected:
238 
239 public:
240  // W = (det(J) - 1)^2.
241  virtual double EvalW(const DenseMatrix &Jpt) const;
242 
243  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
244 
245  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
246  const double weight, DenseMatrix &A) const;
247 
248 };
249 
250 /// Area, ideal barrier metric, 2D
252 {
253 protected:
255 
256 public:
257  // W = 0.5( sqrt(det(J)) - 1 / sqrt(det(J)) )^2
258  // = 0.5( det(J) - 1 )^2 / det(J)
259  // = 0.5( det(J) + 1/det(J) ) - 1.
260  virtual double EvalW(const DenseMatrix &Jpt) const;
261 
262  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
263 
264  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
265  const double weight, DenseMatrix &A) const;
266 
267 };
268 
269 /// Shape, ideal barrier metric, 2D
271 {
272 protected:
274 
275 public:
276  // W = |J^t J|^2 / det(J)^2 - 2|J|^2 / det(J) + 2
277  // = I1b (I1b - 2).
278  virtual double EvalW(const DenseMatrix &Jpt) const;
279 
280  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
281 
282  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
283  const double weight, DenseMatrix &A) const;
284 
285 };
286 
287 /// Area, ideal barrier metric, 2D
289 {
290 protected:
292 
293 public:
294  // W = 0.5(det(J) - 1 / det(J))^2.
295  virtual double EvalW(const DenseMatrix &Jpt) const;
296 
297  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
298 
299  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
300  const double weight, DenseMatrix &A) const;
301 
302 };
303 
304 /// Untangling metric, 2D
306 {
307 protected:
308  const double eps;
310 
311 public:
312  TMOP_Metric_211(double epsilon = 1e-4) : eps(epsilon) { }
313 
314  // W = (det(J) - 1)^2 - det(J) + sqrt(det(J)^2 + eps).
315  virtual double EvalW(const DenseMatrix &Jpt) const;
316 
317  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
318 
319  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
320  const double weight, DenseMatrix &A) const;
321 };
322 
323 /// Shifted barrier form of metric 56 (area, ideal barrier metric), 2D
325 {
326 protected:
327  double &tau0;
329 
330 public:
331  /// Note that @a t0 is stored by reference
332  TMOP_Metric_252(double &t0): tau0(t0) {}
333 
334  // W = 0.5(det(J) - 1)^2 / (det(J) - tau0).
335  virtual double EvalW(const DenseMatrix &Jpt) const;
336 
337  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
338 
339  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
340  const double weight, DenseMatrix &A) const;
341 };
342 
343 /// Shape, ideal barrier metric, 3D
345 {
346 protected:
348 
349 public:
350  // W = |J| |J^-1| / 3 - 1.
351  virtual double EvalW(const DenseMatrix &Jpt) const;
352 
353  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
354 
355  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
356  const double weight, DenseMatrix &A) const;
357 };
358 
359 /// Shape, ideal barrier metric, 3D
361 {
362 protected:
364 
365 public:
366  // W = |J|^2 |J^-1|^2 / 9 - 1.
367  virtual double EvalW(const DenseMatrix &Jpt) const;
368 
369  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
370 
371  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
372  const double weight, DenseMatrix &A) const;
373 };
374 
375 /// Shape, ideal barrier metric, 3D
377 {
378 protected:
380 
381 public:
382  // W = |J|^2 / 3 * det(J)^(2/3) - 1.
383  virtual double EvalW(const DenseMatrix &Jpt) const;
384 
385  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
386 
387  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
388  const double weight, DenseMatrix &A) const;
389 };
390 
391 /// Volume metric, 3D
393 {
394 protected:
396 
397 public:
398  // W = (det(J) - 1)^2.
399  virtual double EvalW(const DenseMatrix &Jpt) const;
400 
401  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
402 
403  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
404  const double weight, DenseMatrix &A) const;
405 };
406 
407 /// Volume, ideal barrier metric, 3D
409 {
410 protected:
412 
413 public:
414  // W = 0.5( sqrt(det(J)) - 1 / sqrt(det(J)) )^2
415  // = 0.5( det(J) - 1 )^2 / det(J)
416  // = 0.5( det(J) + 1/det(J) ) - 1.
417  virtual double EvalW(const DenseMatrix &Jpt) const;
418 
419  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
420 
421  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
422  const double weight, DenseMatrix &A) const;
423 };
424 
425 /// Shape & volume, ideal barrier metric, 3D
427 {
428 protected:
430 
431 public:
432  // W = |J - J^-t|^2.
433  virtual double EvalW(const DenseMatrix &Jpt) const;
434 
435  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
436 
437  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
438  const double weight, DenseMatrix &A) const;
439 };
440 
441 /// Shifted barrier form of 3D metric 16 (volume, ideal barrier metric), 3D
443 {
444 protected:
445  double &tau0;
447 
448 public:
449  TMOP_Metric_352(double &t0): tau0(t0) {}
450 
451  // W = 0.5(det(J) - 1)^2 / (det(J) - tau0).
452  virtual double EvalW(const DenseMatrix &Jpt) const;
453 
454  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
455 
456  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
457  const double weight, DenseMatrix &A) const;
458 };
459 
460 
461 /// Base class for limiting functions to be used in class TMOP_Integrator.
462 /** This class represents a scalar function f(x, x0, d), where x and x0 are
463  positions in physical space, and d is a reference physical distance
464  associated with the point x0. */
466 {
467 public:
468  /// Returns the limiting function, f(x, x0, d).
469  virtual double Eval(const Vector &x, const Vector &x0, double d) const = 0;
470 
471  /** @brief Returns the gradient of the limiting function f(x, x0, d) with
472  respect to x. */
473  virtual void Eval_d1(const Vector &x, const Vector &x0, double dist,
474  Vector &d1) const = 0;
475 
476  /** @brief Returns the Hessian of the limiting function f(x, x0, d) with
477  respect to x. */
478  virtual void Eval_d2(const Vector &x, const Vector &x0, double dist,
479  DenseMatrix &d2) const = 0;
480 
481  /// Virtual destructor.
482  virtual ~TMOP_LimiterFunction() { }
483 };
484 
485 /// Default limiter function in TMOP_Integrator.
487 {
488 public:
489  virtual double Eval(const Vector &x, const Vector &x0, double dist) const
490  {
491  MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
492 
493  return 0.5 * x.DistanceSquaredTo(x0) / (dist * dist);
494  }
495 
496  virtual void Eval_d1(const Vector &x, const Vector &x0, double dist,
497  Vector &d1) const
498  {
499  MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
500 
501  d1.SetSize(x.Size());
502  subtract(1.0 / (dist * dist), x, x0, d1);
503  }
504 
505  virtual void Eval_d2(const Vector &x, const Vector &x0, double dist,
506  DenseMatrix &d2) const
507  {
508  MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
509 
510  d2.Diag(1.0 / (dist * dist), x.Size());
511  }
512 
514 };
515 
516 class FiniteElementCollection;
517 class FiniteElementSpace;
518 class ParFiniteElementSpace;
519 
521 {
522 protected:
523  // Owned.
526 
527 #ifdef MFEM_USE_MPI
528  // Owned.
531 #endif
532 
533 public:
534  AdaptivityEvaluator() : mesh(NULL), fes(NULL)
535  {
536 #ifdef MFEM_USE_MPI
537  pmesh = NULL;
538  pfes = NULL;
539 #endif
540  }
541  virtual ~AdaptivityEvaluator();
542 
543  /** Specifies the Mesh and FiniteElementCollection of the solution that will
544  be evaluated. The given mesh will be copied into the internal object. */
545  void SetSerialMetaInfo(const Mesh &m,
546  const FiniteElementCollection &fec, int num_comp);
547 
548 #ifdef MFEM_USE_MPI
549  /// Parallel version of SetSerialMetaInfo.
550  void SetParMetaInfo(const ParMesh &m,
551  const FiniteElementCollection &fec, int num_comp);
552 #endif
553 
554  // TODO use GridFunctions to make clear it's on the ldofs?
555  virtual void SetInitialField(const Vector &init_nodes,
556  const Vector &init_field) = 0;
557 
558  virtual void ComputeAtNewPosition(const Vector &new_nodes,
559  Vector &new_field) = 0;
560 };
561 
562 /** @brief Base class representing target-matrix construction algorithms for
563  mesh optimization via the target-matrix optimization paradigm (TMOP). */
564 /** This class is used by class TMOP_Integrator to construct the target Jacobian
565  matrices (reference-element to target-element) at quadrature points. It
566  supports a set of algorithms chosen by the #TargetType enumeration.
567 
568  New target-matrix construction algorithms can be defined by deriving new
569  classes and overriding the method ComputeElementTargets(). */
571 {
572 public:
573  /// Target-matrix construction algorithms supported by this class.
575  {
577  Ideal shape, unit size; the nodes are not used. */
579  Ideal shape, equal size/volume; the given nodes define the total target
580  volume; for each mesh element, the target volume is the average volume
581  multiplied by the volume scale, set with SetVolumeScale(). */
583  Ideal shape, given size/volume; the given nodes define the target
584  volume at all quadrature points. */
586  Given shape, given size/volume; the given nodes define the exact target
587  Jacobian matrix at all quadrature points. */
589  Full target tensor is specified at every quadrature point. */
590  };
591 
592 protected:
593  // Nodes that are used in ComputeElementTargets(), depending on target_type.
594  const GridFunction *nodes; // not owned
595  mutable double avg_volume;
596  double volume_scale;
598 
599 #ifdef MFEM_USE_MPI
600  MPI_Comm comm;
601  bool Parallel() const { return (comm != MPI_COMM_NULL); }
602 #else
603  bool Parallel() const { return false; }
604 #endif
605 
606  // should be called only if avg_volume == 0.0, i.e. avg_volume is not
607  // computed yet
608  void ComputeAvgVolume() const;
609 
610 public:
611  /// Constructor for use in serial
613  : nodes(NULL), avg_volume(), volume_scale(1.0), target_type(ttype)
614  {
615 #ifdef MFEM_USE_MPI
616  comm = MPI_COMM_NULL;
617 #endif
618  }
619 #ifdef MFEM_USE_MPI
620  /// Constructor for use in parallel
621  TargetConstructor(TargetType ttype, MPI_Comm mpicomm)
622  : nodes(NULL), avg_volume(), volume_scale(1.0), target_type(ttype),
623  comm(mpicomm) { }
624 #endif
625  virtual ~TargetConstructor() { }
626 
627  /** @brief Set the nodes to be used in the target-matrix construction.
628 
629  This method should be called every time the target nodes are updated
630  externally and recomputation of the target average volume is needed. The
631  nodes are used by all target types except IDEAL_SHAPE_UNIT_SIZE. */
632  void SetNodes(const GridFunction &n) { nodes = &n; avg_volume = 0.0; }
633 
634  /// Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
635  void SetVolumeScale(double vol_scale) { volume_scale = vol_scale; }
636 
637  /** @brief Given an element and quadrature rule, computes ref->target
638  transformation Jacobians for each quadrature point in the element.
639  The physical positions of the element's nodes are given by @a elfun. */
640  virtual void ComputeElementTargets(int e_id, const FiniteElement &fe,
641  const IntegrationRule &ir,
642  const Vector &elfun,
643  DenseTensor &Jtr) const;
644 };
645 
647 {
648 protected:
649  // Analytic target specification.
653 
654 public:
656  : TargetConstructor(ttype),
657  scalar_tspec(NULL), vector_tspec(NULL), matrix_tspec(NULL) { }
658 
659  virtual void SetAnalyticTargetSpec(Coefficient *sspec,
660  VectorCoefficient *vspec,
661  MatrixCoefficient *mspec);
662 
663  /** @brief Given an element and quadrature rule, computes ref->target
664  transformation Jacobians for each quadrature point in the element.
665  The physical positions of the element's nodes are given by @a elfun. */
666  virtual void ComputeElementTargets(int e_id, const FiniteElement &fe,
667  const IntegrationRule &ir,
668  const Vector &elfun,
669  DenseTensor &Jtr) const;
670 };
671 
672 class ParGridFunction;
673 
675 {
676 protected:
677  // Discrete target specification.
678  // Data is owned, updated by UpdateTargetSpecification.
680  // Note: do not use the Nodes of this space as they may not be on the
681  // positions corresponding to the values of tspec.
683 
684  // Evaluation of the discrete target specification on different meshes.
685  // Owned.
687 
688 public:
690  : TargetConstructor(ttype),
691  target_spec(), tspec_fes(NULL), adapt_eval(NULL) { }
692 
693  virtual ~DiscreteAdaptTC() { delete adapt_eval; }
694 
695  virtual void SetSerialDiscreteTargetSpec(GridFunction &tspec);
696 #ifdef MFEM_USE_MPI
697  virtual void SetParDiscreteTargetSpec(ParGridFunction &tspec);
698 #endif
699 
700  /** Used to update the target specification after the mesh has changed. The
701  new mesh positions are given by new_x. */
702  void UpdateTargetSpecification(const Vector &new_x);
703 
705  {
706  if (adapt_eval) { delete adapt_eval; }
707  adapt_eval = ae;
708  }
709 
710  /** @brief Given an element and quadrature rule, computes ref->target
711  transformation Jacobians for each quadrature point in the element.
712  The physical positions of the element's nodes are given by @a elfun.
713  Note that this function assumes that UpdateTargetSpecification() has
714  been called with the position vector corresponding to @a elfun. */
715  virtual void ComputeElementTargets(int e_id, const FiniteElement &fe,
716  const IntegrationRule &ir,
717  const Vector &elfun,
718  DenseTensor &Jtr) const;
719 };
720 
721 /** @brief A TMOP integrator class based on any given TMOP_QualityMetric and
722  TargetConstructor.
723 
724  Represents @f$ \int W(Jpt) dx @f$ over a target zone, where W is the
725  metric's strain energy density function, and Jpt is the Jacobian of the
726  target->physical coordinates transformation. The virtual target zone is
727  defined by the TargetConstructor. */
729 {
730 protected:
731  TMOP_QualityMetric *metric; // not owned
732  const TargetConstructor *targetC; // not owned
733 
734  // Weight Coefficient multiplying the quality metric term.
735  Coefficient *coeff1; // not owned, if NULL -> coeff1 is 1.
736  // Normalization factor for the metric term.
738 
739  // Nodes and weight Coefficient used for "limiting" the TMOP_Integrator.
740  // These are both NULL when there is no limiting.
741  // The class doesn't own nodes0 and coeff0.
744  // Limiting reference distance. Not owned.
746  // Limiting function. Owned.
748  // Normalization factor for the limiting term.
749  double lim_normal;
750 
751  // Jrt: the inverse of the ref->target Jacobian, Jrt = Jtr^{-1}.
752  // Jpr: the ref->physical transformation Jacobian, Jpr = PMatI^t DS.
753  // Jpt: the target->physical transformation Jacobian, Jpt = Jpr Jrt.
754  // P: represents dW_d(Jtp) (dim x dim).
755  // DSh: gradients of reference shape functions (dof x dim).
756  // DS: gradients of the shape functions in the target configuration,
757  // DS = DSh Jrt (dof x dim).
758  // PMatI: current coordinates of the nodes (dof x dim).
759  // PMat0: reshaped view into the local element contribution to the operator
760  // output - the result of AssembleElementVector() (dof x dim).
762 
764  double &metric_energy, double &lim_energy);
765 
766 public:
767  /** @param[in] m TMOP_QualityMetric that will be integrated (not owned).
768  @param[in] tc Target-matrix construction algorithm to use (not owned). */
770  : metric(m), targetC(tc),
771  coeff1(NULL), metric_normal(1.0),
772  nodes0(NULL), coeff0(NULL),
773  lim_dist(NULL), lim_func(NULL), lim_normal(1.0)
774  { }
775 
777 
778  /// Sets a scaling Coefficient for the quality metric term of the integrator.
779  /** With this addition, the integrator becomes
780  @f$ \int w1 W(Jpt) dx @f$.
781 
782  Note that the Coefficient is evaluated in the physical configuration and
783  not in the target configuration which may be undefined. */
784  void SetCoefficient(Coefficient &w1) { coeff1 = &w1; }
785 
786  /// Adds a limiting term to the integrator (general version).
787  /** With this addition, the integrator becomes
788  @f$ \int w1 W(Jpt) + w0 f(x, x_0, d) dx @f$,
789  where the second term measures the change with respect to the original
790  physical positions, @a n0.
791  @param[in] n0 Original mesh node coordinates.
792  @param[in] dist Limiting physical distances.
793  @param[in] w0 Coefficient scaling the limiting term.
794  @param[in] lfunc TMOP_LimiterFunction defining the limiting term f. If
795  NULL, a TMOP_QuadraticLimiter will be used. The
796  TMOP_Integrator assumes ownership of this pointer. */
797  void EnableLimiting(const GridFunction &n0, const GridFunction &dist,
798  Coefficient &w0, TMOP_LimiterFunction *lfunc = NULL);
799 
800  /** @brief Adds a limiting term to the integrator with limiting distance
801  function (@a dist in the general version of the method) equal to 1. */
802  void EnableLimiting(const GridFunction &n0,
803  Coefficient &w0, TMOP_LimiterFunction *lfunc = NULL);
804 
805  /// Update the original/reference nodes used for limiting.
806  void SetLimitingNodes(const GridFunction &n0) { nodes0 = &n0; }
807 
808  /** @brief Computes the integral of W(Jacobian(Trt)) over a target zone.
809  @param[in] el Type of FiniteElement.
810  @param[in] T Mesh element transformation.
811  @param[in] elfun Physical coordinates of the zone. */
812  virtual double GetElementEnergy(const FiniteElement &el,
814  const Vector &elfun);
815 
816  virtual void AssembleElementVector(const FiniteElement &el,
818  const Vector &elfun, Vector &elvect);
819 
820  virtual void AssembleElementGrad(const FiniteElement &el,
822  const Vector &elfun, DenseMatrix &elmat);
823 
824  /** @brief Computes the normalization factors of the metric and limiting
825  integrals using the mesh position given by @a x. */
826  void EnableNormalization(const GridFunction &x);
827 #ifdef MFEM_USE_MPI
829 #endif
830 };
831 
832 
833 /// Interpolates the @a metric's values at the nodes of @a metric_gf.
834 /** Assumes that @a metric_gf's FiniteElementSpace is initialized. */
835 void InterpolateTMOP_QualityMetric(TMOP_QualityMetric &metric,
836  const TargetConstructor &tc,
837  const Mesh &mesh, GridFunction &metric_gf);
838 
839 }
840 
841 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:232
MatrixCoefficient * matrix_tspec
Definition: tmop.hpp:652
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...
Definition: tmop.hpp:145
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:291
virtual void SetParDiscreteTargetSpec(ParGridFunction &tspec)
Definition: tmop.cpp:925
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...
Definition: tmop.cpp:277
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...
Definition: tmop.cpp:34
Shifted barrier form of metric 56 (area, ideal barrier metric), 2D.
Definition: tmop.hpp:324
Shifted barrier form of 3D metric 16 (volume, ideal barrier metric), 3D.
Definition: tmop.hpp:442
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:258
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)=0
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...
void SetSerialMetaInfo(const Mesh &m, const FiniteElementCollection &fec, int num_comp)
Definition: tmop.cpp:1006
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...
Definition: tmop.hpp:115
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
VectorCoefficient * vector_tspec
Definition: tmop.hpp:651
const TargetType target_type
Definition: tmop.hpp:597
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
DenseMatrix PMatO
Definition: tmop.hpp:761
Shape &amp; volume, ideal barrier metric, 3D.
Definition: tmop.hpp:426
const DenseMatrix * Jtr
Definition: tmop.hpp:26
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:429
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:572
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:142
DenseMatrix Jrt
Definition: tmop.hpp:761
const GridFunction * lim_dist
Definition: tmop.hpp:745
Coefficient * coeff0
Definition: tmop.hpp:743
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.
Definition: tmop.hpp:505
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:237
void EnableNormalization(const GridFunction &x)
Computes the normalization factors of the metric and limiting integrals using the mesh position given...
Definition: tmop.cpp:1366
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:407
void InterpolateTMOP_QualityMetric(TMOP_QualityMetric &metric, const TargetConstructor &tc, const Mesh &mesh, GridFunction &metric_gf)
Interpolates the metric&#39;s values at the nodes of metric_gf.
Definition: tmop.cpp:1436
Shape, ideal barrier metric, 3D.
Definition: tmop.hpp:376
TMOP_QualityMetric * metric
Definition: tmop.hpp:731
void SetAdaptivityEvaluator(AdaptivityEvaluator *ae)
Definition: tmop.hpp:704
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:254
virtual ~TMOP_LimiterFunction()
Virtual destructor.
Definition: tmop.hpp:482
Base class for limiting functions to be used in class TMOP_Integrator.
Definition: tmop.hpp:465
virtual double EvalW(const DenseMatrix &Jpt) const =0
Evaluate the strain energy density function, W = W(Jpt).
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.
Skew metric, 3D.
Definition: tmop.hpp:106
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...
Definition: tmop.cpp:699
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:626
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:163
DenseMatrix Jpr
Definition: tmop.hpp:761
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:482
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const =0
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:78
ParFiniteElementSpace * pfes
Definition: tmop.hpp:530
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:562
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Shape, ideal barrier metric, 3D.
Definition: tmop.hpp:360
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
void SetVolumeScale(double vol_scale)
Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
Definition: tmop.hpp:635
Shape, ideal barrier metric, 2D.
Definition: tmop.hpp:151
double metric_normal
Definition: tmop.hpp:737
Volume metric, 3D.
Definition: tmop.hpp:392
Area, ideal barrier metric, 2D.
Definition: tmop.hpp:288
virtual ~DiscreteAdaptTC()
Definition: tmop.hpp:693
virtual void SetSerialDiscreteTargetSpec(GridFunction &tspec)
Definition: tmop.cpp:943
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...
Definition: tmop.cpp:494
Coefficient * coeff1
Definition: tmop.hpp:735
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:127
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:112
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:186
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:273
const GridFunction * nodes
Definition: tmop.hpp:594
Aspect ratio metric, 2D.
Definition: tmop.hpp:121
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:688
Shape &amp; area, ideal barrier metric, 2D.
Definition: tmop.hpp:167
Default limiter function in TMOP_Integrator.
Definition: tmop.hpp:486
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-&gt;target transformation Jacobians for each quadratur...
Definition: tmop.cpp:967
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:365
Shifted barrier form of metric 2 (shape, ideal barrier metric), 2D.
Definition: tmop.hpp:199
TMOP_Metric_352(double &t0)
Definition: tmop.hpp:449
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-&gt;target transformation Jacobians for each quadratur...
Definition: tmop.cpp:891
Coefficient * scalar_tspec
Definition: tmop.hpp:650
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:518
void SetTargetJacobian(const DenseMatrix &_Jtr)
Specify the reference-element -&gt; target-element Jacobian matrix for the point of interest.
Definition: tmop.hpp:43
Shape, ideal barrier metric, 2D.
Definition: tmop.hpp:270
virtual ~TMOP_QualityMetric()
Definition: tmop.hpp:35
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:68
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...
Definition: tmop.cpp:608
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...
Definition: tmop.cpp:381
Abstract class for local mesh quality metrics in the target-matrix optimization paradigm (TMOP) by P...
Definition: tmop.hpp:23
virtual ~AdaptivityEvaluator()
Definition: tmop.cpp:1028
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:395
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:655
void EnableLimiting(const GridFunction &n0, const GridFunction &dist, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
Adds a limiting term to the integrator (general version).
Definition: tmop.cpp:1038
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:452
void SetCoefficient(Coefficient &w1)
Sets a scaling Coefficient for the quality metric term of the integrator.
Definition: tmop.hpp:784
Area metric, 2D.
Definition: tmop.hpp:234
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...
Definition: tmop.cpp:323
TargetConstructor(TargetType ttype)
Constructor for use in serial.
Definition: tmop.hpp:612
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...
Definition: tmop.cpp:352
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:225
Volume, ideal barrier metric, 3D.
Definition: tmop.hpp:408
DiscreteAdaptTC(TargetType ttype)
Definition: tmop.hpp:689
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:309
void SetNodes(const GridFunction &n)
Set the nodes to be used in the target-matrix construction.
Definition: tmop.hpp:632
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:596
const TargetConstructor * targetC
Definition: tmop.hpp:732
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:461
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...
Definition: tmop.cpp:411
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:203
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:108
virtual ~TargetConstructor()
Definition: tmop.hpp:625
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:424
void UpdateTargetSpecification(const Vector &new_x)
Definition: tmop.cpp:960
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:328
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
Definition: tmop.cpp:1168
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...
Definition: tmop.cpp:200
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:28
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:192
void ComputeAvgVolume() const
Definition: tmop.cpp:769
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:677
virtual double Eval(const Vector &x, const Vector &x0, double dist) const
Returns the limiting function, f(x, x0, d).
Definition: tmop.hpp:489
virtual ~TMOP_QuadraticLimiter()
Definition: tmop.hpp:513
Shape, ideal barrier metric, 2D.
Definition: tmop.hpp:218
TMOP_Metric_211(double epsilon=1e-4)
Definition: tmop.hpp:312
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
Definition: tmop.cpp:1262
Untangling metric, 2D.
Definition: tmop.hpp:305
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:97
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
const double eps
Definition: tmop.hpp:308
Area, ideal barrier metric, 2D.
Definition: tmop.hpp:251
void SetLimitingNodes(const GridFunction &n0)
Update the original/reference nodes used for limiting.
Definition: tmop.hpp:806
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:431
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:446
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition: vector.cpp:385
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-&gt;target transformation Jacobians for each quadratur...
Definition: tmop.cpp:812
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:221
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:44
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...
Definition: tmop.cpp:175
double DistanceSquaredTo(const double *p) const
Compute the square of the Euclidean distance to another vector.
Definition: vector.hpp:537
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:304
virtual void SetAnalyticTargetSpec(Coefficient *sspec, VectorCoefficient *vspec, MatrixCoefficient *mspec)
Definition: tmop.cpp:882
TMOP_Metric_252(double &t0)
Note that t0 is stored by reference.
Definition: tmop.hpp:332
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:336
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:474
AnalyticAdaptTC(TargetType ttype)
Definition: tmop.hpp:655
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:395
Skew metric, 2D.
Definition: tmop.hpp:91
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:185
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...
Definition: tmop.cpp:242
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...
Definition: tmop.cpp:634
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:403
DenseMatrix PMatI
Definition: tmop.hpp:761
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:411
Shape &amp; area metric, 2D.
Definition: tmop.hpp:183
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:233
void SetTransformation(ElementTransformation &)
The method SetTransformation() is hidden for TMOP_QualityMetrics, because it is not used...
Definition: tmop.hpp:31
virtual double Eval(const Vector &x, const Vector &x0, double d) const =0
Returns the limiting function, f(x, x0, d).
AdaptivityEvaluator * adapt_eval
Definition: tmop.hpp:686
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:347
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...
Definition: tmop.cpp:580
void ComputeNormalizationEnergies(const GridFunction &x, double &metric_energy, double &lim_energy)
Definition: tmop.cpp:1384
void Diag(double c, int n)
Creates n x n diagonal matrix with diagonal elements c.
Definition: densemat.cpp:1346
void SetParMetaInfo(const ParMesh &m, const FiniteElementCollection &fec, int num_comp)
Parallel version of SetSerialMetaInfo.
Definition: tmop.cpp:1017
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:373
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
Computes the integral of W(Jacobian(Trt)) over a target zone.
Definition: tmop.cpp:1074
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:22
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.
Aspect ratio metric, 3D.
Definition: tmop.hpp:136
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:647
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:154
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...
Definition: tmop.cpp:527
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:170
Abstract class for hyperelastic models.
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:128
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:267
TargetConstructor(TargetType ttype, MPI_Comm mpicomm)
Constructor for use in parallel.
Definition: tmop.hpp:621
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...
Definition: tmop.cpp:663
double epsilon(const Vector &x)
Definition: maxwell.cpp:90
TMOP_Integrator(TMOP_QualityMetric *m, TargetConstructor *tc)
Definition: tmop.hpp:769
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:512
Vector data type.
Definition: vector.hpp:48
DenseMatrix DSh
Definition: tmop.hpp:761
DenseMatrix DS
Definition: tmop.hpp:761
TMOP_Metric_022(double &t0)
Definition: tmop.hpp:206
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:379
const GridFunction * nodes0
Definition: tmop.hpp:742
const FiniteElementSpace * tspec_fes
Definition: tmop.hpp:682
TargetType
Target-matrix construction algorithms supported by this class.
Definition: tmop.hpp:574
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:344
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...
Definition: tmop.cpp:466
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:315
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.
Definition: tmop.hpp:496
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...
Definition: tmop.hpp:130
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Base class representing target-matrix construction algorithms for mesh optimization via the target-ma...
Definition: tmop.hpp:570
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:169
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:725
Class for parallel meshes.
Definition: pmesh.hpp:32
TMOP_LimiterFunction * lim_func
Definition: tmop.hpp:747
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:618
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:602
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:363
DenseMatrix Jpt
Definition: tmop.hpp:761
void ParEnableNormalization(const ParGridFunction &x)
Definition: tmop.cpp:1374
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:733
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...
Definition: tmop.hpp:100
DenseMatrix P
Definition: tmop.hpp:761
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...
Definition: tmop.cpp:745
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field)=0
bool Parallel() const
Definition: tmop.hpp:601
Shape, ideal barrier metric, 3D.
Definition: tmop.hpp:344
Metric without a type, 2D.
Definition: tmop.hpp:75
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:725
FiniteElementSpace * fes
Definition: tmop.hpp:525
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...
Definition: tmop.cpp:440
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:728