MFEM  v4.4.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-2022, 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  virtual 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  /** @brief Return the metric ID.
73  */
74  virtual int Id() const { return 0; }
75 };
76 
77 /// Abstract class used to define combination of metrics with constant coefficients.
79 {
80 protected:
83 
84 public:
85  virtual void AddQualityMetric(TMOP_QualityMetric *tq, double wt = 1.0)
86  {
87  tmop_q_arr.Append(tq);
88  wt_arr.Append(wt);
89  }
90 
91  virtual void SetTargetJacobian(const DenseMatrix &Jtr_)
92  {
93  for (int i = 0; i < tmop_q_arr.Size(); i++)
94  {
95  tmop_q_arr[i]->SetTargetJacobian(Jtr_);
96  }
97  }
98 
99  virtual double EvalW(const DenseMatrix &Jpt) const;
100 
101  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
102 
103  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
104  const double weight, DenseMatrix &A) const;
105 };
106 
107 /// 2D non-barrier metric without a type.
109 {
110 protected:
112 
113 public:
114  // W = |J|^2.
115  virtual double EvalW(const DenseMatrix &Jpt) const;
116 
117  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
118 
119  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
120  const double weight, DenseMatrix &A) const;
121 
122  virtual int Id() const { return 1; }
123 };
124 
125 /// 2D non-barrier Skew metric.
127 {
128 public:
129  // W = 0.5 (1 - cos(angle_Jpr - angle_Jtr)).
130  virtual double EvalW(const DenseMatrix &Jpt) const;
131 
132  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
133  { MFEM_ABORT("Not implemented"); }
134 
135  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
136  const double weight, DenseMatrix &A) const
137  { MFEM_ABORT("Not implemented"); }
138 };
139 
140 /// 3D non-barrier Skew metric.
142 {
143 public:
144  // W = 1/6 (3 - sum_i cos(angle_Jpr_i - angle_Jtr_i)), i = 1..3.
145  virtual double EvalW(const DenseMatrix &Jpt) const;
146 
147  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
148  { MFEM_ABORT("Not implemented"); }
149 
150  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
151  const double weight, DenseMatrix &A) const
152  { MFEM_ABORT("Not implemented"); }
153 };
154 
155 /// 2D non-barrier Aspect ratio metric.
157 {
158 public:
159  // W = 0.5 (ar_Jpr/ar_Jtr + ar_Jtr/ar_Jpr) - 1.
160  virtual double EvalW(const DenseMatrix &Jpt) const;
161 
162  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
163  { MFEM_ABORT("Not implemented"); }
164 
165  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
166  const double weight, DenseMatrix &A) const
167  { MFEM_ABORT("Not implemented"); }
168 };
169 
170 /// 3D non-barrier Aspect ratio metric.
172 {
173 public:
174  // W = 1/3 sum [0.5 (ar_Jpr_i/ar_Jtr_i + ar_Jtr_i/ar_Jpr_i) - 1], i = 1..3.
175  virtual double EvalW(const DenseMatrix &Jpt) const;
176 
177  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
178  { MFEM_ABORT("Not implemented"); }
179 
180  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
181  const double weight, DenseMatrix &A) const
182  { MFEM_ABORT("Not implemented"); }
183 };
184 
185 /// 2D barrier shape (S) metric (polyconvex).
187 {
188 protected:
190 
191 public:
192  // W = 0.5|J|^2 / det(J) - 1.
193  virtual double EvalW(const DenseMatrix &Jpt) const;
194 
195  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
196 
197  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
198  const double weight, DenseMatrix &A) const;
199 
200  virtual int Id() const { return 2; }
201 };
202 
203 /// 2D barrier Shape+Size (VS) metric (not polyconvex).
205 {
206 protected:
208 
209 public:
210  // W = |J - J^-t|^2.
211  virtual double EvalW(const DenseMatrix &Jpt) const;
212 
213  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
214 
215  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
216  const double weight, DenseMatrix &A) const;
217 
218  virtual int Id() const { return 7; }
219 };
220 
221 /// 2D barrier Shape+Size (VS) metric (not polyconvex).
223 {
224 protected:
226 
227 public:
228  // W = det(J) * |J - J^-t|^2.
229  virtual double EvalW(const DenseMatrix &Jpt) const;
230 
231  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
232 
233  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
234  const double weight, DenseMatrix &A) const;
235 };
236 
237 /// 2D non-barrier Shape+Size+Orientation (VOS) metric (polyconvex).
239 {
240 public:
241  // W = |T-I|^2.
242  virtual double EvalW(const DenseMatrix &Jpt) const;
243 
244  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
245  { MFEM_ABORT("Not implemented"); }
246 
247  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
248  const double weight, DenseMatrix &A) const
249  { MFEM_ABORT("Not implemented"); }
250 };
251 
252 /// 2D Shifted barrier form of shape metric (mu_2).
254 {
255 protected:
256  double &min_detT;
258 
259 public:
260  TMOP_Metric_022(double &t0): min_detT(t0) {}
261 
262  // W = 0.5(|J|^2 - 2det(J)) / (det(J) - tau0).
263  virtual double EvalW(const DenseMatrix &Jpt) const;
264 
265  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
266 
267  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
268  const double weight, DenseMatrix &A) const;
269 };
270 
271 /// 2D barrier (not a shape) metric (polyconvex).
273 {
274 protected:
276 
277 public:
278  // W = 0.5|J^t J|^2 / det(J)^2 - 1.
279  virtual double EvalW(const DenseMatrix &Jpt) const;
280 
281  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
282 
283  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
284  const double weight, DenseMatrix &A) const;
285 };
286 
287 /// 2D non-barrier size (V) metric (not polyconvex).
289 {
290 protected:
292 
293 public:
294  // W = (det(J) - 1)^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 /// 2D barrier size (V) metric (polyconvex).
306 {
307 protected:
309 
310 public:
311  // W = 0.5( sqrt(det(J)) - 1 / sqrt(det(J)) )^2
312  // = 0.5( det(J) - 1 )^2 / det(J)
313  // = 0.5( det(J) + 1/det(J) ) - 1.
314  virtual double EvalW(const DenseMatrix &Jpt) const;
315 
316  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
317 
318  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
319  const double weight, DenseMatrix &A) const;
320 };
321 
322 /// 2D barrier shape (S) metric (not polyconvex).
324 {
325 protected:
327 
328 public:
329  // W = |J^t J|^2 / det(J)^2 - 2|J|^2 / det(J) + 2
330  // = I1b (I1b - 2).
331  virtual double EvalW(const DenseMatrix &Jpt) const;
332 
333  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
334 
335  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
336  const double weight, DenseMatrix &A) const;
337 };
338 
339 /// 2D barrier size (V) metric (polyconvex).
341 {
342 protected:
344 
345 public:
346  // W = 0.5(det(J) - 1 / det(J))^2.
347  virtual double EvalW(const DenseMatrix &Jpt) const;
348 
349  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
350 
351  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
352  const double weight, DenseMatrix &A) const;
353 
354  virtual int Id() const { return 77; }
355 };
356 
357 /// 2D barrier Shape+Size (VS) metric (polyconvex).
359 {
360 protected:
362  double gamma;
364 
365 public:
366  TMOP_Metric_080(double gamma_) : gamma(gamma_),
369  {
370  // (1-gamma) mu_2 + gamma mu_77
371  AddQualityMetric(sh_metric, 1.-gamma_);
372  AddQualityMetric(sz_metric, gamma_);
373  }
374  virtual int Id() const { return 80; }
375  double GetGamma() const { return gamma; }
376 
377  virtual ~TMOP_Metric_080() { delete sh_metric; delete sz_metric; }
378 };
379 
380 /// 2D barrier Shape+Orientation (OS) metric (polyconvex).
382 {
383 public:
384  // W = |T-T'|^2, where T'= |T|*I/sqrt(2).
385  virtual double EvalW(const DenseMatrix &Jpt) const;
386 
387  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
388  { MFEM_ABORT("Not implemented"); }
389 
390  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
391  const double weight, DenseMatrix &A) const
392  { MFEM_ABORT("Not implemented"); }
393 };
394 
395 /// 2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
397 {
398 public:
399  // W = 1/tau |T-I|^2.
400  virtual double EvalW(const DenseMatrix &Jpt) const;
401 
402  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
403  { MFEM_ABORT("Not implemented"); }
404 
405  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
406  const double weight, DenseMatrix &A) const
407  { MFEM_ABORT("Not implemented"); }
408 };
409 
410 /// 2D untangling metric.
412 {
413 protected:
414  const double eps;
416 
417 public:
418  TMOP_Metric_211(double epsilon = 1e-4) : eps(epsilon) { }
419 
420  // W = (det(J) - 1)^2 - det(J) + sqrt(det(J)^2 + eps).
421  virtual double EvalW(const DenseMatrix &Jpt) const;
422 
423  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
424 
425  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
426  const double weight, DenseMatrix &A) const;
427 };
428 
429 /// Shifted barrier form of metric 56 (area, ideal barrier metric), 2D
431 {
432 protected:
433  double &tau0;
435 
436 public:
437  /// Note that @a t0 is stored by reference
438  TMOP_Metric_252(double &t0): tau0(t0) {}
439 
440  // W = 0.5(det(J) - 1)^2 / (det(J) - tau0).
441  virtual double EvalW(const DenseMatrix &Jpt) const;
442 
443  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
444 
445  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
446  const double weight, DenseMatrix &A) const;
447 };
448 
449 /// 3D barrier Shape (S) metric.
451 {
452 protected:
454 
455 public:
456  // W = |J| |J^-1| / 3 - 1.
457  virtual double EvalW(const DenseMatrix &Jpt) const;
458 
459  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
460 
461  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
462  const double weight, DenseMatrix &A) const;
463 };
464 
465 /// 3D barrier Shape (S) metric.
467 {
468 protected:
470 
471 public:
472  // W = |J|^2 |J^-1|^2 / 9 - 1.
473  virtual double EvalW(const DenseMatrix &Jpt) const;
474 
475  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
476 
477  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
478  const double weight, DenseMatrix &A) const;
479 
480  virtual int Id() const { return 302; }
481 };
482 
483 /// 3D barrier Shape (S) metric.
485 {
486 protected:
488 
489 public:
490  // W = |J|^2 / 3 * det(J)^(-2/3) - 1.
491  virtual double EvalW(const DenseMatrix &Jpt) const;
492 
493  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
494 
495  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
496  const double weight, DenseMatrix &A) const;
497 
498  virtual int Id() const { return 303; }
499 };
500 
501 /// 3D Size (V) untangling metric.
503 {
504 protected:
505  const double eps;
507 
508 public:
509  TMOP_Metric_311(double epsilon = 1e-4) : eps(epsilon) { }
510 
511  // W = (det(J) - 1)^2 - det(J) + (det(J)^2 + eps)^(1/2).
512  virtual double EvalW(const DenseMatrix &Jpt) const;
513 
514  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
515 
516  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
517  const double weight, DenseMatrix &A) const;
518 };
519 
520 /// 3D Shape (S) metric, untangling version of 303.
522 {
523 protected:
524  double &min_detT;
526 
527 public:
528  TMOP_Metric_313(double &mindet) : min_detT(mindet) { }
529 
530  // W = 1/3 |J|^2 / [det(J)-tau0]^(-2/3).
531  virtual double EvalW(const DenseMatrix &Jpt) const;
532 
533  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
534 
535  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
536  const double weight, DenseMatrix &A) const;
537 
538  virtual int Id() const { return 313; }
539 };
540 
541 /// 3D non-barrier Size (V) metric.
543 {
544 protected:
546 
547 public:
548  // W = (det(J) - 1)^2.
549  virtual double EvalW(const DenseMatrix &Jpt) const;
550 
551  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
552 
553  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
554  const double weight, DenseMatrix &A) const;
555 
556  virtual int Id() const { return 315; }
557 };
558 
559 /// 3D barrier Size (V) metric.
561 {
562 protected:
564 
565 public:
566  // W = 0.5( sqrt(det(J)) - 1 / sqrt(det(J)) )^2
567  // = 0.5( det(J) - 1 )^2 / det(J)
568  // = 0.5( det(J) + 1/det(J) ) - 1.
569  virtual double EvalW(const DenseMatrix &Jpt) const;
570 
571  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
572 
573  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
574  const double weight, DenseMatrix &A) const;
575 };
576 
577 /// 3D barrier Shape+Size (VS) metric.
579 {
580 protected:
582 
583 public:
584  // W = |J - J^-t|^2.
585  virtual double EvalW(const DenseMatrix &Jpt) const;
586 
587  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
588 
589  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
590  const double weight, DenseMatrix &A) const;
591 
592  virtual int Id() const { return 321; }
593 };
594 
595 /// 3D barrier Shape+Size (VS) metric (polyconvex).
597 {
598 protected:
600  double gamma;
602 
603 public:
604  TMOP_Metric_328(double gamma_) : gamma(gamma_),
607  {
608  // (1-gamma) mu_301 + gamma mu_316
609  AddQualityMetric(sh_metric, 1.-gamma_);
610  AddQualityMetric(sz_metric, gamma_);
611  }
612 
613  virtual ~TMOP_Metric_328() { delete sh_metric; delete sz_metric; }
614 };
615 
616 /// 3D barrier Shape+Size (VS) metric (polyconvex).
618 {
619 protected:
620  double gamma;
622 
623 public:
624  TMOP_Metric_332(double gamma_) : gamma(gamma_),
627  {
628  // (1-gamma) mu_302 + gamma mu_315
629  AddQualityMetric(sh_metric, 1.-gamma_);
630  AddQualityMetric(sz_metric, gamma_);
631  }
632 
633  virtual int Id() const { return 332; }
634  double GetGamma() const { return gamma; }
635 
636  virtual ~TMOP_Metric_332() { delete sh_metric; delete sz_metric; }
637 };
638 
639 /// 3D barrier Shape+Size (VS) metric (polyconvex).
641 {
642 protected:
644  double gamma;
646 
647 public:
648  TMOP_Metric_333(double gamma_) : gamma(gamma_),
651  {
652  // (1-gamma) mu_302 + gamma mu_316
653  AddQualityMetric(sh_metric, 1.-gamma_);
654  AddQualityMetric(sz_metric, gamma_);
655  }
656 
657  virtual ~TMOP_Metric_333() { delete sh_metric; delete sz_metric; }
658 };
659 
660 /// 3D barrier Shape+Size (VS) metric (polyconvex).
662 {
663 protected:
665  double gamma;
667 
668 public:
669  TMOP_Metric_334(double gamma_) : gamma(gamma_),
672  {
673  // (1-gamma) mu_303 + gamma mu_316
674  AddQualityMetric(sh_metric, 1.-gamma_);
675  AddQualityMetric(sz_metric, gamma_);
676  }
677 
678  virtual ~TMOP_Metric_334() { delete sh_metric; delete sz_metric; }
679 };
680 
681 /// Shifted barrier form of 3D metric 16 (volume, ideal barrier metric), 3D
683 {
684 protected:
685  double &tau0;
687 
688 public:
689  TMOP_Metric_352(double &t0): tau0(t0) {}
690 
691  // W = 0.5(det(J) - 1)^2 / (det(J) - tau0).
692  virtual double EvalW(const DenseMatrix &Jpt) const;
693 
694  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
695 
696  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
697  const double weight, DenseMatrix &A) const;
698 };
699 
700 /// A-metrics
701 /// 2D barrier Shape (S) metric (polyconvex).
703 {
704 protected:
706 
707 public:
708  // (1/4 alpha) | A - (adj A)^t W^t W / omega |^2
709  virtual double EvalW(const DenseMatrix &Jpt) const;
710 
711  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
712  { MFEM_ABORT("Not implemented"); }
713 
714  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
715  const double weight, DenseMatrix &A) const
716  { MFEM_ABORT("Not implemented"); }
717 };
718 
719 /// 2D barrier Size (V) metric (polyconvex).
721 {
722 protected:
724 
725 public:
726  // 0.5 * ( sqrt(alpha/omega) - sqrt(omega/alpha) )^2
727  virtual double EvalW(const DenseMatrix &Jpt) const;
728 
729  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
730  { MFEM_ABORT("Not implemented"); }
731 
732  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
733  const double weight, DenseMatrix &A) const
734  { MFEM_ABORT("Not implemented"); }
735 };
736 
737 /// 2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
739 {
740 protected:
742 
743 public:
744  // (1/alpha) | A - W |^2
745  virtual double EvalW(const DenseMatrix &Jpt) const;
746 
747  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
748  { MFEM_ABORT("Not implemented"); }
749 
750  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
751  const double weight, DenseMatrix &A) const
752  { MFEM_ABORT("Not implemented"); }
753 };
754 
755 /// 2D barrier Shape+Orientation (OS) metric (polyconvex).
757 {
758 protected:
760 
761 public:
762  // (1/2 alpha) | A - (|A|/|W|) W |^2
763  virtual double EvalW(const DenseMatrix &Jpt) const;
764 
765  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
766  { MFEM_ABORT("Not implemented"); }
767 
768  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
769  const double weight, DenseMatrix &A) const
770  { MFEM_ABORT("Not implemented"); }
771 };
772 
773 /// 2D barrier Shape+Size (VS) metric (polyconvex).
775 {
776 protected:
778  double gamma;
780 
781 public:
782  TMOP_AMetric_126(double gamma_) : gamma(gamma_),
785  {
786  // (1-gamma) nu_11 + gamma nu_14
787  AddQualityMetric(sh_metric, 1.-gamma_);
788  AddQualityMetric(sz_metric, gamma_);
789  }
790 
791  virtual ~TMOP_AMetric_126() { delete sh_metric; delete sz_metric; }
792 };
793 
794 /// Base class for limiting functions to be used in class TMOP_Integrator.
795 /** This class represents a scalar function f(x, x0, d), where x and x0 are
796  positions in physical space, and d is a reference physical distance
797  associated with the point x0. */
799 {
800 public:
801  /// Returns the limiting function, f(x, x0, d).
802  virtual double Eval(const Vector &x, const Vector &x0, double d) const = 0;
803 
804  /** @brief Returns the gradient of the limiting function f(x, x0, d) with
805  respect to x. */
806  virtual void Eval_d1(const Vector &x, const Vector &x0, double dist,
807  Vector &d1) const = 0;
808 
809  /** @brief Returns the Hessian of the limiting function f(x, x0, d) with
810  respect to x. */
811  virtual void Eval_d2(const Vector &x, const Vector &x0, double dist,
812  DenseMatrix &d2) const = 0;
813 
814  /// Virtual destructor.
815  virtual ~TMOP_LimiterFunction() { }
816 };
817 
818 /// Default limiter function in TMOP_Integrator.
820 {
821 public:
822  virtual double Eval(const Vector &x, const Vector &x0, double dist) const
823  {
824  MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
825 
826  return 0.5 * x.DistanceSquaredTo(x0) / (dist * dist);
827  }
828 
829  virtual void Eval_d1(const Vector &x, const Vector &x0, double dist,
830  Vector &d1) const
831  {
832  MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
833 
834  d1.SetSize(x.Size());
835  subtract(1.0 / (dist * dist), x, x0, d1);
836  }
837 
838  virtual void Eval_d2(const Vector &x, const Vector &x0, double dist,
839  DenseMatrix &d2) const
840  {
841  MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
842 
843  d2.Diag(1.0 / (dist * dist), x.Size());
844  }
845 
847 };
848 
849 class FiniteElementCollection;
850 class FiniteElementSpace;
851 class ParFiniteElementSpace;
852 
854 {
855 protected:
856  // Owned.
859 
860 #ifdef MFEM_USE_MPI
861  // Owned.
864 #endif
865 
866  int dim, ncomp;
867 
868 public:
869  AdaptivityEvaluator() : mesh(NULL), fes(NULL)
870  {
871 #ifdef MFEM_USE_MPI
872  pmesh = NULL;
873  pfes = NULL;
874 #endif
875  }
876  virtual ~AdaptivityEvaluator();
877 
878  /** Specifies the Mesh and FiniteElementCollection of the solution that will
879  be evaluated. The given mesh will be copied into the internal object. */
880  void SetSerialMetaInfo(const Mesh &m,
881  const FiniteElementCollection &fec, int num_comp);
882 
883 #ifdef MFEM_USE_MPI
884  /// Parallel version of SetSerialMetaInfo.
885  void SetParMetaInfo(const ParMesh &m,
886  const FiniteElementCollection &fec, int num_comp);
887 #endif
888 
889  // TODO use GridFunctions to make clear it's on the ldofs?
890  virtual void SetInitialField(const Vector &init_nodes,
891  const Vector &init_field) = 0;
892 
893  virtual void ComputeAtNewPosition(const Vector &new_nodes,
894  Vector &new_field) = 0;
895 
896  void ClearGeometricFactors();
897 };
898 
899 /** @brief Base class representing target-matrix construction algorithms for
900  mesh optimization via the target-matrix optimization paradigm (TMOP). */
901 /** This class is used by class TMOP_Integrator to construct the target Jacobian
902  matrices (reference-element to target-element) at quadrature points. It
903  supports a set of algorithms chosen by the #TargetType enumeration.
904 
905  New target-matrix construction algorithms can be defined by deriving new
906  classes and overriding the methods ComputeElementTargets() and
907  ContainsVolumeInfo(). */
909 {
910 public:
911  /// Target-matrix construction algorithms supported by this class.
913  {
915  Ideal shape, unit size; the nodes are not used. */
917  Ideal shape, equal size/volume; the given nodes define the total target
918  volume; for each mesh element, the target volume is the average volume
919  multiplied by the volume scale, set with SetVolumeScale(). */
921  Ideal shape, given size/volume; the given nodes define the target
922  volume at all quadrature points. */
924  Given shape, given size/volume; the given nodes define the exact target
925  Jacobian matrix at all quadrature points. */
927  Full target tensor is specified at every quadrature point. */
928  };
929 
930 protected:
931  // Nodes that are used in ComputeElementTargets(), depending on target_type.
932  const GridFunction *nodes; // not owned
933  mutable double avg_volume;
934  double volume_scale;
936  bool uses_phys_coords; // see UsesPhysicalCoordinates()
937 
938 #ifdef MFEM_USE_MPI
939  MPI_Comm comm;
940 #endif
941 
942  // should be called only if avg_volume == 0.0, i.e. avg_volume is not
943  // computed yet
944  void ComputeAvgVolume() const;
945 
946  template<int DIM>
948  const IntegrationRule &ir,
949  const Vector &xe,
950  DenseTensor &Jtr) const;
951 
952  // CPU fallback that uses ComputeElementTargets()
954  const IntegrationRule &ir,
955  const Vector &xe,
956  DenseTensor &Jtr) const;
957 
958 public:
959  /// Constructor for use in serial
961  : nodes(NULL), avg_volume(), volume_scale(1.0), target_type(ttype),
962  uses_phys_coords(false)
963  {
964 #ifdef MFEM_USE_MPI
965  comm = MPI_COMM_NULL;
966 #endif
967  }
968 #ifdef MFEM_USE_MPI
969  /// Constructor for use in parallel
970  TargetConstructor(TargetType ttype, MPI_Comm mpicomm)
971  : nodes(NULL), avg_volume(), volume_scale(1.0), target_type(ttype),
972  uses_phys_coords(false), comm(mpicomm) { }
973 #endif
974  virtual ~TargetConstructor() { }
975 
976 #ifdef MFEM_USE_MPI
977  bool Parallel() const { return (comm != MPI_COMM_NULL); }
978  MPI_Comm GetComm() const { return comm; }
979 #else
980  bool Parallel() const { return false; }
981 #endif
982 
983  /** @brief Set the nodes to be used in the target-matrix construction.
984 
985  This method should be called every time the target nodes are updated
986  externally and recomputation of the target average volume is needed. The
987  nodes are used by all target types except IDEAL_SHAPE_UNIT_SIZE. */
988  void SetNodes(const GridFunction &n) { nodes = &n; avg_volume = 0.0; }
989 
990  /** @brief Get the nodes to be used in the target-matrix construction. */
991  const GridFunction *GetNodes() const { return nodes; }
992 
993  /// Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
994  void SetVolumeScale(double vol_scale) { volume_scale = vol_scale; }
995 
997 
998  /** @brief Return true if the methods ComputeElementTargets(),
999  ComputeAllElementTargets(), and ComputeElementTargetsGradient() use the
1000  physical node coordinates provided by the parameters 'elfun', or 'xe'. */
1002 
1003  /// Checks if the target matrices contain non-trivial size specification.
1004  virtual bool ContainsVolumeInfo() const;
1005 
1006  /** @brief Given an element and quadrature rule, computes ref->target
1007  transformation Jacobians for each quadrature point in the element.
1008  The physical positions of the element's nodes are given by @a elfun. */
1009  virtual void ComputeElementTargets(int e_id, const FiniteElement &fe,
1010  const IntegrationRule &ir,
1011  const Vector &elfun,
1012  DenseTensor &Jtr) const;
1013 
1014  /** @brief Computes reference-to-target transformation Jacobians for all
1015  quadrature points in all elements.
1016 
1017  @param[in] fes The nodal FE space
1018  @param[in] ir The quadrature rule to use for all elements
1019  @param[in] xe E-vector with the current physical coordinates/positions;
1020  this parameter is used only when needed by the target
1021  constructor, see UsesPhysicalCoordinates()
1022  @param[out] Jtr The computed ref->target Jacobian matrices. */
1023  virtual void ComputeAllElementTargets(const FiniteElementSpace &fes,
1024  const IntegrationRule &ir,
1025  const Vector &xe,
1026  DenseTensor &Jtr) const;
1027 
1028  virtual void ComputeElementTargetsGradient(const IntegrationRule &ir,
1029  const Vector &elfun,
1031  DenseTensor &dJtr) const;
1032 };
1033 
1035 {
1036 public:
1037  explicit TMOPMatrixCoefficient(int dim) : MatrixCoefficient(dim, dim) { }
1038 
1039  /** @brief Evaluate the derivative of the matrix coefficient with respect to
1040  @a comp in the element described by @a T at the point @a ip, storing the
1041  result in @a K. */
1042  virtual void EvalGrad(DenseMatrix &K, ElementTransformation &T,
1043  const IntegrationPoint &ip, int comp) = 0;
1044 
1046 };
1047 
1049 {
1050 protected:
1051  // Analytic target specification.
1055 
1056 public:
1058  : TargetConstructor(ttype),
1059  scalar_tspec(NULL), vector_tspec(NULL), matrix_tspec(NULL)
1060  { uses_phys_coords = true; }
1061 
1062  virtual void SetAnalyticTargetSpec(Coefficient *sspec,
1063  VectorCoefficient *vspec,
1064  TMOPMatrixCoefficient *mspec);
1065 
1066  /** @brief Given an element and quadrature rule, computes ref->target
1067  transformation Jacobians for each quadrature point in the element.
1068  The physical positions of the element's nodes are given by @a elfun. */
1069  virtual void ComputeElementTargets(int e_id, const FiniteElement &fe,
1070  const IntegrationRule &ir,
1071  const Vector &elfun,
1072  DenseTensor &Jtr) const;
1073 
1074  virtual void ComputeAllElementTargets(const FiniteElementSpace &fes,
1075  const IntegrationRule &ir,
1076  const Vector &xe,
1077  DenseTensor &Jtr) const;
1078 
1079  virtual void ComputeElementTargetsGradient(const IntegrationRule &ir,
1080  const Vector &elfun,
1082  DenseTensor &dJtr) const;
1083 };
1084 
1085 #ifdef MFEM_USE_MPI
1086 class ParGridFunction;
1087 #endif
1088 
1090 {
1091 protected:
1092  // Discrete target specification.
1093  // Data is owned, updated by UpdateTargetSpecification.
1095  Vector tspec; //eta(x)
1097  Vector tspec_pert1h; //eta(x+h)
1098  Vector tspec_pert2h; //eta(x+2*h)
1099  Vector tspec_pertmix; //eta(x+h,y+h)
1100  // The order inside these perturbation vectors (e.g. in 2D) is
1101  // eta1(x+h,y), eta2(x+h,y) ... etan(x+h,y), eta1(x,y+h), eta2(x,y+h) ...
1102  // same for tspec_pert2h and tspec_pertmix.
1103 
1104  // DenseMatrix to hold target_spec values for the (children of the)
1105  // element being refined to consider for h-refinement.
1107  // Vector to hold the target_spec values for the coarse version of the
1108  // current mesh. Used for derefinement decision with hr-adaptivity.
1110 
1111  // Components of Target Jacobian at each quadrature point of an element. This
1112  // is required for computation of the derivative using chain rule.
1114 
1115  // Note: do not use the Nodes of this space as they may not be on the
1116  // positions corresponding to the values of tspec.
1118  FiniteElementSpace *coarse_tspec_fesv; //not owned, derefinement FESpace
1119  GridFunction *tspec_gf; //owned, uses tspec and tspec_fes
1120  // discrete adaptivity
1121 #ifdef MFEM_USE_MPI
1122  ParFiniteElementSpace *ptspec_fesv; //owned, needed for derefinement to
1123  // get update operator.
1124  ParGridFunction *tspec_pgf; // similar to tspec_gf
1125 #endif
1126 
1127  int amr_el;
1129 
1130  // These flags can be used by outside functions to avoid recomputing the
1131  // tspec and tspec_perth fields again on the same mesh.
1133 
1134  // Evaluation of the discrete target specification on different meshes.
1135  // Owned.
1137 
1138  void SetDiscreteTargetBase(const GridFunction &tspec_);
1139  void SetTspecAtIndex(int idx, const GridFunction &tspec_);
1140  void FinalizeSerialDiscreteTargetSpec(const GridFunction &tspec_);
1141 #ifdef MFEM_USE_MPI
1142  void SetTspecAtIndex(int idx, const ParGridFunction &tspec_);
1143  void FinalizeParDiscreteTargetSpec(const ParGridFunction &tspec_);
1144 #endif
1145 
1146 public:
1148  : TargetConstructor(ttype),
1149  ncomp(0),
1150  sizeidx(-1), skewidx(-1), aspectratioidx(-1), orientationidx(-1),
1153  tspec_fesv(NULL), coarse_tspec_fesv(NULL), tspec_gf(NULL),
1154 #ifdef MFEM_USE_MPI
1155  ptspec_fesv(NULL), tspec_pgf(NULL),
1156 #endif
1157  amr_el(-1), lim_min_size(-0.1),
1158  good_tspec(false), good_tspec_grad(false), good_tspec_hess(false),
1159  adapt_eval(NULL) { }
1160 
1161  virtual ~DiscreteAdaptTC();
1162 
1163  /** @name Target specification methods.
1164  The following methods are used to specify geometric parameters of the
1165  targets when these parameters are given by discrete FE functions.
1166  Note that every GridFunction given to the Set methods must use a
1167  H1_FECollection of the same order. The number of components must
1168  correspond to the type of geometric parameter and dimension.
1169 
1170  @param[in] tspec_ Input values of a geometric parameter. Note that
1171  the methods in this class support only functions that
1172  use H1_FECollection collection of the same order. */
1173  ///@{
1174  virtual void SetSerialDiscreteTargetSpec(const GridFunction &tspec_);
1175  virtual void SetSerialDiscreteTargetSize(const GridFunction &tspec_);
1176  virtual void SetSerialDiscreteTargetSkew(const GridFunction &tspec_);
1177  virtual void SetSerialDiscreteTargetAspectRatio(const GridFunction &tspec_);
1178  virtual void SetSerialDiscreteTargetOrientation(const GridFunction &tspec_);
1179 #ifdef MFEM_USE_MPI
1180  virtual void SetParDiscreteTargetSpec(const ParGridFunction &tspec_);
1181  virtual void SetParDiscreteTargetSize(const ParGridFunction &tspec_);
1182  virtual void SetParDiscreteTargetSkew(const ParGridFunction &tspec_);
1183  virtual void SetParDiscreteTargetAspectRatio(const ParGridFunction &tspec_);
1184  virtual void SetParDiscreteTargetOrientation(const ParGridFunction &tspec_);
1185 #endif
1186  ///@}
1187 
1188  /// Used in combination with the Update methods to avoid extra computations.
1190  { good_tspec = good_tspec_grad = good_tspec_hess = false; }
1191 
1192  /// Get one of the discrete fields from tspec.
1193  void GetDiscreteTargetSpec(GridFunction &tspec_, int idx);
1194  /// Get the FESpace associated with tspec.
1196  /// Get the entire tspec.
1198  /// Update all discrete fields based on tspec and update for AMR
1200 
1201 #ifdef MFEM_USE_MPI
1204 #endif
1205 
1206  /** Used to update the target specification after the mesh has changed. The
1207  new mesh positions are given by new_x. If @a use_flags is true, repeated
1208  calls won't do anything until ResetUpdateFlags() is called. */
1209  void UpdateTargetSpecification(const Vector &new_x, bool use_flag = false);
1210 
1211  void UpdateTargetSpecification(Vector &new_x, Vector &IntData);
1212 
1215  int nodenum, int idir,
1216  const Vector &IntData);
1217 
1219 
1220  /** Used for finite-difference based computations. Computes the target
1221  specifications after a mesh perturbation in x or y direction.
1222  If @a use_flags is true, repeated calls won't do anything until
1223  ResetUpdateFlags() is called. */
1224  void UpdateGradientTargetSpecification(const Vector &x, double dx,
1225  bool use_flag = false);
1226  /** Used for finite-difference based computations. Computes the target
1227  specifications after two mesh perturbations in x and/or y direction.
1228  If @a use_flags is true, repeated calls won't do anything until
1229  ResetUpdateFlags() is called. */
1230  void UpdateHessianTargetSpecification(const Vector &x, double dx,
1231  bool use_flag = false);
1232 
1234  {
1235  if (adapt_eval) { delete adapt_eval; }
1236  adapt_eval = ae;
1237  }
1238 
1239  const Vector &GetTspecPert1H() { return tspec_pert1h; }
1240  const Vector &GetTspecPert2H() { return tspec_pert2h; }
1242 
1243  /** @brief Given an element and quadrature rule, computes ref->target
1244  transformation Jacobians for each quadrature point in the element.
1245  The physical positions of the element's nodes are given by @a elfun.
1246  Note that this function assumes that UpdateTargetSpecification() has
1247  been called with the position vector corresponding to @a elfun. */
1248  virtual void ComputeElementTargets(int e_id, const FiniteElement &fe,
1249  const IntegrationRule &ir,
1250  const Vector &elfun,
1251  DenseTensor &Jtr) const;
1252 
1253  virtual void ComputeAllElementTargets(const FiniteElementSpace &fes,
1254  const IntegrationRule &ir,
1255  const Vector &xe,
1256  DenseTensor &Jtr) const;
1257 
1258  virtual void ComputeElementTargetsGradient(const IntegrationRule &ir,
1259  const Vector &elfun,
1261  DenseTensor &dJtr) const;
1262 
1263  // Generates tspec_vals for target construction using intrule
1264  // Used for the refinement component in hr-adaptivity.
1265  void SetTspecFromIntRule(int e_id, const IntegrationRule &intrule);
1266 
1267  // Targets based on discrete functions can result in invalid (negative)
1268  // size at the quadrature points. This method can be used to set a
1269  // minimum target size.
1270  void SetMinSizeForTargets(double min_size_) { lim_min_size = min_size_; }
1271 
1272  /// Computes target specification data with respect to the coarse FE space.
1274 
1275  // Reset refinement data associated with h-adaptivity component.
1277  {
1278  tspec_refine.Clear();
1279  amr_el = -1;
1280  }
1281 
1282  // Reset derefinement data associated with h-adaptivity component.
1284  {
1286  coarse_tspec_fesv = NULL;
1287  }
1288 
1289  // Used to specify the fine element for determining energy of children of a
1290  // parent element.
1291  void SetRefinementSubElement(int amr_el_) { amr_el = amr_el_; }
1292 };
1293 
1294 class TMOPNewtonSolver;
1295 
1296 /** @brief A TMOP integrator class based on any given TMOP_QualityMetric and
1297  TargetConstructor.
1298 
1299  Represents @f$ \int W(Jpt) dx @f$ over a target zone, where W is the
1300  metric's strain energy density function, and Jpt is the Jacobian of the
1301  target->physical coordinates transformation. The virtual target zone is
1302  defined by the TargetConstructor. */
1304 {
1305 protected:
1306  friend class TMOPNewtonSolver;
1307  friend class TMOPComboIntegrator;
1308 
1311  const TargetConstructor *targetC; // not owned
1312 
1313  // Custom integration rules.
1316 
1317  // Weight Coefficient multiplying the quality metric term.
1318  Coefficient *metric_coeff; // not owned, if NULL -> metric_coeff is 1.
1319  // Normalization factor for the metric term.
1321 
1322  // Nodes and weight Coefficient used for "limiting" the TMOP_Integrator.
1323  // These are both NULL when there is no limiting.
1324  // The class doesn't own lim_nodes0 and lim_coeff.
1327  // Limiting reference distance. Not owned.
1329  // Limiting function. Owned.
1331  // Normalization factor for the limiting term.
1332  double lim_normal;
1333 
1334  // Adaptive limiting.
1335  const GridFunction *adapt_lim_gf0; // Not owned.
1336 #ifdef MFEM_USE_MPI
1338 #endif
1339  GridFunction *adapt_lim_gf; // Owned. Updated by adapt_lim_eval.
1342 
1343  // Surface fitting.
1344  GridFunction *surf_fit_gf; // Owned, Updated by surf_fit_eval.
1345  const Array<bool> *surf_fit_marker; // Not owned.
1346  Coefficient *surf_fit_coeff; // Not owned.
1349 
1351 
1352  // Parameters for FD-based Gradient & Hessian calculation.
1353  bool fdflag;
1354  double dx;
1355  double dxscale;
1356  // Specifies that ComputeElementTargets is being called by a FD function.
1357  // It's used to skip terms that have exact derivative calculations.
1359  // Compute the exact action of the Integrator (includes derivative of the
1360  // target with respect to spatial position)
1362 
1365 
1366  // Jrt: the inverse of the ref->target Jacobian, Jrt = Jtr^{-1}.
1367  // Jpr: the ref->physical transformation Jacobian, Jpr = PMatI^t DS.
1368  // Jpt: the target->physical transformation Jacobian, Jpt = Jpr Jrt.
1369  // P: represents dW_d(Jtp) (dim x dim).
1370  // DSh: gradients of reference shape functions (dof x dim).
1371  // DS: gradients of the shape functions in the target configuration,
1372  // DS = DSh Jrt (dof x dim).
1373  // PMatI: current coordinates of the nodes (dof x dim).
1374  // PMat0: reshaped view into the local element contribution to the operator
1375  // output - the result of AssembleElementVector() (dof x dim).
1377 
1378  // PA extension
1379  // ------------
1380  // E: Q-vector for TMOP-energy
1381  // O: Q-Vector of 1.0, used to compute sums using the dot product kernel.
1382  // X0: E-vector for initial nodal coordinates used for limiting.
1383  // H: Q-Vector for Hessian associated with the metric term.
1384  // C0: Q-Vector for spatial weight used for the limiting term.
1385  // LD: E-Vector constructed using limiting distance grid function (delta).
1386  // H0: Q-Vector for Hessian associated with the limiting term.
1387  //
1388  // maps: Dof2Quad map for fespace associate with nodal coordinates.
1389  // maps_lim: Dof2Quad map for fespace associated with the limiting distance
1390  // grid function.
1391  //
1392  // Jtr_debug_grad
1393  // We keep track if Jtr was set by AssembleGradPA() in Jtr_debug_grad: it
1394  // is set to true by AssembleGradPA(); any other call to
1395  // ComputeAllElementTargets() will set the flag to false. This flag will
1396  // be used to check that Jtr is the one set by AssembleGradPA() when
1397  // performing operations with the gradient like AddMultGradPA() and
1398  // AssembleGradDiagonalPA().
1399  //
1400  // TODO:
1401  // * Merge LD, C0, H0 into one scalar Q-vector
1402  struct
1403  {
1404  bool enabled;
1405  int dim, ne, nq;
1406  mutable DenseTensor Jtr;
1407  mutable bool Jtr_needs_update;
1408  mutable bool Jtr_debug_grad;
1409  mutable Vector E, O, X0, H, C0, LD, H0;
1410  const DofToQuad *maps;
1411  const DofToQuad *maps_lim = nullptr;
1415  } PA;
1416 
1418  double &metric_energy, double &lim_energy,
1419  double &surf_fit_gf_energy);
1420 
1423  const Vector &elfun, Vector &elvect);
1424 
1425  void AssembleElementGradExact(const FiniteElement &el,
1427  const Vector &elfun, DenseMatrix &elmat);
1428 
1429  void AssembleElementVectorFD(const FiniteElement &el,
1431  const Vector &elfun, Vector &elvect);
1432 
1433  // Assumes that AssembleElementVectorFD has been called.
1434  void AssembleElementGradFD(const FiniteElement &el,
1436  const Vector &elfun, DenseMatrix &elmat);
1437 
1438  void AssembleElemVecAdaptLim(const FiniteElement &el,
1440  const IntegrationRule &ir,
1441  const Vector &weights, DenseMatrix &mat);
1442  void AssembleElemGradAdaptLim(const FiniteElement &el,
1444  const IntegrationRule &ir,
1445  const Vector &weights, DenseMatrix &m);
1446 
1447  // First derivative of the surface fitting term.
1448  void AssembleElemVecSurfFit(const FiniteElement &el_x,
1450  DenseMatrix &mat);
1451 
1452  // Second derivative of the surface fitting term.
1453  void AssembleElemGradSurfFit(const FiniteElement &el_x,
1455  DenseMatrix &mat);
1456 
1457  double GetFDDerivative(const FiniteElement &el,
1459  Vector &elfun, const int nodenum,const int idir,
1460  const double baseenergy, bool update_stored);
1461 
1462  /** @brief Determines the perturbation, h, for FD-based approximation. */
1463  void ComputeFDh(const Vector &x, const FiniteElementSpace &fes);
1464 #ifdef MFEM_USE_MPI
1465  void ComputeFDh(const Vector &x, const ParFiniteElementSpace &pfes);
1466 #endif
1467  void ComputeMinJac(const Vector &x, const FiniteElementSpace &fes);
1468 
1469  void UpdateAfterMeshPositionChange(const Vector &new_x);
1470 
1472  {
1473  lim_nodes0 = NULL; lim_coeff = NULL; lim_dist = NULL;
1474  delete lim_func; lim_func = NULL;
1475  }
1476 
1478  {
1479  if (IntegRules)
1480  {
1481  return IntegRules->Get(el.GetGeomType(), integ_order);
1482  }
1483  return (IntRule) ? *IntRule
1484  /* */ : IntRules.Get(el.GetGeomType(), 2*el.GetOrder() + 3);
1485  }
1487  {
1488  // TODO the energy most likely needs less integration points.
1489  return EnergyIntegrationRule(el);
1490  }
1492  {
1493  // TODO the action and energy most likely need less integration points.
1494  return EnergyIntegrationRule(el);
1495  }
1496 
1497  // Auxiliary PA methods
1498  void AssembleGradPA_2D(const Vector&) const;
1499  void AssembleGradPA_3D(const Vector&) const;
1500  void AssembleGradPA_C0_2D(const Vector&) const;
1501  void AssembleGradPA_C0_3D(const Vector&) const;
1502 
1503  double GetLocalStateEnergyPA_2D(const Vector&) const;
1504  double GetLocalStateEnergyPA_C0_2D(const Vector&) const;
1505  double GetLocalStateEnergyPA_3D(const Vector&) const;
1506  double GetLocalStateEnergyPA_C0_3D(const Vector&) const;
1507 
1508  void AddMultPA_2D(const Vector&, Vector&) const;
1509  void AddMultPA_3D(const Vector&, Vector&) const;
1510  void AddMultPA_C0_2D(const Vector&, Vector&) const;
1511  void AddMultPA_C0_3D(const Vector&, Vector&) const;
1512 
1513  void AddMultGradPA_2D(const Vector&, Vector&) const;
1514  void AddMultGradPA_3D(const Vector&, Vector&) const;
1515  void AddMultGradPA_C0_2D(const Vector&, Vector&) const;
1516  void AddMultGradPA_C0_3D(const Vector&, Vector&) const;
1517 
1518  void AssembleDiagonalPA_2D(Vector&) const;
1519  void AssembleDiagonalPA_3D(Vector&) const;
1520  void AssembleDiagonalPA_C0_2D(Vector&) const;
1521  void AssembleDiagonalPA_C0_3D(Vector&) const;
1522 
1523  void AssemblePA_Limiting();
1524  void ComputeAllElementTargets(const Vector &xe = Vector()) const;
1525 
1526 public:
1527  /** @param[in] m TMOP_QualityMetric for r-adaptivity (not owned).
1528  @param[in] tc Target-matrix construction algorithm to use (not owned).
1529  @param[in] hm TMOP_QualityMetric for h-adaptivity (not owned). */
1531  TMOP_QualityMetric *hm)
1532  : h_metric(hm), metric(m), targetC(tc), IntegRules(NULL),
1533  integ_order(-1), metric_coeff(NULL), metric_normal(1.0),
1534  lim_nodes0(NULL), lim_coeff(NULL),
1535  lim_dist(NULL), lim_func(NULL), lim_normal(1.0),
1536  adapt_lim_gf0(NULL), adapt_lim_gf(NULL), adapt_lim_coeff(NULL),
1537  adapt_lim_eval(NULL),
1538  surf_fit_gf(NULL), surf_fit_marker(NULL),
1539  surf_fit_coeff(NULL),
1540  surf_fit_eval(NULL), surf_fit_normal(1.0),
1541  discr_tc(dynamic_cast<DiscreteAdaptTC *>(tc)),
1542  fdflag(false), dxscale(1.0e3), fd_call_flag(false), exact_action(false)
1543  { PA.enabled = false; }
1544 
1546  : TMOP_Integrator(m, tc, m) { }
1547 
1548  ~TMOP_Integrator();
1549 
1550  /// Release the device memory of large PA allocations. This will copy device
1551  /// memory back to the host before releasing.
1552  void ReleasePADeviceMemory(bool copy_to_host = true);
1553 
1554  /// Prescribe a set of integration rules; relevant for mixed meshes.
1555  /** This function has priority over SetIntRule(), if both are called. */
1556  void SetIntegrationRules(IntegrationRules &irules, int order)
1557  {
1558  IntegRules = &irules;
1559  integ_order = order;
1560  }
1561 
1562  /// Sets a scaling Coefficient for the quality metric term of the integrator.
1563  /** With this addition, the integrator becomes
1564  @f$ \int w1 W(Jpt) dx @f$.
1565 
1566  Note that the Coefficient is evaluated in the physical configuration and
1567  not in the target configuration which may be undefined. */
1569 
1570  /** @brief Limiting of the mesh displacements (general version).
1571 
1572  Adds the term @f$ \int w_0 f(x, x_0, d) dx @f$, where f is a measure of
1573  the displacement between x and x_0, given the max allowed displacement d.
1574 
1575  @param[in] n0 Original mesh node coordinates (x0 above).
1576  @param[in] dist Allowed displacement in physical space (d above).
1577  @param[in] w0 Coefficient scaling the limiting integral.
1578  @param[in] lfunc TMOP_LimiterFunction defining the function f. If
1579  NULL, a TMOP_QuadraticLimiter will be used. The
1580  TMOP_Integrator assumes ownership of this pointer. */
1581  void EnableLimiting(const GridFunction &n0, const GridFunction &dist,
1582  Coefficient &w0, TMOP_LimiterFunction *lfunc = NULL);
1583 
1584  /** @brief Adds a limiting term to the integrator with limiting distance
1585  function (@a dist in the general version of the method) equal to 1. */
1586  void EnableLimiting(const GridFunction &n0, Coefficient &w0,
1587  TMOP_LimiterFunction *lfunc = NULL);
1588 
1589  /** @brief Restriction of the node positions to certain regions.
1590 
1591  Adds the term @f$ \int c (z(x) - z_0(x_0))^2 @f$, where z0(x0) is a given
1592  function on the starting mesh, and z(x) is its image on the new mesh.
1593  Minimizing this term means that a node at x0 is allowed to move to a
1594  position x(x0) only if z(x) ~ z0(x0).
1595  Such term can be used for tangential mesh relaxation.
1596 
1597  @param[in] z0 Function z0 that controls the adaptive limiting.
1598  @param[in] coeff Coefficient c for the above integral.
1599  @param[in] ae AdaptivityEvaluator to compute z(x) from z0(x0). */
1600  void EnableAdaptiveLimiting(const GridFunction &z0, Coefficient &coeff,
1601  AdaptivityEvaluator &ae);
1602 #ifdef MFEM_USE_MPI
1603  /// Parallel support for adaptive limiting.
1604  void EnableAdaptiveLimiting(const ParGridFunction &z0, Coefficient &coeff,
1605  AdaptivityEvaluator &ae);
1606 #endif
1607 
1608  /** @brief Fitting of certain DOFs to the zero level set of a function.
1609 
1610  Having a level set function s0(x0) on the starting mesh, and a set of
1611  marked nodes (or DOFs), we move these nodes to the zero level set of s0.
1612  If s(x) is the image of s0(x0) on the current mesh, this function adds to
1613  the TMOP functional the term @f$ \int c \bar{s}(x))^2 @f$, where
1614  @f$\bar{s}(x)@f$ is the restriction of s(x) on the aligned DOFs.
1615  Minimizing this term means that a marked node at x0 is allowed to move to
1616  a position x(x0) only if s(x) ~ 0.
1617  Such term can be used for surface fitting and tangential relaxation.
1618 
1619  @param[in] s0 The level set function on the initial mesh.
1620  @param[in] smarker Indicates which DOFs will be aligned.
1621  @param[in] coeff Coefficient c for the above integral.
1622  @param[in] ae AdaptivityEvaluator to compute s(x) from s0(x0). */
1623  void EnableSurfaceFitting(const GridFunction &s0,
1624  const Array<bool> &smarker, Coefficient &coeff,
1625  AdaptivityEvaluator &ae);
1626 #ifdef MFEM_USE_MPI
1627  /// Parallel support for surface fitting.
1628  void EnableSurfaceFitting(const ParGridFunction &s0,
1629  const Array<bool> &smarker, Coefficient &coeff,
1630  AdaptivityEvaluator &ae);
1631 #endif
1632  void GetSurfaceFittingErrors(double &err_avg, double &err_max);
1633  bool IsSurfaceFittingEnabled() { return (surf_fit_gf != NULL); }
1634 
1635  /// Update the original/reference nodes used for limiting.
1636  void SetLimitingNodes(const GridFunction &n0) { lim_nodes0 = &n0; }
1637 
1638  /** @brief Computes the integral of W(Jacobian(Trt)) over a target zone.
1639  @param[in] el Type of FiniteElement.
1640  @param[in] T Mesh element transformation.
1641  @param[in] elfun Physical coordinates of the zone. */
1642  virtual double GetElementEnergy(const FiniteElement &el,
1644  const Vector &elfun);
1645 
1646  /** @brief Computes the mean of the energies of the given element's children.
1647 
1648  In addition to the inputs for GetElementEnergy, this function requires an
1649  IntegrationRule to be specified that will give the decomposition of the
1650  given element based on the refinement type being considered. */
1651  virtual double GetRefinementElementEnergy(const FiniteElement &el,
1653  const Vector &elfun,
1654  const IntegrationRule &irule);
1655 
1656  /// This function is similar to GetElementEnergy, but ignores components
1657  /// such as limiting etc. to compute the element energy.
1658  virtual double GetDerefinementElementEnergy(const FiniteElement &el,
1660  const Vector &elfun);
1661 
1662  virtual void AssembleElementVector(const FiniteElement &el,
1664  const Vector &elfun, Vector &elvect);
1665 
1666  virtual void AssembleElementGrad(const FiniteElement &el,
1668  const Vector &elfun, DenseMatrix &elmat);
1669 
1671 
1673 #ifdef MFEM_USE_MPI
1675 #endif
1676 
1677  // PA extension
1679  virtual void AssemblePA(const FiniteElementSpace&);
1680 
1681  virtual void AssembleGradPA(const Vector&, const FiniteElementSpace&);
1682 
1683  virtual double GetLocalStateEnergyPA(const Vector&) const;
1684 
1685  virtual void AddMultPA(const Vector&, Vector&) const;
1686 
1687  virtual void AddMultGradPA(const Vector&, Vector&) const;
1688 
1689  virtual void AssembleGradDiagonalPA(Vector&) const;
1690 
1692 
1693  /** @brief Computes the normalization factors of the metric and limiting
1694  integrals using the mesh position given by @a x. */
1695  void EnableNormalization(const GridFunction &x);
1696 #ifdef MFEM_USE_MPI
1697  void ParEnableNormalization(const ParGridFunction &x);
1698 #endif
1699 
1700  /** @brief Enables FD-based approximation and computes dx. */
1701  void EnableFiniteDifferences(const GridFunction &x);
1702 #ifdef MFEM_USE_MPI
1704 #endif
1705 
1706  void SetFDhScale(double dxscale_) { dxscale = dxscale_; }
1707  bool GetFDFlag() const { return fdflag; }
1708  double GetFDh() const { return dx; }
1709 
1710  /** @brief Flag to control if exact action of Integration is effected. */
1711  void SetExactActionFlag(bool flag_) { exact_action = flag_; }
1712 
1713  /// Update the surface fitting weight as surf_fit_coeff *= factor;
1714  void UpdateSurfaceFittingWeight(double factor);
1715 
1716  /// Get the surface fitting weight.
1717  double GetSurfaceFittingWeight();
1718 };
1719 
1721 {
1722 protected:
1723  // Integrators in the combination. Owned.
1725 
1726 public:
1728 
1730  {
1731  for (int i = 0; i < tmopi.Size(); i++) { delete tmopi[i]; }
1732  }
1733 
1734  /// Adds a new TMOP_Integrator to the combination.
1735  void AddTMOPIntegrator(TMOP_Integrator *ti) { tmopi.Append(ti); }
1736 
1738 
1739  /// Adds the limiting term to the first integrator. Disables it for the rest.
1740  void EnableLimiting(const GridFunction &n0, const GridFunction &dist,
1741  Coefficient &w0, TMOP_LimiterFunction *lfunc = NULL);
1742 
1743  /** @brief Adds the limiting term to the first integrator. Disables it for
1744  the rest (@a dist in the general version of the method) equal to 1. */
1745  void EnableLimiting(const GridFunction &n0, Coefficient &w0,
1746  TMOP_LimiterFunction *lfunc = NULL);
1747 
1748  /// Update the original/reference nodes used for limiting.
1749  void SetLimitingNodes(const GridFunction &n0);
1750 
1751  virtual double GetElementEnergy(const FiniteElement &el,
1753  const Vector &elfun);
1754  virtual void AssembleElementVector(const FiniteElement &el,
1756  const Vector &elfun, Vector &elvect);
1757  virtual void AssembleElementGrad(const FiniteElement &el,
1759  const Vector &elfun, DenseMatrix &elmat);
1760 
1761  virtual double GetRefinementElementEnergy(const FiniteElement &el,
1763  const Vector &elfun,
1764  const IntegrationRule &irule);
1765 
1766  virtual double GetDerefinementElementEnergy(const FiniteElement &el,
1768  const Vector &elfun);
1769 
1770  /// Normalization factor that considers all integrators in the combination.
1771  void EnableNormalization(const GridFunction &x);
1772 #ifdef MFEM_USE_MPI
1773  void ParEnableNormalization(const ParGridFunction &x);
1774 #endif
1775 
1776  // PA extension
1778  virtual void AssemblePA(const FiniteElementSpace&);
1779  virtual void AssembleGradPA(const Vector&, const FiniteElementSpace&);
1780  virtual double GetLocalStateEnergyPA(const Vector&) const;
1781  virtual void AddMultPA(const Vector&, Vector&) const;
1782  virtual void AddMultGradPA(const Vector&, Vector&) const;
1783  virtual void AssembleGradDiagonalPA(Vector&) const;
1784 };
1785 
1786 /// Interpolates the @a metric's values at the nodes of @a metric_gf.
1787 /** Assumes that @a metric_gf's FiniteElementSpace is initialized. */
1788 void InterpolateTMOP_QualityMetric(TMOP_QualityMetric &metric,
1789  const TargetConstructor &tc,
1790  const Mesh &mesh, GridFunction &metric_gf);
1791 }
1792 
1793 #endif
Abstract class for all finite elements.
Definition: fe_base.hpp:235
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:180
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:343
virtual double GetDerefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
Definition: tmop.cpp:2714
void AssembleGradPA_2D(const Vector &) 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...
Definition: tmop.cpp:344
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:74
virtual void SetParDiscreteTargetSize(const ParGridFunction &tspec_)
Definition: tmop.cpp:1375
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:73
Shifted barrier form of metric 56 (area, ideal barrier metric), 2D.
Definition: tmop.hpp:430
Shifted barrier form of 3D metric 16 (volume, ideal barrier metric), 3D.
Definition: tmop.hpp:682
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:961
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:314
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:506
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)=0
void AddMultGradPA_C0_3D(const Vector &, Vector &) const
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 AddMultPA_C0_2D(const Vector &, Vector &) const
void AssembleElemVecSurfFit(const FiniteElement &el_x, IsoparametricTransformation &Tpr, DenseMatrix &mat)
Definition: tmop.cpp:3142
void SetSerialMetaInfo(const Mesh &m, const FiniteElementCollection &fec, int num_comp)
Definition: tmop.cpp:2269
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:150
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
VectorCoefficient * vector_tspec
Definition: tmop.hpp:1053
virtual void AddMultGradPA(const Vector &, Vector &) const
Method for partially assembled gradient action.
void EnableFiniteDifferences(const GridFunction &x)
Enables FD-based approximation and computes dx.
Definition: tmop.cpp:3642
const TargetType target_type
Definition: tmop.hpp:935
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:768
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:765
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
DenseMatrix PMatO
Definition: tmop.hpp:1376
virtual void AddMultPA(const Vector &, Vector &) const
Method for partially assembled action.
Definition: tmop.cpp:3826
DenseTensor Jtr
Definition: tmop.hpp:1406
double & min_detT
Definition: tmop.hpp:256
3D barrier Shape+Size (VS) metric.
Definition: tmop.hpp:578
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:525
const DenseMatrix * Jtr
Definition: tmop.hpp:26
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:581
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:664
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:674
double GetSurfaceFittingWeight()
Get the surface fitting weight.
Definition: tmop.cpp:3469
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:402
3D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:661
Base class for vector Coefficients that optionally depend on time and space.
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).
Definition: tmop.hpp:177
virtual ~TMOP_Metric_334()
Definition: tmop.hpp:678
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:732
DenseMatrix Jrt
Definition: tmop.hpp:1376
const GridFunction * lim_dist
Definition: tmop.hpp:1328
virtual void SetSerialDiscreteTargetOrientation(const GridFunction &tspec_)
Definition: tmop.cpp:1474
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
Compute the local energy.
Definition: tmop.cpp:3699
double epsilon
Definition: ex25.cpp:140
virtual ~TMOP_Metric_333()
Definition: tmop.hpp:657
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:838
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.
Definition: tmop.hpp:354
double GetGamma() const
Definition: tmop.hpp:634
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:218
TMOP_AMetric_126(double gamma_)
Definition: tmop.hpp:782
virtual void AssemblePA(const FiniteElementSpace &)
Method defining partial assembly.
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:291
const double eps
Definition: tmop.hpp:505
void EnableNormalization(const GridFunction &x)
Computes the normalization factors of the metric and limiting integrals using the mesh position given...
Definition: tmop.cpp:3480
TMOP_QualityMetric & GetAMRQualityMetric()
Definition: tmop.hpp:1670
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
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:3852
3D barrier Shape (S) metric.
Definition: tmop.hpp:484
TMOP_QualityMetric * metric
Definition: tmop.hpp:1310
void SetAdaptivityEvaluator(AdaptivityEvaluator *ae)
Definition: tmop.hpp:1233
virtual ~TMOP_Metric_332()
Definition: tmop.hpp:636
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:308
virtual ~TMOP_LimiterFunction()
Virtual destructor.
Definition: tmop.hpp:815
Base class for limiting functions to be used in class TMOP_Integrator.
Definition: tmop.hpp:798
virtual double EvalW(const DenseMatrix &Jpt) const =0
Evaluate the strain energy density function, W = W(Jpt).
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:935
const GridFunction * GetNodes() const
Get the nodes to be used in the target-matrix construction.
Definition: tmop.hpp:991
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 void ComputeElementTargetsGradient(const IntegrationRule &ir, const Vector &elfun, IsoparametricTransformation &Tpr, DenseTensor &dJtr) const
Definition: tmop.cpp:1832
3D non-barrier Skew metric.
Definition: tmop.hpp:141
2D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:358
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:866
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:793
void ParUpdateAfterMeshTopologyChange()
Definition: tmop.cpp:1343
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:202
DenseMatrix Jpr
Definition: tmop.hpp:1376
IntegrationRules * IntegRules
Definition: tmop.hpp:1314
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...
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:363
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:779
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:584
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:993
void SetFDhScale(double dxscale_)
Definition: tmop.hpp:1706
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:111
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe_base.hpp:327
void SetMinSizeForTargets(double min_size_)
Definition: tmop.hpp:1270
virtual void ComputeElementTargetsGradient(const IntegrationRule &ir, const Vector &elfun, IsoparametricTransformation &Tpr, DenseTensor &dJtr) const
Definition: tmop.cpp:1270
Container class for integration rules.
Definition: intrules.hpp:311
ParFiniteElementSpace * pfes
Definition: tmop.hpp:863
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:664
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
3D barrier Shape (S) metric.
Definition: tmop.hpp:466
void GetSurfaceFittingErrors(double &err_avg, double &err_max)
Definition: tmop.cpp:2436
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
void SetVolumeScale(double vol_scale)
Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
Definition: tmop.hpp:994
2D barrier shape (S) metric (polyconvex).
Definition: tmop.hpp:186
void EnableAdaptiveLimiting(const GridFunction &z0, Coefficient &coeff, AdaptivityEvaluator &ae)
Restriction of the node positions to certain regions.
Definition: tmop.cpp:2365
virtual void AddMultPA(const Vector &, Vector &) const
Method for partially assembled action.
3D non-barrier Size (V) metric.
Definition: tmop.hpp:542
2D barrier size (V) metric (polyconvex).
Definition: tmop.hpp:340
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:737
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition: tmop.hpp:1556
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:759
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:122
void EnableNormalization(const GridFunction &x)
Normalization factor that considers all integrators in the combination.
Definition: tmop.cpp:3769
2D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:774
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:596
bool IsSurfaceFittingEnabled()
Definition: tmop.hpp:1633
double & min_detT
Definition: tmop.hpp:524
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:162
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:520
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:147
virtual void AssemblePA(const FiniteElementSpace &)
Method defining partial assembly.
Definition: tmop.cpp:3801
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:200
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:225
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:244
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
Definition: mesh.hpp:1814
const IntegrationRule * ir
Definition: tmop.hpp:1414
void SetExactActionFlag(bool flag_)
Flag to control if exact action of Integration is effected.
Definition: tmop.hpp:1711
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:326
void AssembleElementVectorFD(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
Definition: tmop.cpp:3294
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:480
const GridFunction * nodes
Definition: tmop.hpp:932
virtual void SetSerialDiscreteTargetSkew(const GridFunction &tspec_)
Definition: tmop.cpp:1458
void UpdateSurfaceFittingWeight(double factor)
Update the surface fitting weight as surf_fit_coeff *= factor;.
Definition: tmop.cpp:3457
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:556
struct mfem::TMOP_Integrator::@13 PA
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:390
void UpdateGradientTargetSpecification(const Vector &x, double dx, bool use_flag=false)
Definition: tmop.cpp:2179
void SetRefinementSubElement(int amr_el_)
Definition: tmop.hpp:1291
2D non-barrier Aspect ratio metric.
Definition: tmop.hpp:156
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:643
TMOP_Metric_080(double gamma_)
Definition: tmop.hpp:366
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:855
2D barrier Shape+Size (VS) metric (not polyconvex).
Definition: tmop.hpp:204
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:645
Default limiter function in TMOP_Integrator.
Definition: tmop.hpp:819
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:1617
void AssembleElementGradExact(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
Definition: tmop.cpp:2937
Array< Vector * > ElemDer
Definition: tmop.hpp:1363
virtual ~TMOP_AMetric_126()
Definition: tmop.hpp:791
const Vector & GetTspecPertMixH()
Definition: tmop.hpp:1241
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.
Definition: tmop.cpp:3670
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:599
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:432
double GetFDh() const
Definition: tmop.hpp:1708
2D Shifted barrier form of shape metric (mu_2).
Definition: tmop.hpp:253
TMOP_Metric_313(double &mindet)
Definition: tmop.hpp:528
TMOP_Metric_352(double &t0)
Definition: tmop.hpp:689
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:1236
Coefficient * scalar_tspec
Definition: tmop.hpp:1052
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:620
virtual void SetSerialDiscreteTargetAspectRatio(const GridFunction &tspec_)
Definition: tmop.cpp:1466
void UpdateTargetSpecification(const Vector &new_x, bool use_flag=false)
Definition: tmop.cpp:1534
TMOPMatrixCoefficient(int dim)
Definition: tmop.hpp:1037
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:720
2D barrier shape (S) metric (not polyconvex).
Definition: tmop.hpp:323
const DofToQuad * maps
Definition: tmop.hpp:1410
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:107
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:710
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:448
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:33
const GeometricFactors * geom
Definition: tmop.hpp:1412
virtual void AddMultGradPA(const Vector &, Vector &) const
Method for partially assembled gradient action.
Definition: tmop.cpp:3834
virtual void ComputeElementTargetsGradient(const IntegrationRule &ir, const Vector &elfun, IsoparametricTransformation &Tpr, DenseTensor &dJtr) const
Definition: tmop.cpp:1213
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
2D barrier Size (V) metric (polyconvex).
Definition: tmop.hpp:720
Abstract class for local mesh quality metrics in the target-matrix optimization paradigm (TMOP) by P...
Definition: tmop.hpp:23
void FinalizeSerialDiscreteTargetSpec(const GridFunction &tspec_)
Definition: tmop.cpp:1482
Coefficient * adapt_lim_coeff
Definition: tmop.hpp:1340
virtual ~TMOPMatrixCoefficient()
Definition: tmop.hpp:1045
virtual ~AdaptivityEvaluator()
Definition: tmop.cpp:2304
void UpdateAfterMeshPositionChange(const Vector &new_x)
Definition: tmop.cpp:3608
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:545
virtual void AddQualityMetric(TMOP_QualityMetric *tq, double wt=1.0)
Definition: tmop.hpp:85
void SetLimitingNodes(const GridFunction &n0)
Update the original/reference nodes used for limiting.
Definition: tmop.cpp:3691
virtual ~TMOP_Metric_080()
Definition: tmop.hpp:377
void SetDiscreteTargetBase(const GridFunction &tspec_)
Definition: tmop.cpp:1413
AdaptivityEvaluator * adapt_lim_eval
Definition: tmop.hpp:1341
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:822
DiscreteAdaptTC * discr_tc
Definition: tmop.hpp:1350
void EnableLimiting(const GridFunction &n0, const GridFunction &dist, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
Limiting of the mesh displacements (general version).
Definition: tmop.cpp:2340
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:554
void AddMultGradPA_C0_2D(const Vector &, Vector &) const
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:705
void SetCoefficient(Coefficient &w1)
Sets a scaling Coefficient for the quality metric term of the integrator.
Definition: tmop.hpp:1568
2D non-barrier size (V) metric (not polyconvex).
Definition: tmop.hpp:288
void UpdateTargetSpecificationAtNode(const FiniteElement &el, ElementTransformation &T, int nodenum, int idir, const Vector &IntData)
Definition: tmop.cpp:1552
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:666
void AssembleDiagonalPA_C0_2D(Vector &) const
void SetTspecFromIntRule(int e_id, const IntegrationRule &intrule)
Definition: tmop.cpp:1583
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:390
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:751
TargetConstructor(TargetType ttype)
Constructor for use in serial.
Definition: tmop.hpp:960
TargetType GetTargetType() const
Definition: tmop.hpp:996
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:46
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:419
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:264
3D barrier Size (V) metric.
Definition: tmop.hpp:560
DiscreteAdaptTC(TargetType ttype)
Definition: tmop.hpp:1147
void ResetRefinementTspecData()
Definition: tmop.hpp:1276
void AssembleDiagonalPA_3D(Vector &) const
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:415
void SetNodes(const GridFunction &n)
Set the nodes to be used in the target-matrix construction.
Definition: tmop.hpp:988
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:698
void AssembleElementVectorExact(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
Definition: tmop.cpp:2797
void AssembleDiagonalPA_2D(Vector &) const
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:601
virtual void AssembleGradDiagonalPA(Vector &) const
Method for computing the diagonal of the gradient with partial assembly.
Definition: tmop.cpp:3818
const TargetConstructor * targetC
Definition: tmop.hpp:1311
Array< Vector * > ElemPertEnergy
Definition: tmop.hpp:1364
const GridFunction * lim_nodes0
Definition: tmop.hpp:1325
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:563
const ParGridFunction * adapt_lim_pgf0
Definition: tmop.hpp:1337
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:645
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:621
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
Definition: nonlininteg.cpp:25
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:711
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:478
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:633
Abstract class used to define combination of metrics with constant coefficients.
Definition: tmop.hpp:78
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:257
virtual void AssembleGradPA(const Vector &, const FiniteElementSpace &)
Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at th...
Definition: tmop.cpp:3809
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:147
virtual void SetParDiscreteTargetAspectRatio(const ParGridFunction &tspec_)
Definition: tmop.cpp:1391
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 ~TargetConstructor()
Definition: tmop.hpp:974
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:491
void ParEnableNormalization(const ParGridFunction &x)
Definition: tmop.cpp:3785
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:434
double GetLocalStateEnergyPA_2D(const Vector &) const
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
Definition: tmop.cpp:2768
virtual double GetLocalStateEnergyPA(const Vector &) const
Compute the local (to the MPI rank) energy with partial assembly.
Definition: tmop.cpp:3842
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:396
GridFunction * GetTSpecData()
Get the entire tspec.
Definition: tmop.hpp:1197
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:777
FiniteElementSpace * tspec_fesv
Definition: tmop.hpp:1117
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:239
Coefficient * metric_coeff
Definition: tmop.hpp:1318
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:67
void AssembleElemGradSurfFit(const FiniteElement &el_x, IsoparametricTransformation &Tpr, DenseMatrix &mat)
Definition: tmop.cpp:3196
virtual bool ContainsVolumeInfo() const
Checks if the target matrices contain non-trivial size specification.
Definition: tmop.cpp:1128
TMOP_QualityMetric * h_metric
Definition: tmop.hpp:1309
Array< TMOP_QualityMetric * > tmop_q_arr
Definition: tmop.hpp:81
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:729
void RestoreTargetSpecificationAtNode(ElementTransformation &T, int nodenum)
Definition: tmop.cpp:1569
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:231
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:363
2D barrier Shape+Orientation (OS) metric (polyconvex).
Definition: tmop.hpp:756
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:387
void SetTspecAtIndex(int idx, const GridFunction &tspec_)
Definition: tmop.cpp:1437
virtual double GetRefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, const IntegrationRule &irule)
Computes the mean of the energies of the given element&#39;s children.
Definition: tmop.cpp:2628
void ComputeAvgVolume() const
Definition: tmop.cpp:1014
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:844
Coefficient * lim_coeff
Definition: tmop.hpp:1326
3D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:596
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:23
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:361
TMOP_Metric_328(double gamma_)
Definition: tmop.hpp:604
virtual double Eval(const Vector &x, const Vector &x0, double dist) const
Returns the limiting function, f(x, x0, d).
Definition: tmop.hpp:822
AdaptivityEvaluator * surf_fit_eval
Definition: tmop.hpp:1347
virtual ~TMOP_QuadraticLimiter()
Definition: tmop.hpp:846
FiniteElementSpace * GetTSpecFESpace()
Get the FESpace associated with tspec.
Definition: tmop.hpp:1195
2D barrier (not a shape) metric (polyconvex).
Definition: tmop.hpp:272
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:779
virtual void SetTargetJacobian(const DenseMatrix &Jtr_)
Specify the reference-element -&gt; target-element Jacobian matrix for the point of interest.
Definition: tmop.hpp:43
TMOP_Metric_211(double epsilon=1e-4)
Definition: tmop.hpp:418
void UpdateHessianTargetSpecification(const Vector &x, double dx, bool use_flag=false)
Definition: tmop.cpp:2205
virtual void AssembleGradDiagonalPA(Vector &) const
Method for computing the diagonal of the gradient with partial assembly.
void AssembleElemGradAdaptLim(const FiniteElement &el, IsoparametricTransformation &Tpr, const IntegrationRule &ir, const Vector &weights, DenseMatrix &m)
Definition: tmop.cpp:3081
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
Definition: tmop.cpp:2782
const DofToQuad * maps_lim
Definition: tmop.hpp:1411
2D untangling metric.
Definition: tmop.hpp:411
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:132
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
3D Size (V) untangling metric.
Definition: tmop.hpp:502
const double eps
Definition: tmop.hpp:414
virtual void SetTargetJacobian(const DenseMatrix &Jtr_)
Specify the reference-element -&gt; target-element Jacobian matrix for the point of interest.
Definition: tmop.hpp:91
2D barrier size (V) metric (polyconvex).
Definition: tmop.hpp:305
void SetLimitingNodes(const GridFunction &n0)
Update the original/reference nodes used for limiting.
Definition: tmop.hpp:1636
DenseMatrix tspec_refine
Definition: tmop.hpp:1106
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:498
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:686
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition: vector.cpp:447
TMOPMatrixCoefficient * matrix_tspec
Definition: tmop.hpp:1054
virtual double GetLocalStateEnergyPA(const Vector &) const
Compute the local (to the MPI rank) energy with partial assembly.
TMOP_Metric_333(double gamma_)
Definition: tmop.hpp:648
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:1142
const IntegrationRule & ActionIntegrationRule(const FiniteElement &el) const
Definition: tmop.hpp:1486
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...
Base class for Matrix Coefficients that optionally depend on time and space.
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:275
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:83
void ReleasePADeviceMemory(bool copy_to_host=true)
Definition: tmop.cpp:2314
bool UsesPhysicalCoordinates() const
Return true if the methods ComputeElementTargets(), ComputeAllElementTargets(), and ComputeElementTar...
Definition: tmop.hpp:1001
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:214
void Clear()
Delete the matrix data array (if owned) and reset the matrix state.
Definition: densemat.hpp:98
void DisableLimiting()
Definition: tmop.hpp:1471
double GetLocalStateEnergyPA_C0_2D(const Vector &) const
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
Definition: tmop.cpp:3727
double DistanceSquaredTo(const double *p) const
Compute the square of the Euclidean distance to another vector.
Definition: vector.hpp:657
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:371
DiscreteAdaptTC * GetDiscreteAdaptTC() const
Definition: tmop.hpp:1691
void AddMultGradPA_2D(const Vector &, Vector &) const
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
Definition: fe_base.hpp:136
TMOP_Metric_252(double &t0)
Note that t0 is stored by reference.
Definition: tmop.hpp:438
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:403
DenseTensor Jtrcomp
Definition: tmop.hpp:1113
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:576
const IntegrationRule * IntRule
Definition: nonlininteg.hpp:30
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:374
void ResetUpdateFlags()
Used in combination with the Update methods to avoid extra computations.
Definition: tmop.hpp:1189
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:592
AnalyticAdaptTC(TargetType ttype)
Definition: tmop.hpp:1057
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:741
void AddTMOPIntegrator(TMOP_Integrator *ti)
Adds a new TMOP_Integrator to the combination.
Definition: tmop.hpp:1735
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:462
GridFunction * tspec_gf
Definition: tmop.hpp:1119
2D non-barrier Skew metric.
Definition: tmop.hpp:126
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:752
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:298
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:224
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:281
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:801
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:470
void UpdateAfterMeshTopologyChange()
Definition: tmop.cpp:2465
void AddMultPA_3D(const Vector &, Vector &) const
DenseMatrix PMatI
Definition: tmop.hpp:1376
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:563
2D barrier Shape+Size (VS) metric (not polyconvex).
Definition: tmop.hpp:222
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:272
void AddMultPA_2D(const Vector &, Vector &) const
void SetTransformation(ElementTransformation &)
The method SetTransformation() is hidden for TMOP_QualityMetrics, because it is not used...
Definition: tmop.hpp:31
void ComputeFDh(const Vector &x, const FiniteElementSpace &fes)
Determines the perturbation, h, for FD-based approximation.
Definition: tmop.cpp:3624
void AssembleGradPA_C0_2D(const Vector &) const
void ComputeNormalizationEnergies(const GridFunction &x, double &metric_energy, double &lim_energy, double &surf_fit_gf_energy)
Definition: tmop.cpp:3503
TMOP_Metric_334(double gamma_)
Definition: tmop.hpp:669
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
Definition: tmop.hpp:1346
Class for integration point with weight.
Definition: intrules.hpp:25
AdaptivityEvaluator * adapt_eval
Definition: tmop.hpp:1136
ParFiniteElementSpace * ptspec_fesv
Definition: tmop.hpp:1122
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:453
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:405
void AssembleDiagonalPA_C0_3D(Vector &) const
const FiniteElementSpace * fes
Definition: tmop.hpp:1413
Array< TMOP_Integrator * > tmopi
Definition: tmop.hpp:1724
A standard isoparametric element transformation.
Definition: eltrans.hpp:361
void ComputeMinJac(const Vector &x, const FiniteElementSpace &fes)
Definition: tmop.cpp:3574
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:682
TMOP_Integrator(TMOP_QualityMetric *m, TargetConstructor *tc, TMOP_QualityMetric *hm)
Definition: tmop.hpp:1530
void Diag(double c, int n)
Creates n x n diagonal matrix with diagonal elements c.
Definition: densemat.cpp:1329
This class is used to express the local action of a general nonlinear finite element operator...
Definition: nonlininteg.hpp:27
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:750
const IntegrationRule & EnergyIntegrationRule(const FiniteElement &el) const
Definition: tmop.hpp:1477
const IntegrationRule & GradientIntegrationRule(const FiniteElement &el) const
Definition: tmop.hpp:1491
double GetFDDerivative(const FiniteElement &el, ElementTransformation &T, Vector &elfun, const int nodenum, const int idir, const double baseenergy, bool update_stored)
Definition: tmop.cpp:3272
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:747
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:738
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:714
void SetParMetaInfo(const ParMesh &m, const FiniteElementCollection &fec, int num_comp)
Parallel version of SetSerialMetaInfo.
Definition: tmop.cpp:2282
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:440
double surf_fit_normal
Definition: tmop.hpp:1348
3D Shape (S) metric, untangling version of 303.
Definition: tmop.hpp:521
int dim
Definition: ex24.cpp:53
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:2491
void AssembleElemVecAdaptLim(const FiniteElement &el, IsoparametricTransformation &Tpr, const IntegrationRule &ir, const Vector &weights, DenseMatrix &mat)
Definition: tmop.cpp:3045
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:61
void ParUpdateAfterMeshTopologyChange()
Definition: tmop.cpp:2478
virtual double GetDerefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
Definition: tmop.cpp:3756
void AssembleGradPA_3D(const Vector &) const
void Destroy()
Destroy a vector.
Definition: vector.hpp:598
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.
void AddMultGradPA_3D(const Vector &, Vector &) const
void ComputeAllElementTargets_Fallback(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
Definition: tmop.cpp:1056
3D non-barrier Aspect ratio metric.
Definition: tmop.hpp:171
const Vector & GetTspecPert2H()
Definition: tmop.hpp:1240
const Vector & GetTspecPert1H()
Definition: tmop.hpp:1239
void AssembleElementGradFD(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
Definition: tmop.cpp:3360
ParFiniteElementSpace * GetTSpecParFESpace()
Definition: tmop.hpp:1202
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:666
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:814
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:189
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:629
virtual void SetParDiscreteTargetSpec(const ParGridFunction &tspec_)
Definition: tmop.cpp:1407
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:207
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:601
Abstract class for hyperelastic models.
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:167
bool ComputeAllElementTargets(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:334
GridFunction * adapt_lim_gf
Definition: tmop.hpp:1339
TargetConstructor(TargetType ttype, MPI_Comm mpicomm)
Constructor for use in parallel.
Definition: tmop.hpp:970
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:830
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:247
GridFunction * surf_fit_gf
Definition: tmop.hpp:1344
void ComputeAllElementTargets(const Vector &xe=Vector()) const
TMOP_Integrator(TMOP_QualityMetric *m, TargetConstructor *tc)
Definition: tmop.hpp:1545
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:614
void FinalizeParDiscreteTargetSpec(const ParGridFunction &tspec_)
Definition: tmop.cpp:1317
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:538
ParGridFunction * tspec_pgf
Definition: tmop.hpp:1124
bool GetFDFlag() const
Definition: tmop.hpp:1707
Vector data type.
Definition: vector.hpp:60
void GetDiscreteTargetSpec(GridFunction &tspec_, int idx)
Get one of the discrete fields from tspec.
Definition: tmop.cpp:1502
DenseMatrix DSh
Definition: tmop.hpp:1376
const Array< bool > * surf_fit_marker
Definition: tmop.hpp:1345
DenseMatrix DS
Definition: tmop.hpp:1376
TMOP_Metric_022(double &t0)
Definition: tmop.hpp:260
3D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:640
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:487
void SetTspecDataForDerefinement(FiniteElementSpace *fes)
Computes target specification data with respect to the coarse FE space.
Definition: tmop.cpp:1609
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:538
TargetType
Target-matrix construction algorithms supported by this class.
Definition: tmop.hpp:912
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:411
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.
Definition: tmop.cpp:2400
double GetLocalStateEnergyPA_C0_3D(const Vector &) 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...
Definition: tmop.cpp:568
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:729
TMOP_Metric_311(double epsilon=1e-4)
Definition: tmop.hpp:509
virtual void SetSerialDiscreteTargetSpec(const GridFunction &tspec_)
Definition: tmop.cpp:1528
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:723
virtual void SetParDiscreteTargetSkew(const ParGridFunction &tspec_)
Definition: tmop.cpp:1383
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:382
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:829
virtual double GetRefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, const IntegrationRule &irule)
Definition: tmop.cpp:3743
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:165
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:777
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:498
virtual void SetParDiscreteTargetOrientation(const ParGridFunction &tspec_)
Definition: tmop.cpp:1399
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:908
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:977
void AddMultPA_C0_3D(const Vector &, Vector &) const
virtual ~DiscreteAdaptTC()
Definition: tmop.cpp:2259
TMOP_Metric_332(double gamma_)
Definition: tmop.hpp:624
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:208
virtual void SetAnalyticTargetSpec(Coefficient *sspec, VectorCoefficient *vspec, TMOPMatrixCoefficient *mspec)
Definition: tmop.cpp:1227
3D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:617
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:786
Class for parallel meshes.
Definition: pmesh.hpp:32
virtual void SetSerialDiscreteTargetSize(const GridFunction &tspec_)
Definition: tmop.cpp:1450
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:621
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:379
TMOP_LimiterFunction * lim_func
Definition: tmop.hpp:1330
const Array< TMOP_Integrator * > & GetTMOPIntegrators() const
Definition: tmop.hpp:1737
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:785
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:704
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:469
const GridFunction * adapt_lim_gf0
Definition: tmop.hpp:1335
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
Definition: tmop.cpp:3711
DenseMatrix Jpt
Definition: tmop.hpp:1376
double GetGamma() const
Definition: tmop.hpp:375
2D barrier Shape+Orientation (OS) metric (polyconvex).
Definition: tmop.hpp:381
Array< double > wt_arr
Definition: tmop.hpp:82
void ParEnableNormalization(const ParGridFunction &x)
Definition: tmop.cpp:3490
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:900
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:135
DenseMatrix P
Definition: tmop.hpp:1376
MPI_Comm GetComm() const
Definition: tmop.hpp:978
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:912
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field)=0
void UpdateAfterMeshTopologyChange()
Update all discrete fields based on tspec and update for AMR.
Definition: tmop.cpp:1516
void ResetDerefinementTspecData()
Definition: tmop.hpp:1283
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:772
bool Parallel() const
Definition: tmop.hpp:977
3D barrier Shape (S) metric.
Definition: tmop.hpp:450
void AssembleGradPA_C0_3D(const Vector &) const
2D non-barrier metric without a type.
Definition: tmop.hpp:108
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt).
Definition: tmop.cpp:892
FiniteElementSpace * fes
Definition: tmop.hpp:858
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:507
virtual ~TMOP_Metric_328()
Definition: tmop.hpp:613
FiniteElementSpace * coarse_tspec_fesv
Definition: tmop.hpp:1118
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:1303
2D non-barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:238