MFEM  v4.5.1
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 HyperelasticModel::SetTransformation() is hidden
30  for TMOP_QualityMetric%s, 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 Evaluates the metric in matrix form (opposed to invariant form).
46  Used for validating the invariant evaluations. */
47  virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
48  { return -1.0; /* not implemented -> checks would fail. */ }
49 
50  /** @brief Evaluate the strain energy density function, W = W(Jpt), by using
51  the 2D or 3D matrix invariants, see linalg/invariants.hpp.
52  @param[in] Jpt Represents the target->physical transformation
53  Jacobian matrix. */
54  virtual double EvalW(const DenseMatrix &Jpt) const = 0;
55 
56  /** @brief Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
57  @param[in] Jpt Represents the target->physical transformation
58  Jacobian matrix.
59  @param[out] P The evaluated 1st Piola-Kirchhoff stress tensor. */
60  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const = 0;
61 
62  /** @brief Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor
63  and assemble its contribution to the local gradient matrix 'A'.
64  @param[in] Jpt Represents the target->physical transformation
65  Jacobian matrix.
66  @param[in] DS Gradient of the basis matrix (dof x dim).
67  @param[in] weight Quadrature weight coefficient for the point.
68  @param[in,out] A Local gradient matrix where the contribution from this
69  point will be added.
70 
71  Computes weight * d(dW_dxi)_d(xj) at the current point, for all i and j,
72  where x1 ... xn are the FE dofs. This function is usually defined using
73  the matrix invariants and their derivatives. */
74  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
75  const double weight, DenseMatrix &A) const = 0;
76 
77  /** @brief Return the metric ID. */
78  virtual int Id() const { return 0; }
79 };
80 
81 /// Abstract class used to define combination of metrics with constant coefficients.
83 {
84 protected:
87 
88 public:
89  virtual void AddQualityMetric(TMOP_QualityMetric *tq, double wt = 1.0)
90  {
91  tmop_q_arr.Append(tq);
92  wt_arr.Append(wt);
93  }
94 
95  virtual void SetTargetJacobian(const DenseMatrix &Jtr_)
96  {
97  for (int i = 0; i < tmop_q_arr.Size(); i++)
98  {
99  tmop_q_arr[i]->SetTargetJacobian(Jtr_);
100  }
101  }
102 
103  virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const;
104 
105  virtual double EvalW(const DenseMatrix &Jpt) const;
106 
107  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
108 
109  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
110  const double weight, DenseMatrix &A) const;
111 };
112 
113 /// Simultaneous Untangler + Worst Case Improvement Metric
114 /// Uses a base metric mu and is defined as:
115 /// mu_tilde = mu_hat, when WorstCaseType = None,
116 /// = mu_hat/(beta - mu_hat), when WorstCaseType = Beta,
117 /// = mu_hat^p, when WorstCaseType = PMean,
118 /// where beta = max(mu_hat) + muT_ep,
119 /// and mu_hat = (mu/2phi(tau,ep)) where
120 /// 2phi(tau,ep) = 1, when when BarrierType = None,
121 /// = 2*(tau - min(alpha*min(tau)-detT_ep,0)), when BarrierType = Shifted
122 /// = tau^2 + sqrt(tau^2 + ep^2), when BarrierType = Pseudo
123 /// where tau = det(T), and max(mu_hat) and min(tau) are computed over the
124 /// entire mesh.
125 /// Ultimately, this metric can be used for mesh untangling with the BarrierType
126 /// option and for worst case quality improvement with the WorstCaseType option.
128 {
129 public:
130  enum class BarrierType
131  {
132  None,
133  Shifted,
134  Pseudo
135  };
136  enum class WorstCaseType
137  {
138  None,
139  Beta,
140  PMean
141  };
142 
143 protected:
144  TMOP_QualityMetric &tmop_metric; // non-barrier metric to use
145  double min_detT; // minimum Jacobian in the mesh
146  double max_muT; // max mu_k/phi(tau,ep) in the mesh
147  int exponent; // used for p-mean metrics
148  double alpha; // scaling factor for min(det(T))
149  double detT_ep; // small constant subtracted from min(detT)
150  double muT_ep; // small constant added to muT term
153 
154 public:
156  int exponent_ = 1,
157  double alpha_ = 1.5,
158  double detT_ep_ = 0.0001,
159  double muT_ep_ = 0.0001,
162  tmop_metric(tmop_metric_), exponent(exponent_), alpha(alpha_),
163  detT_ep(detT_ep_), muT_ep(muT_ep_), btype(btype_), wctype(wctype_)
164  {
165  MFEM_VERIFY(wctype == WorstCaseType::None,
166  "Worst-case optimization has not been fully developed!");
167  if (btype != BarrierType::None)
168  {
169  const int m_id = tmop_metric.Id();
170  MFEM_VERIFY(m_id == 4 || m_id == 14 || m_id == 66,
171  "Incorrect input barrier metric -- must be 4 / 14 / 66");
172  }
173  }
174 
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  // Compute mu_hat.
185  virtual double EvalWBarrier(const DenseMatrix &Jpt) const;
186 
187  virtual void SetMinDetT(double min_detT_) { min_detT = min_detT_; }
188 
189  virtual void SetMaxMuT(double max_muT_) { max_muT = max_muT_; }
190 
191  virtual BarrierType GetBarrierType() { return btype; }
192 
193  virtual WorstCaseType GetWorstCaseType() { return wctype; }
194 };
195 
196 /// 2D non-barrier metric without a type.
198 {
199 protected:
201 
202 public:
203  // W = |J|^2.
204  virtual double EvalW(const DenseMatrix &Jpt) const;
205 
206  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
207 
208  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
209  const double weight, DenseMatrix &A) const;
210 
211  virtual int Id() const { return 1; }
212 };
213 
214 /// 2D non-barrier Skew metric.
216 {
217 public:
218  // W = 0.5 (1 - cos(angle_Jpr - angle_Jtr)).
219  virtual double EvalW(const DenseMatrix &Jpt) const;
220 
221  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
222  { MFEM_ABORT("Not implemented"); }
223 
224  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
225  const double weight, DenseMatrix &A) const
226  { MFEM_ABORT("Not implemented"); }
227 };
228 
229 /// 3D non-barrier Skew metric.
231 {
232 public:
233  // W = 1/6 (3 - sum_i cos(angle_Jpr_i - angle_Jtr_i)), i = 1..3.
234  virtual double EvalW(const DenseMatrix &Jpt) const;
235 
236  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
237  { MFEM_ABORT("Not implemented"); }
238 
239  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
240  const double weight, DenseMatrix &A) const
241  { MFEM_ABORT("Not implemented"); }
242 };
243 
244 /// 2D non-barrier Aspect ratio metric.
246 {
247 public:
248  // W = 0.5 (ar_Jpr/ar_Jtr + ar_Jtr/ar_Jpr) - 1.
249  virtual double EvalW(const DenseMatrix &Jpt) const;
250 
251  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
252  { MFEM_ABORT("Not implemented"); }
253 
254  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
255  const double weight, DenseMatrix &A) const
256  { MFEM_ABORT("Not implemented"); }
257 };
258 
259 /// 3D non-barrier Aspect ratio metric.
261 {
262 public:
263  // W = 1/3 sum [0.5 (ar_Jpr_i/ar_Jtr_i + ar_Jtr_i/ar_Jpr_i) - 1], i = 1..3.
264  virtual double EvalW(const DenseMatrix &Jpt) const;
265 
266  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
267  { MFEM_ABORT("Not implemented"); }
268 
269  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
270  const double weight, DenseMatrix &A) const
271  { MFEM_ABORT("Not implemented"); }
272 };
273 
274 /// 2D barrier shape (S) metric (polyconvex).
276 {
277 protected:
279 
280 public:
281  // W = 0.5 |J|^2 / det(J) - 1.
282  virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const;
283 
284  // W = 0.5 I1b - 1.
285  virtual double EvalW(const DenseMatrix &Jpt) const;
286 
287  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
288 
289  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
290  const double weight, DenseMatrix &A) const;
291 
292  virtual int Id() const { return 2; }
293 };
294 
295 /// 2D non-barrier shape (S) metric.
297 {
298 protected:
300 
301 public:
302  // W = |J|^2 - 2*det(J)
303  virtual double EvalW(const DenseMatrix &Jpt) const;
304 
305  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
306 
307  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
308  const double weight, DenseMatrix &A) const;
309 
310  virtual int Id() const { return 4; }
311 };
312 
313 /// 2D barrier Shape+Size (VS) metric (not polyconvex).
315 {
316 protected:
318 
319 public:
320  // W = |J - J^-t|^2.
321  virtual double EvalW(const DenseMatrix &Jpt) const;
322 
323  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
324 
325  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
326  const double weight, DenseMatrix &A) const;
327 
328  virtual int Id() const { return 7; }
329 };
330 
331 /// 2D barrier Shape+Size (VS) metric (not polyconvex).
333 {
334 protected:
336 
337 public:
338  // W = det(J) * |J - J^-t|^2.
339  virtual double EvalW(const DenseMatrix &Jpt) const;
340 
341  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
342 
343  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
344  const double weight, DenseMatrix &A) const;
345 };
346 
347 /// 2D non-barrier Shape+Size+Orientation (VOS) metric (polyconvex).
349 {
350 public:
351  // W = |T-I|^2.
352  virtual double EvalW(const DenseMatrix &Jpt) const;
353 
354  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
355  { MFEM_ABORT("Not implemented"); }
356 
357  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
358  const double weight, DenseMatrix &A) const
359  { MFEM_ABORT("Not implemented"); }
360 };
361 
362 /// 2D Shifted barrier form of shape metric (mu_2).
364 {
365 protected:
366  double &min_detT;
368 
369 public:
370  TMOP_Metric_022(double &t0): min_detT(t0) {}
371 
372  // W = 0.5(|J|^2 - 2det(J)) / (det(J) - tau0).
373  virtual double EvalW(const DenseMatrix &Jpt) const;
374 
375  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
376 
377  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
378  const double weight, DenseMatrix &A) const;
379 };
380 
381 /// 2D barrier (not a shape) metric (polyconvex).
383 {
384 protected:
386 
387 public:
388  // W = 0.5|J^t J|^2 / det(J)^2 - 1.
389  virtual double EvalW(const DenseMatrix &Jpt) const;
390 
391  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
392 
393  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
394  const double weight, DenseMatrix &A) const;
395 };
396 
397 /// 2D non-barrier size (V) metric (not polyconvex).
399 {
400 protected:
402 
403 public:
404  // W = (det(J) - 1)^2.
405  virtual double EvalW(const DenseMatrix &Jpt) const;
406 
407  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
408 
409  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
410  const double weight, DenseMatrix &A) const;
411 
412 };
413 
414 /// 2D barrier size (V) metric (polyconvex).
416 {
417 protected:
419 
420 public:
421  // W = 0.5( sqrt(det(J)) - 1 / sqrt(det(J)) )^2
422  // = 0.5( det(J) - 1 )^2 / det(J)
423  // = 0.5( det(J) + 1/det(J) ) - 1.
424  virtual double EvalW(const DenseMatrix &Jpt) const;
425 
426  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
427 
428  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
429  const double weight, DenseMatrix &A) const;
430 };
431 
432 /// 2D barrier shape (S) metric (not polyconvex).
434 {
435 protected:
437 
438 public:
439  // W = |J^t J|^2 / det(J)^2 - 2|J|^2 / det(J) + 2
440  virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const;
441 
442  // W = I1b (I1b - 2).
443  virtual double EvalW(const DenseMatrix &Jpt) const;
444 
445  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
446 
447  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
448  const double weight, DenseMatrix &A) const;
449 };
450 
451 /// 2D non-barrier Shape+Size (VS) metric.
453 {
454 protected:
456  double gamma;
458 
459 public:
460  TMOP_Metric_066(double gamma_) : gamma(gamma_),
463  {
464  // (1-gamma) mu_4 + gamma mu_55
465  AddQualityMetric(sh_metric, 1.-gamma_);
466  AddQualityMetric(sz_metric, gamma_);
467  }
468  virtual int Id() const { return 66; }
469  double GetGamma() const { return gamma; }
470 
471  virtual ~TMOP_Metric_066() { delete sh_metric; delete sz_metric; }
472 };
473 
474 /// 2D barrier size (V) metric (polyconvex).
476 {
477 protected:
479 
480 public:
481  // W = 0.5(det(J) - 1 / det(J))^2.
482  virtual double EvalW(const DenseMatrix &Jpt) const;
483 
484  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
485 
486  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
487  const double weight, DenseMatrix &A) const;
488 
489  virtual int Id() const { return 77; }
490 };
491 
492 /// 2D barrier Shape+Size (VS) metric (polyconvex).
494 {
495 protected:
497  double gamma;
499 
500 public:
501  TMOP_Metric_080(double gamma_) : gamma(gamma_),
504  {
505  // (1-gamma) mu_2 + gamma mu_77
506  AddQualityMetric(sh_metric, 1.-gamma_);
507  AddQualityMetric(sz_metric, gamma_);
508  }
509  virtual int Id() const { return 80; }
510  double GetGamma() const { return gamma; }
511 
512  virtual ~TMOP_Metric_080() { delete sh_metric; delete sz_metric; }
513 };
514 
515 /// 2D barrier Shape+Orientation (OS) metric (polyconvex).
517 {
518 public:
519  // W = |T-T'|^2, where T'= |T|*I/sqrt(2).
520  virtual double EvalW(const DenseMatrix &Jpt) const;
521 
522  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
523  { MFEM_ABORT("Not implemented"); }
524 
525  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
526  const double weight, DenseMatrix &A) const
527  { MFEM_ABORT("Not implemented"); }
528 };
529 
530 /// 2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
532 {
533 public:
534  // W = 1/tau |T-I|^2.
535  virtual double EvalW(const DenseMatrix &Jpt) const;
536 
537  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
538  { MFEM_ABORT("Not implemented"); }
539 
540  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
541  const double weight, DenseMatrix &A) const
542  { MFEM_ABORT("Not implemented"); }
543 };
544 
545 /// 2D untangling metric.
547 {
548 protected:
549  const double eps;
551 
552 public:
553  TMOP_Metric_211(double epsilon = 1e-4) : eps(epsilon) { }
554 
555  // W = (det(J) - 1)^2 - det(J) + sqrt(det(J)^2 + eps).
556  virtual double EvalW(const DenseMatrix &Jpt) const;
557 
558  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
559 
560  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
561  const double weight, DenseMatrix &A) const;
562 };
563 
564 /// Shifted barrier form of metric 56 (area, ideal barrier metric), 2D
566 {
567 protected:
568  double &tau0;
570 
571 public:
572  /// Note that @a t0 is stored by reference
573  TMOP_Metric_252(double &t0): tau0(t0) {}
574 
575  // W = 0.5(det(J) - 1)^2 / (det(J) - tau0).
576  virtual double EvalW(const DenseMatrix &Jpt) const;
577 
578  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
579 
580  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
581  const double weight, DenseMatrix &A) const;
582 };
583 
584 /// 3D barrier Shape (S) metric, well-posed (polyconvex & invex).
586 {
587 protected:
589 
590 public:
591  // W = 1/3 |J| |J^-1| - 1.
592  virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const;
593 
594  // W = 1/3 sqrt(I1b * I2b) - 1
595  virtual double EvalW(const DenseMatrix &Jpt) const;
596 
597  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
598 
599  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
600  const double weight, DenseMatrix &A) const;
601 };
602 
603 /// 3D barrier Shape (S) metric, well-posed (polyconvex & invex).
605 {
606 protected:
608 
609 public:
610  // W = |J|^2 |J^{-1}|^2 / 9 - 1.
611  virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const;
612 
613  // W = I1b * I2b / 9 - 1.
614  virtual double EvalW(const DenseMatrix &Jpt) const;
615 
616  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
617 
618  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
619  const double weight, DenseMatrix &A) const;
620 
621  virtual int Id() const { return 302; }
622 };
623 
624 /// 3D barrier Shape (S) metric, well-posed (polyconvex & invex).
626 {
627 protected:
629 
630 public:
631  // W = |J|^2 / 3 / det(J)^(2/3) - 1.
632  virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const;
633 
634  // W = I1b / 3 - 1.
635  virtual double EvalW(const DenseMatrix &Jpt) const;
636 
637  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
638 
639  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
640  const double weight, DenseMatrix &A) const;
641 
642  virtual int Id() const { return 303; }
643 };
644 
645 /// 3D barrier Shape (S) metric, well-posed (polyconvex & invex).
647 {
648 protected:
650 
651 public:
652  // W = |J|^3 / 3^(3/2) / det(J) - 1.
653  virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const;
654 
655  // W = (I1b/3)^3/2 - 1.
656  virtual double EvalW(const DenseMatrix &Jpt) const;
657 
658  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
659 
660  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
661  const double weight, DenseMatrix &A) const;
662 
663  virtual int Id() const { return 304; }
664 };
665 
666 /// 3D Size (V) untangling metric.
668 {
669 protected:
670  const double eps;
672 
673 public:
674  TMOP_Metric_311(double epsilon = 1e-4) : eps(epsilon) { }
675 
676  // W = (det(J) - 1)^2 - det(J) + (det(J)^2 + eps)^(1/2).
677  virtual double EvalW(const DenseMatrix &Jpt) const;
678 
679  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
680 
681  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
682  const double weight, DenseMatrix &A) const;
683 };
684 
685 /// 3D Shape (S) metric, untangling version of 303.
687 {
688 protected:
689  double &min_detT;
691 
692 public:
693  TMOP_Metric_313(double &mindet) : min_detT(mindet) { }
694 
695  // W = 1/3 |J|^2 / [det(J)-tau0]^(-2/3).
696  virtual double EvalW(const DenseMatrix &Jpt) const;
697 
698  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
699 
700  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
701  const double weight, DenseMatrix &A) const;
702 
703  virtual int Id() const { return 313; }
704 };
705 
706 /// 3D non-barrier metric without a type.
708 {
709 protected:
711 
712 public:
713  // W = (det(J) - 1)^2.
714  virtual double EvalW(const DenseMatrix &Jpt) const;
715 
716  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
717 
718  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
719  const double weight, DenseMatrix &A) const;
720 
721  virtual int Id() const { return 315; }
722 };
723 
724 /// 3D barrier metric without a type.
726 {
727 protected:
729 
730 public:
731  // W = 0.5 (det(J) + 1/det(J)) - 1.
732  virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const;
733 
734  // W = 0.5 (I3b + 1/I3b) - 1.
735  virtual double EvalW(const DenseMatrix &Jpt) const;
736 
737  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
738 
739  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
740  const double weight, DenseMatrix &A) const;
741 };
742 
743 /// 3D barrier Shape+Size (VS) metric, well-posed (invex).
745 {
746 protected:
748 
749 public:
750  // W = |J - J^-t|^2.
751  virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const;
752 
753  // W = I1 + I2/I3 - 6.
754  virtual double EvalW(const DenseMatrix &Jpt) const;
755 
756  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
757 
758  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
759  const double weight, DenseMatrix &A) const;
760 
761  virtual int Id() const { return 321; }
762 };
763 
764 /// 3D barrier Shape+Size (VS) metric, well-posed (invex).
766 {
767 protected:
769 
770 public:
771  // W = |J - adjJ^-t|^2.
772  virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const;
773 
774  // W = I1b / (I3b^-1/3) / 6 + I2b (I3b^1/3) / 6 - 1
775  virtual double EvalW(const DenseMatrix &Jpt) const;
776 
777  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
778 
779  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
780  const double weight, DenseMatrix &A) const;
781 
782  virtual int Id() const { return 322; }
783 };
784 
785 /// 3D barrier Shape+Size (VS) metric, well-posed (invex).
787 {
788 protected:
790 
791 public:
792  // W = |J|^3 - 3 sqrt(3) ln(det(J)) - 3 sqrt(3).
793  virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const;
794 
795  // W = I1^3/2 - 3 sqrt(3) ln(I3b) - 3 sqrt(3).
796  virtual double EvalW(const DenseMatrix &Jpt) const;
797 
798  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
799 
800  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
801  const double weight, DenseMatrix &A) const;
802 
803  virtual int Id() const { return 323; }
804 };
805 
806 /// 3D barrier Shape+Size (VS) metric (polyconvex).
808 {
809 protected:
811  double gamma;
813 
814 public:
815  TMOP_Metric_328(double gamma_) : gamma(gamma_),
818  {
819  // (1-gamma) mu_301 + gamma mu_316
820  AddQualityMetric(sh_metric, 1.-gamma_);
821  AddQualityMetric(sz_metric, gamma_);
822  }
823 
824  virtual ~TMOP_Metric_328() { delete sh_metric; delete sz_metric; }
825 };
826 
827 /// 3D barrier Shape+Size (VS) metric (polyconvex).
829 {
830 protected:
831  double gamma;
833 
834 public:
835  TMOP_Metric_332(double gamma_) : gamma(gamma_),
838  {
839  // (1-gamma) mu_302 + gamma mu_315
840  AddQualityMetric(sh_metric, 1.-gamma_);
841  AddQualityMetric(sz_metric, gamma_);
842  }
843 
844  virtual int Id() const { return 332; }
845  double GetGamma() const { return gamma; }
846 
847  virtual ~TMOP_Metric_332() { delete sh_metric; delete sz_metric; }
848 };
849 
850 /// 3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
852 {
853 protected:
855  double gamma;
857 
858 public:
859  TMOP_Metric_333(double gamma_) : gamma(gamma_),
862  {
863  // (1-gamma) mu_302 + gamma mu_316
864  AddQualityMetric(sh_metric, 1.-gamma_);
865  AddQualityMetric(sz_metric, gamma_);
866  }
867 
868  virtual ~TMOP_Metric_333() { delete sh_metric; delete sz_metric; }
869 };
870 
871 /// 3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
873 {
874 protected:
876  double gamma;
878 
879 public:
880  TMOP_Metric_334(double gamma_) : gamma(gamma_),
883  {
884  // (1-gamma) mu_303 + gamma mu_316
885  AddQualityMetric(sh_metric, 1.-gamma_);
886  AddQualityMetric(sz_metric, gamma_);
887  }
888 
889  virtual int Id() const { return 334; }
890  double GetGamma() const { return gamma; }
891 
892  virtual ~TMOP_Metric_334() { delete sh_metric; delete sz_metric; }
893 };
894 
895 /// 3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
897 {
898 protected:
900  double gamma;
902 
903 public:
904  TMOP_Metric_347(double gamma_) : gamma(gamma_),
907  {
908  // (1-gamma) mu_304 + gamma mu_316
909  AddQualityMetric(sh_metric, 1.-gamma_);
910  AddQualityMetric(sz_metric, gamma_);
911  }
912 
913  virtual int Id() const { return 347; }
914  double GetGamma() const { return gamma; }
915 
916  virtual ~TMOP_Metric_347() { delete sh_metric; delete sz_metric; }
917 };
918 
919 /// 3D shifted barrier form of metric 316 (not typed).
921 {
922 protected:
923  double &tau0;
925 
926 public:
927  TMOP_Metric_352(double &t0): tau0(t0) {}
928 
929  // W = 0.5(det(J) - 1)^2 / (det(J) - tau0).
930  virtual double EvalW(const DenseMatrix &Jpt) const;
931 
932  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
933 
934  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
935  const double weight, DenseMatrix &A) const;
936 };
937 
938 /// 3D non-barrier Shape (S) metric.
940 {
941 protected:
943 
944 public:
945  // W = |J|^3 / 3^(3/2) - det(J).
946  virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const;
947 
948  // W = (I1b/3)^3/2 - 1.
949  virtual double EvalW(const DenseMatrix &Jpt) const;
950 
951  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const;
952 
953  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
954  const double weight, DenseMatrix &A) const;
955 
956  virtual int Id() const { return 360; }
957 };
958 
959 /// A-metrics
960 /// 2D barrier Shape (S) metric (polyconvex).
962 {
963 protected:
965 
966 public:
967  // (1/4 alpha) | A - (adj A)^t W^t W / omega |^2
968  virtual double EvalW(const DenseMatrix &Jpt) const;
969 
970  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
971  { MFEM_ABORT("Not implemented"); }
972 
973  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
974  const double weight, DenseMatrix &A) const
975  { MFEM_ABORT("Not implemented"); }
976 };
977 
978 /// 2D barrier Size (V) metric (polyconvex).
980 {
981 protected:
983 
984 public:
985  // 0.5 * ( sqrt(alpha/omega) - sqrt(omega/alpha) )^2
986  virtual double EvalW(const DenseMatrix &Jpt) const;
987 
988  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
989  { MFEM_ABORT("Not implemented"); }
990 
991  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
992  const double weight, DenseMatrix &A) const
993  { MFEM_ABORT("Not implemented"); }
994 };
995 
996 /// 2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
998 {
999 protected:
1001 
1002 public:
1003  // (1/alpha) | A - W |^2
1004  virtual double EvalW(const DenseMatrix &Jpt) const;
1005 
1006  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
1007  { MFEM_ABORT("Not implemented"); }
1008 
1009  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
1010  const double weight, DenseMatrix &A) const
1011  { MFEM_ABORT("Not implemented"); }
1012 };
1013 
1014 /// 2D barrier Shape+Orientation (OS) metric (polyconvex).
1016 {
1017 protected:
1019 
1020 public:
1021  // (1/2 alpha) | A - (|A|/|W|) W |^2
1022  virtual double EvalW(const DenseMatrix &Jpt) const;
1023 
1024  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
1025  { MFEM_ABORT("Not implemented"); }
1026 
1027  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
1028  const double weight, DenseMatrix &A) const
1029  { MFEM_ABORT("Not implemented"); }
1030 };
1031 
1032 /// 2D barrier Shape+Size (VS) metric (polyconvex).
1034 {
1035 protected:
1037  double gamma;
1039 
1040 public:
1041  TMOP_AMetric_126(double gamma_) : gamma(gamma_),
1044  {
1045  // (1-gamma) nu_11 + gamma nu_14
1046  AddQualityMetric(sh_metric, 1.-gamma_);
1047  AddQualityMetric(sz_metric, gamma_);
1048  }
1049 
1050  virtual ~TMOP_AMetric_126() { delete sh_metric; delete sz_metric; }
1051 };
1052 
1053 /// Base class for limiting functions to be used in class TMOP_Integrator.
1054 /** This class represents a scalar function f(x, x0, d), where x and x0 are
1055  positions in physical space, and d is a reference physical distance
1056  associated with the point x0. */
1058 {
1059 public:
1060  /// Returns the limiting function, f(x, x0, d).
1061  virtual double Eval(const Vector &x, const Vector &x0, double d) const = 0;
1062 
1063  /** @brief Returns the gradient of the limiting function f(x, x0, d) with
1064  respect to x. */
1065  virtual void Eval_d1(const Vector &x, const Vector &x0, double dist,
1066  Vector &d1) const = 0;
1067 
1068  /** @brief Returns the Hessian of the limiting function f(x, x0, d) with
1069  respect to x. */
1070  virtual void Eval_d2(const Vector &x, const Vector &x0, double dist,
1071  DenseMatrix &d2) const = 0;
1072 
1073  /// Virtual destructor.
1075 };
1076 
1077 /// Default limiter function in TMOP_Integrator.
1079 {
1080 public:
1081  virtual double Eval(const Vector &x, const Vector &x0, double dist) const
1082  {
1083  MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
1084 
1085  return 0.5 * x.DistanceSquaredTo(x0) / (dist * dist);
1086  }
1087 
1088  virtual void Eval_d1(const Vector &x, const Vector &x0, double dist,
1089  Vector &d1) const
1090  {
1091  MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
1092 
1093  d1.SetSize(x.Size());
1094  subtract(1.0 / (dist * dist), x, x0, d1);
1095  }
1096 
1097  virtual void Eval_d2(const Vector &x, const Vector &x0, double dist,
1098  DenseMatrix &d2) const
1099  {
1100  MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
1101 
1102  d2.Diag(1.0 / (dist * dist), x.Size());
1103  }
1104 
1106 };
1107 
1108 class FiniteElementCollection;
1109 class FiniteElementSpace;
1110 class ParFiniteElementSpace;
1111 
1113 {
1114 protected:
1115  // Owned.
1118 
1119 #ifdef MFEM_USE_MPI
1120  // Owned.
1123 #endif
1124 
1125 public:
1126  AdaptivityEvaluator() : mesh(NULL), fes(NULL)
1127  {
1128 #ifdef MFEM_USE_MPI
1129  pmesh = NULL;
1130  pfes = NULL;
1131 #endif
1132  }
1133  virtual ~AdaptivityEvaluator();
1134 
1135  /** Specifies the Mesh and FiniteElementSpace of the solution that will
1136  be evaluated. The given mesh will be copied into the internal object. */
1137  void SetSerialMetaInfo(const Mesh &m,
1138  const FiniteElementSpace &f);
1139 
1140 #ifdef MFEM_USE_MPI
1141  /// Parallel version of SetSerialMetaInfo.
1142  void SetParMetaInfo(const ParMesh &m,
1143  const ParFiniteElementSpace &f);
1144 #endif
1145 
1146  // TODO use GridFunctions to make clear it's on the ldofs?
1147  virtual void SetInitialField(const Vector &init_nodes,
1148  const Vector &init_field) = 0;
1149 
1150  virtual void ComputeAtNewPosition(const Vector &new_nodes,
1151  Vector &new_field,
1152  int new_nodes_ordering = Ordering::byNODES) = 0;
1153 
1154  void ClearGeometricFactors();
1155 };
1156 
1157 /** @brief Base class representing target-matrix construction algorithms for
1158  mesh optimization via the target-matrix optimization paradigm (TMOP). */
1159 /** This class is used by class TMOP_Integrator to construct the target Jacobian
1160  matrices (reference-element to target-element) at quadrature points. It
1161  supports a set of algorithms chosen by the #TargetType enumeration.
1162 
1163  New target-matrix construction algorithms can be defined by deriving new
1164  classes and overriding the methods ComputeElementTargets() and
1165  ContainsVolumeInfo(). */
1167 {
1168 public:
1169  /// Target-matrix construction algorithms supported by this class.
1171  {
1173  Ideal shape, unit size; the nodes are not used. */
1175  Ideal shape, equal size/volume; the given nodes define the total target
1176  volume; for each mesh element, the target volume is the average volume
1177  multiplied by the volume scale, set with SetVolumeScale(). */
1179  Ideal shape, given size/volume; the given nodes define the target
1180  volume at all quadrature points. */
1182  Given shape, given size/volume; the given nodes define the exact target
1183  Jacobian matrix at all quadrature points. */
1185  Full target tensor is specified at every quadrature point. */
1186  };
1187 
1188 protected:
1189  // Nodes that are used in ComputeElementTargets(), depending on target_type.
1190  const GridFunction *nodes; // not owned
1191  mutable double avg_volume;
1194  bool uses_phys_coords; // see UsesPhysicalCoordinates()
1195 
1196 #ifdef MFEM_USE_MPI
1197  MPI_Comm comm;
1198 #endif
1199 
1200  // should be called only if avg_volume == 0.0, i.e. avg_volume is not
1201  // computed yet
1202  void ComputeAvgVolume() const;
1203 
1204  template<int DIM>
1206  const IntegrationRule &ir,
1207  const Vector &xe,
1208  DenseTensor &Jtr) const;
1209 
1210  // CPU fallback that uses ComputeElementTargets()
1212  const IntegrationRule &ir,
1213  const Vector &xe,
1214  DenseTensor &Jtr) const;
1215 
1216 public:
1217  /// Constructor for use in serial
1219  : nodes(NULL), avg_volume(), volume_scale(1.0), target_type(ttype),
1220  uses_phys_coords(false)
1221  {
1222 #ifdef MFEM_USE_MPI
1223  comm = MPI_COMM_NULL;
1224 #endif
1225  }
1226 #ifdef MFEM_USE_MPI
1227  /// Constructor for use in parallel
1228  TargetConstructor(TargetType ttype, MPI_Comm mpicomm)
1229  : nodes(NULL), avg_volume(), volume_scale(1.0), target_type(ttype),
1230  uses_phys_coords(false), comm(mpicomm) { }
1231 #endif
1232  virtual ~TargetConstructor() { }
1233 
1234 #ifdef MFEM_USE_MPI
1235  bool Parallel() const { return (comm != MPI_COMM_NULL); }
1236  MPI_Comm GetComm() const { return comm; }
1237 #else
1238  bool Parallel() const { return false; }
1239 #endif
1240 
1241  /** @brief Set the nodes to be used in the target-matrix construction.
1242 
1243  This method should be called every time the target nodes are updated
1244  externally and recomputation of the target average volume is needed. The
1245  nodes are used by all target types except IDEAL_SHAPE_UNIT_SIZE. */
1246  void SetNodes(const GridFunction &n) { nodes = &n; avg_volume = 0.0; }
1247 
1248  /** @brief Get the nodes to be used in the target-matrix construction. */
1249  const GridFunction *GetNodes() const { return nodes; }
1250 
1251  /// Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
1252  void SetVolumeScale(double vol_scale) { volume_scale = vol_scale; }
1253 
1255 
1256  /** @brief Return true if the methods ComputeElementTargets(),
1257  ComputeAllElementTargets(), and ComputeElementTargetsGradient() use the
1258  physical node coordinates provided by the parameters 'elfun', or 'xe'. */
1260 
1261  /// Checks if the target matrices contain non-trivial size specification.
1262  virtual bool ContainsVolumeInfo() const;
1263 
1264  /** @brief Given an element and quadrature rule, computes ref->target
1265  transformation Jacobians for each quadrature point in the element.
1266  The physical positions of the element's nodes are given by @a elfun. */
1267  virtual void ComputeElementTargets(int e_id, const FiniteElement &fe,
1268  const IntegrationRule &ir,
1269  const Vector &elfun,
1270  DenseTensor &Jtr) const;
1271 
1272  /** @brief Computes reference-to-target transformation Jacobians for all
1273  quadrature points in all elements.
1274 
1275  @param[in] fes The nodal FE space
1276  @param[in] ir The quadrature rule to use for all elements
1277  @param[in] xe E-vector with the current physical coordinates/positions;
1278  this parameter is used only when needed by the target
1279  constructor, see UsesPhysicalCoordinates()
1280  @param[out] Jtr The computed ref->target Jacobian matrices. */
1281  virtual void ComputeAllElementTargets(const FiniteElementSpace &fes,
1282  const IntegrationRule &ir,
1283  const Vector &xe,
1284  DenseTensor &Jtr) const;
1285 
1286  virtual void ComputeElementTargetsGradient(const IntegrationRule &ir,
1287  const Vector &elfun,
1289  DenseTensor &dJtr) const;
1290 };
1291 
1293 {
1294 public:
1295  explicit TMOPMatrixCoefficient(int dim) : MatrixCoefficient(dim, dim) { }
1296 
1297  /** @brief Evaluate the derivative of the matrix coefficient with respect to
1298  @a comp in the element described by @a T at the point @a ip, storing the
1299  result in @a K. */
1300  virtual void EvalGrad(DenseMatrix &K, ElementTransformation &T,
1301  const IntegrationPoint &ip, int comp) = 0;
1302 
1304 };
1305 
1307 {
1308 protected:
1309  // Analytic target specification.
1313 
1314 public:
1316  : TargetConstructor(ttype),
1317  scalar_tspec(NULL), vector_tspec(NULL), matrix_tspec(NULL)
1318  { uses_phys_coords = true; }
1319 
1320  virtual void SetAnalyticTargetSpec(Coefficient *sspec,
1321  VectorCoefficient *vspec,
1322  TMOPMatrixCoefficient *mspec);
1323 
1324  /** @brief Given an element and quadrature rule, computes ref->target
1325  transformation Jacobians for each quadrature point in the element.
1326  The physical positions of the element's nodes are given by @a elfun. */
1327  virtual void ComputeElementTargets(int e_id, const FiniteElement &fe,
1328  const IntegrationRule &ir,
1329  const Vector &elfun,
1330  DenseTensor &Jtr) const;
1331 
1332  virtual void ComputeAllElementTargets(const FiniteElementSpace &fes,
1333  const IntegrationRule &ir,
1334  const Vector &xe,
1335  DenseTensor &Jtr) const;
1336 
1337  virtual void ComputeElementTargetsGradient(const IntegrationRule &ir,
1338  const Vector &elfun,
1340  DenseTensor &dJtr) const;
1341 };
1342 
1343 #ifdef MFEM_USE_MPI
1344 class ParGridFunction;
1345 #endif
1346 
1348 {
1349 protected:
1350  // Discrete target specification.
1351  // Data is owned, updated by UpdateTargetSpecification.
1353  Vector tspec; //eta(x) - we enforce Ordering::byNODES
1355  Vector tspec_pert1h; //eta(x+h)
1356  Vector tspec_pert2h; //eta(x+2*h)
1357  Vector tspec_pertmix; //eta(x+h,y+h)
1358  // The order inside these perturbation vectors (e.g. in 2D) is
1359  // eta1(x+h,y), eta2(x+h,y) ... etan(x+h,y), eta1(x,y+h), eta2(x,y+h) ...
1360  // same for tspec_pert2h and tspec_pertmix.
1361 
1362  // DenseMatrix to hold target_spec values for the (children of the)
1363  // element being refined to consider for h-refinement.
1365  // Vector to hold the target_spec values for the coarse version of the
1366  // current mesh. Used for derefinement decision with hr-adaptivity.
1368 
1369  // Components of Target Jacobian at each quadrature point of an element. This
1370  // is required for computation of the derivative using chain rule.
1372 
1373  // Note: do not use the Nodes of this space as they may not be on the
1374  // positions corresponding to the values of tspec.
1376  FiniteElementSpace *coarse_tspec_fesv; //not owned, derefinement FESpace
1377  GridFunction *tspec_gf; //owned, uses tspec and tspec_fes
1378  // discrete adaptivity
1379 #ifdef MFEM_USE_MPI
1380  ParFiniteElementSpace *ptspec_fesv; //owned, needed for derefinement to
1381  // get update operator.
1382  ParGridFunction *tspec_pgf; // similar to tspec_gf
1383 #endif
1384 
1385  int amr_el;
1387 
1388  // These flags can be used by outside functions to avoid recomputing the
1389  // tspec and tspec_perth fields again on the same mesh.
1391 
1392  // Evaluation of the discrete target specification on different meshes.
1393  // Owned.
1395 
1396  void SetDiscreteTargetBase(const GridFunction &tspec_);
1397  void SetTspecAtIndex(int idx, const GridFunction &tspec_);
1398  void FinalizeSerialDiscreteTargetSpec(const GridFunction &tspec_);
1399 #ifdef MFEM_USE_MPI
1400  void SetTspecAtIndex(int idx, const ParGridFunction &tspec_);
1401  void FinalizeParDiscreteTargetSpec(const ParGridFunction &tspec_);
1402 #endif
1403 
1404 public:
1406  : TargetConstructor(ttype),
1407  ncomp(0),
1408  sizeidx(-1), skewidx(-1), aspectratioidx(-1), orientationidx(-1),
1411  tspec_fesv(NULL), coarse_tspec_fesv(NULL), tspec_gf(NULL),
1412 #ifdef MFEM_USE_MPI
1413  ptspec_fesv(NULL), tspec_pgf(NULL),
1414 #endif
1415  amr_el(-1), lim_min_size(-0.1),
1416  good_tspec(false), good_tspec_grad(false), good_tspec_hess(false),
1417  adapt_eval(NULL) { }
1418 
1419  virtual ~DiscreteAdaptTC();
1420 
1421  /** @name Target specification methods.
1422  The following methods are used to specify geometric parameters of the
1423  targets when these parameters are given by discrete FE functions.
1424  Note that every GridFunction given to the Set methods must use a
1425  H1_FECollection of the same order. The number of components must
1426  correspond to the type of geometric parameter and dimension.
1427 
1428  @param[in] tspec_ Input values of a geometric parameter. Note that
1429  the methods in this class support only functions that
1430  use H1_FECollection collection of the same order. */
1431  ///@{
1432  virtual void SetSerialDiscreteTargetSpec(const GridFunction &tspec_);
1433  virtual void SetSerialDiscreteTargetSize(const GridFunction &tspec_);
1434  virtual void SetSerialDiscreteTargetSkew(const GridFunction &tspec_);
1435  virtual void SetSerialDiscreteTargetAspectRatio(const GridFunction &tspec_);
1436  virtual void SetSerialDiscreteTargetOrientation(const GridFunction &tspec_);
1437 #ifdef MFEM_USE_MPI
1438  virtual void SetParDiscreteTargetSpec(const ParGridFunction &tspec_);
1439  virtual void SetParDiscreteTargetSize(const ParGridFunction &tspec_);
1440  virtual void SetParDiscreteTargetSkew(const ParGridFunction &tspec_);
1441  virtual void SetParDiscreteTargetAspectRatio(const ParGridFunction &tspec_);
1442  virtual void SetParDiscreteTargetOrientation(const ParGridFunction &tspec_);
1443 #endif
1444  ///@}
1445 
1446  /// Used in combination with the Update methods to avoid extra computations.
1448  { good_tspec = good_tspec_grad = good_tspec_hess = false; }
1449 
1450  /// Get one of the discrete fields from tspec.
1451  void GetDiscreteTargetSpec(GridFunction &tspec_, int idx);
1452  /// Get the FESpace associated with tspec.
1454  /// Get the entire tspec.
1456  /// Update all discrete fields based on tspec and update for AMR
1458 
1459 #ifdef MFEM_USE_MPI
1462 #endif
1463 
1464  /** Used to update the target specification after the mesh has changed. The
1465  new mesh positions are given by new_x. If @a use_flags is true, repeated
1466  calls won't do anything until ResetUpdateFlags() is called. */
1467  void UpdateTargetSpecification(const Vector &new_x, bool use_flag = false,
1468  int new_x_ordering=Ordering::byNODES);
1469 
1470  void UpdateTargetSpecification(Vector &new_x, Vector &IntData,
1471  int new_x_ordering=Ordering::byNODES);
1472 
1475  int nodenum, int idir,
1476  const Vector &IntData);
1477 
1479 
1480  /** Used for finite-difference based computations. Computes the target
1481  specifications after a mesh perturbation in x or y direction.
1482  If @a use_flags is true, repeated calls won't do anything until
1483  ResetUpdateFlags() is called. */
1484  void UpdateGradientTargetSpecification(const Vector &x, double dx,
1485  bool use_flag = false,
1486  int x_ordering = Ordering::byNODES);
1487  /** Used for finite-difference based computations. Computes the target
1488  specifications after two mesh perturbations in x and/or y direction.
1489  If @a use_flags is true, repeated calls won't do anything until
1490  ResetUpdateFlags() is called. */
1491  void UpdateHessianTargetSpecification(const Vector &x, double dx,
1492  bool use_flag = false,
1493  int x_ordering = Ordering::byNODES);
1494 
1496  {
1497  if (adapt_eval) { delete adapt_eval; }
1498  adapt_eval = ae;
1499  }
1500 
1501  const Vector &GetTspecPert1H() { return tspec_pert1h; }
1502  const Vector &GetTspecPert2H() { return tspec_pert2h; }
1504 
1505  /** @brief Given an element and quadrature rule, computes ref->target
1506  transformation Jacobians for each quadrature point in the element.
1507  The physical positions of the element's nodes are given by @a elfun.
1508  Note that this function assumes that UpdateTargetSpecification() has
1509  been called with the position vector corresponding to @a elfun. */
1510  virtual void ComputeElementTargets(int e_id, const FiniteElement &fe,
1511  const IntegrationRule &ir,
1512  const Vector &elfun,
1513  DenseTensor &Jtr) const;
1514 
1515  virtual void ComputeAllElementTargets(const FiniteElementSpace &fes,
1516  const IntegrationRule &ir,
1517  const Vector &xe,
1518  DenseTensor &Jtr) const;
1519 
1520  virtual void ComputeElementTargetsGradient(const IntegrationRule &ir,
1521  const Vector &elfun,
1523  DenseTensor &dJtr) const;
1524 
1525  // Generates tspec_vals for target construction using intrule
1526  // Used for the refinement component in hr-adaptivity.
1527  void SetTspecFromIntRule(int e_id, const IntegrationRule &intrule);
1528 
1529  // Targets based on discrete functions can result in invalid (negative)
1530  // size at the quadrature points. This method can be used to set a
1531  // minimum target size.
1532  void SetMinSizeForTargets(double min_size_) { lim_min_size = min_size_; }
1533 
1534  /// Computes target specification data with respect to the coarse FE space.
1536 
1537  // Reset refinement data associated with h-adaptivity component.
1539  {
1540  tspec_refine.Clear();
1541  amr_el = -1;
1542  }
1543 
1544  // Reset derefinement data associated with h-adaptivity component.
1546  {
1548  coarse_tspec_fesv = NULL;
1549  }
1550 
1551  // Used to specify the fine element for determining energy of children of a
1552  // parent element.
1553  void SetRefinementSubElement(int amr_el_) { amr_el = amr_el_; }
1554 };
1555 
1556 class TMOPNewtonSolver;
1557 
1558 /** @brief A TMOP integrator class based on any given TMOP_QualityMetric and
1559  TargetConstructor.
1560 
1561  Represents @f$ \int W(Jpt) dx @f$ over a target zone, where W is the
1562  metric's strain energy density function, and Jpt is the Jacobian of the
1563  target->physical coordinates transformation. The virtual target zone is
1564  defined by the TargetConstructor. */
1566 {
1567 protected:
1568  friend class TMOPNewtonSolver;
1569  friend class TMOPComboIntegrator;
1570 
1573  const TargetConstructor *targetC; // not owned
1574 
1575  // Custom integration rules.
1578 
1579  // Weight Coefficient multiplying the quality metric term.
1580  Coefficient *metric_coeff; // not owned, if NULL -> metric_coeff is 1.
1581  // Normalization factor for the metric term.
1583 
1584  // Nodes and weight Coefficient used for "limiting" the TMOP_Integrator.
1585  // These are both NULL when there is no limiting.
1586  // The class doesn't own lim_nodes0 and lim_coeff.
1589  // Limiting reference distance. Not owned.
1591  // Limiting function. Owned.
1593  // Normalization factor for the limiting term.
1594  double lim_normal;
1595 
1596  // Adaptive limiting.
1597  const GridFunction *adapt_lim_gf0; // Not owned.
1598 #ifdef MFEM_USE_MPI
1600 #endif
1601  GridFunction *adapt_lim_gf; // Owned. Updated by adapt_lim_eval.
1604 
1605  // Surface fitting.
1606  GridFunction *surf_fit_gf; // Owned, Updated by surf_fit_eval.
1607  const Array<bool> *surf_fit_marker; // Not owned.
1608  Coefficient *surf_fit_coeff; // Not owned.
1611 
1613 
1614  // Parameters for FD-based Gradient & Hessian calculation.
1615  bool fdflag;
1616  double dx;
1617  double dxscale;
1618  // Specifies that ComputeElementTargets is being called by a FD function.
1619  // It's used to skip terms that have exact derivative calculations.
1621  // Compute the exact action of the Integrator (includes derivative of the
1622  // target with respect to spatial position)
1624 
1627 
1628  // Jrt: the inverse of the ref->target Jacobian, Jrt = Jtr^{-1}.
1629  // Jpr: the ref->physical transformation Jacobian, Jpr = PMatI^t DS.
1630  // Jpt: the target->physical transformation Jacobian, Jpt = Jpr Jrt.
1631  // P: represents dW_d(Jtp) (dim x dim).
1632  // DSh: gradients of reference shape functions (dof x dim).
1633  // DS: gradients of the shape functions in the target configuration,
1634  // DS = DSh Jrt (dof x dim).
1635  // PMatI: current coordinates of the nodes (dof x dim).
1636  // PMat0: reshaped view into the local element contribution to the operator
1637  // output - the result of AssembleElementVector() (dof x dim).
1639 
1640  // PA extension
1641  // ------------
1642  // E: Q-vector for TMOP-energy
1643  // O: Q-Vector of 1.0, used to compute sums using the dot product kernel.
1644  // X0: E-vector for initial nodal coordinates used for limiting.
1645  // H: Q-Vector for Hessian associated with the metric term.
1646  // C0: Q-Vector for spatial weight used for the limiting term.
1647  // LD: E-Vector constructed using limiting distance grid function (delta).
1648  // H0: Q-Vector for Hessian associated with the limiting term.
1649  //
1650  // maps: Dof2Quad map for fespace associate with nodal coordinates.
1651  // maps_lim: Dof2Quad map for fespace associated with the limiting distance
1652  // grid function.
1653  //
1654  // Jtr_debug_grad
1655  // We keep track if Jtr was set by AssembleGradPA() in Jtr_debug_grad: it
1656  // is set to true by AssembleGradPA(); any other call to
1657  // ComputeAllElementTargets() will set the flag to false. This flag will
1658  // be used to check that Jtr is the one set by AssembleGradPA() when
1659  // performing operations with the gradient like AddMultGradPA() and
1660  // AssembleGradDiagonalPA().
1661  //
1662  // TODO:
1663  // * Merge LD, C0, H0 into one scalar Q-vector
1664  struct
1665  {
1666  bool enabled;
1667  int dim, ne, nq;
1668  mutable DenseTensor Jtr;
1669  mutable bool Jtr_needs_update;
1670  mutable bool Jtr_debug_grad;
1671  mutable Vector E, O, X0, H, C0, LD, H0;
1672  const DofToQuad *maps;
1673  const DofToQuad *maps_lim = nullptr;
1677  } PA;
1678 
1680  double &metric_energy, double &lim_energy,
1681  double &surf_fit_gf_energy);
1682 
1685  const Vector &elfun, Vector &elvect);
1686 
1687  void AssembleElementGradExact(const FiniteElement &el,
1689  const Vector &elfun, DenseMatrix &elmat);
1690 
1691  void AssembleElementVectorFD(const FiniteElement &el,
1693  const Vector &elfun, Vector &elvect);
1694 
1695  // Assumes that AssembleElementVectorFD has been called.
1696  void AssembleElementGradFD(const FiniteElement &el,
1698  const Vector &elfun, DenseMatrix &elmat);
1699 
1700  void AssembleElemVecAdaptLim(const FiniteElement &el,
1702  const IntegrationRule &ir,
1703  const Vector &weights, DenseMatrix &mat);
1704  void AssembleElemGradAdaptLim(const FiniteElement &el,
1706  const IntegrationRule &ir,
1707  const Vector &weights, DenseMatrix &m);
1708 
1709  // First derivative of the surface fitting term.
1710  void AssembleElemVecSurfFit(const FiniteElement &el_x,
1712  DenseMatrix &mat);
1713 
1714  // Second derivative of the surface fitting term.
1715  void AssembleElemGradSurfFit(const FiniteElement &el_x,
1717  DenseMatrix &mat);
1718 
1719  double GetFDDerivative(const FiniteElement &el,
1721  Vector &elfun, const int nodenum,const int idir,
1722  const double baseenergy, bool update_stored);
1723 
1724  /** @brief Determines the perturbation, h, for FD-based approximation. */
1725  void ComputeFDh(const Vector &x, const FiniteElementSpace &fes);
1726  void ComputeMinJac(const Vector &x, const FiniteElementSpace &fes);
1727 
1728  void UpdateAfterMeshPositionChange(const Vector &new_x,
1729  int new_x_ordering = Ordering::byNODES);
1730 
1732  {
1733  lim_nodes0 = NULL; lim_coeff = NULL; lim_dist = NULL;
1734  delete lim_func; lim_func = NULL;
1735  }
1736 
1738  {
1739  if (IntegRules)
1740  {
1741  return IntegRules->Get(el.GetGeomType(), integ_order);
1742  }
1743  return (IntRule) ? *IntRule
1744  /* */ : IntRules.Get(el.GetGeomType(), 2*el.GetOrder() + 3);
1745  }
1747  {
1748  // TODO the energy most likely needs less integration points.
1749  return EnergyIntegrationRule(el);
1750  }
1752  {
1753  // TODO the action and energy most likely need less integration points.
1754  return EnergyIntegrationRule(el);
1755  }
1756 
1757  // Auxiliary PA methods
1758  void AssembleGradPA_2D(const Vector&) const;
1759  void AssembleGradPA_3D(const Vector&) const;
1760  void AssembleGradPA_C0_2D(const Vector&) const;
1761  void AssembleGradPA_C0_3D(const Vector&) const;
1762 
1763  double GetLocalStateEnergyPA_2D(const Vector&) const;
1764  double GetLocalStateEnergyPA_C0_2D(const Vector&) const;
1765  double GetLocalStateEnergyPA_3D(const Vector&) const;
1766  double GetLocalStateEnergyPA_C0_3D(const Vector&) const;
1767 
1768  void AddMultPA_2D(const Vector&, Vector&) const;
1769  void AddMultPA_3D(const Vector&, Vector&) const;
1770  void AddMultPA_C0_2D(const Vector&, Vector&) const;
1771  void AddMultPA_C0_3D(const Vector&, Vector&) const;
1772 
1773  void AddMultGradPA_2D(const Vector&, Vector&) const;
1774  void AddMultGradPA_3D(const Vector&, Vector&) const;
1775  void AddMultGradPA_C0_2D(const Vector&, Vector&) const;
1776  void AddMultGradPA_C0_3D(const Vector&, Vector&) const;
1777 
1778  void AssembleDiagonalPA_2D(Vector&) const;
1779  void AssembleDiagonalPA_3D(Vector&) const;
1780  void AssembleDiagonalPA_C0_2D(Vector&) const;
1781  void AssembleDiagonalPA_C0_3D(Vector&) const;
1782 
1783  void AssemblePA_Limiting();
1784  void ComputeAllElementTargets(const Vector &xe = Vector()) const;
1785 
1786  // Compute Min(Det(Jpt)) in the mesh, does not reduce over MPI.
1787  double ComputeMinDetT(const Vector &x, const FiniteElementSpace &fes);
1788  // Compute Max(mu_hat) for the TMOP_WorstCaseUntangleOptimizer_Metric,
1789  // does not reduce over MPI.
1790  double ComputeUntanglerMaxMuBarrier(const Vector &x,
1791  const FiniteElementSpace &fes);
1792 
1793 public:
1794  /** @param[in] m TMOP_QualityMetric for r-adaptivity (not owned).
1795  @param[in] tc Target-matrix construction algorithm to use (not owned).
1796  @param[in] hm TMOP_QualityMetric for h-adaptivity (not owned). */
1798  TMOP_QualityMetric *hm)
1799  : h_metric(hm), metric(m), targetC(tc), IntegRules(NULL),
1800  integ_order(-1), metric_coeff(NULL), metric_normal(1.0),
1801  lim_nodes0(NULL), lim_coeff(NULL),
1802  lim_dist(NULL), lim_func(NULL), lim_normal(1.0),
1803  adapt_lim_gf0(NULL), adapt_lim_gf(NULL), adapt_lim_coeff(NULL),
1804  adapt_lim_eval(NULL),
1805  surf_fit_gf(NULL), surf_fit_marker(NULL),
1806  surf_fit_coeff(NULL),
1807  surf_fit_eval(NULL), surf_fit_normal(1.0),
1808  discr_tc(dynamic_cast<DiscreteAdaptTC *>(tc)),
1809  fdflag(false), dxscale(1.0e3), fd_call_flag(false), exact_action(false)
1810  { PA.enabled = false; }
1811 
1813  : TMOP_Integrator(m, tc, m) { }
1814 
1815  ~TMOP_Integrator();
1816 
1817  /// Release the device memory of large PA allocations. This will copy device
1818  /// memory back to the host before releasing.
1819  void ReleasePADeviceMemory(bool copy_to_host = true);
1820 
1821  /// Prescribe a set of integration rules; relevant for mixed meshes.
1822  /** This function has priority over SetIntRule(), if both are called. */
1823  void SetIntegrationRules(IntegrationRules &irules, int order)
1824  {
1825  IntegRules = &irules;
1826  integ_order = order;
1827  }
1828 
1829  /// Sets a scaling Coefficient for the quality metric term of the integrator.
1830  /** With this addition, the integrator becomes
1831  @f$ \int w1 W(Jpt) dx @f$.
1832 
1833  Note that the Coefficient is evaluated in the physical configuration and
1834  not in the target configuration which may be undefined. */
1836 
1837  /** @brief Limiting of the mesh displacements (general version).
1838 
1839  Adds the term @f$ \int w_0 f(x, x_0, d) dx @f$, where f is a measure of
1840  the displacement between x and x_0, given the max allowed displacement d.
1841 
1842  @param[in] n0 Original mesh node coordinates (x0 above).
1843  @param[in] dist Allowed displacement in physical space (d above).
1844  @param[in] w0 Coefficient scaling the limiting integral.
1845  @param[in] lfunc TMOP_LimiterFunction defining the function f. If
1846  NULL, a TMOP_QuadraticLimiter will be used. The
1847  TMOP_Integrator assumes ownership of this pointer. */
1848  void EnableLimiting(const GridFunction &n0, const GridFunction &dist,
1849  Coefficient &w0, TMOP_LimiterFunction *lfunc = NULL);
1850 
1851  /** @brief Adds a limiting term to the integrator with limiting distance
1852  function (@a dist in the general version of the method) equal to 1. */
1853  void EnableLimiting(const GridFunction &n0, Coefficient &w0,
1854  TMOP_LimiterFunction *lfunc = NULL);
1855 
1856  /** @brief Restriction of the node positions to certain regions.
1857 
1858  Adds the term @f$ \int c (z(x) - z_0(x_0))^2 @f$, where z0(x0) is a given
1859  function on the starting mesh, and z(x) is its image on the new mesh.
1860  Minimizing this term means that a node at x0 is allowed to move to a
1861  position x(x0) only if z(x) ~ z0(x0).
1862  Such term can be used for tangential mesh relaxation.
1863 
1864  @param[in] z0 Function z0 that controls the adaptive limiting.
1865  @param[in] coeff Coefficient c for the above integral.
1866  @param[in] ae AdaptivityEvaluator to compute z(x) from z0(x0). */
1867  void EnableAdaptiveLimiting(const GridFunction &z0, Coefficient &coeff,
1868  AdaptivityEvaluator &ae);
1869 #ifdef MFEM_USE_MPI
1870  /// Parallel support for adaptive limiting.
1871  void EnableAdaptiveLimiting(const ParGridFunction &z0, Coefficient &coeff,
1872  AdaptivityEvaluator &ae);
1873 #endif
1874 
1875  /** @brief Fitting of certain DOFs to the zero level set of a function.
1876 
1877  Having a level set function s0(x0) on the starting mesh, and a set of
1878  marked nodes (or DOFs), we move these nodes to the zero level set of s0.
1879  If s(x) is the image of s0(x0) on the current mesh, this function adds to
1880  the TMOP functional the term @f$ \int c \bar{s}(x))^2 @f$, where
1881  @f$\bar{s}(x)@f$ is the restriction of s(x) on the aligned DOFs.
1882  Minimizing this term means that a marked node at x0 is allowed to move to
1883  a position x(x0) only if s(x) ~ 0.
1884  Such term can be used for surface fitting and tangential relaxation.
1885 
1886  @param[in] s0 The level set function on the initial mesh.
1887  @param[in] smarker Indicates which DOFs will be aligned.
1888  @param[in] coeff Coefficient c for the above integral.
1889  @param[in] ae AdaptivityEvaluator to compute s(x) from s0(x0). */
1890  void EnableSurfaceFitting(const GridFunction &s0,
1891  const Array<bool> &smarker, Coefficient &coeff,
1892  AdaptivityEvaluator &ae);
1893 #ifdef MFEM_USE_MPI
1894  /// Parallel support for surface fitting.
1895  void EnableSurfaceFitting(const ParGridFunction &s0,
1896  const Array<bool> &smarker, Coefficient &coeff,
1897  AdaptivityEvaluator &ae);
1898 #endif
1899  void GetSurfaceFittingErrors(double &err_avg, double &err_max);
1900  bool IsSurfaceFittingEnabled() { return (surf_fit_gf != NULL); }
1901 
1902  /// Update the original/reference nodes used for limiting.
1903  void SetLimitingNodes(const GridFunction &n0) { lim_nodes0 = &n0; }
1904 
1905  /** @brief Computes the integral of W(Jacobian(Trt)) over a target zone.
1906  @param[in] el Type of FiniteElement.
1907  @param[in] T Mesh element transformation.
1908  @param[in] elfun Physical coordinates of the zone. */
1909  virtual double GetElementEnergy(const FiniteElement &el,
1911  const Vector &elfun);
1912 
1913  /** @brief Computes the mean of the energies of the given element's children.
1914 
1915  In addition to the inputs for GetElementEnergy, this function requires an
1916  IntegrationRule to be specified that will give the decomposition of the
1917  given element based on the refinement type being considered. */
1918  virtual double GetRefinementElementEnergy(const FiniteElement &el,
1920  const Vector &elfun,
1921  const IntegrationRule &irule);
1922 
1923  /// This function is similar to GetElementEnergy, but ignores components
1924  /// such as limiting etc. to compute the element energy.
1925  virtual double GetDerefinementElementEnergy(const FiniteElement &el,
1927  const Vector &elfun);
1928 
1929  virtual void AssembleElementVector(const FiniteElement &el,
1931  const Vector &elfun, Vector &elvect);
1932 
1933  virtual void AssembleElementGrad(const FiniteElement &el,
1935  const Vector &elfun, DenseMatrix &elmat);
1936 
1938 
1940 #ifdef MFEM_USE_MPI
1942 #endif
1943 
1944  // PA extension
1946  virtual void AssemblePA(const FiniteElementSpace&);
1947 
1948  virtual void AssembleGradPA(const Vector&, const FiniteElementSpace&);
1949 
1950  virtual double GetLocalStateEnergyPA(const Vector&) const;
1951 
1952  virtual void AddMultPA(const Vector&, Vector&) const;
1953 
1954  virtual void AddMultGradPA(const Vector&, Vector&) const;
1955 
1956  virtual void AssembleGradDiagonalPA(Vector&) const;
1957 
1959 
1960  /** @brief Computes the normalization factors of the metric and limiting
1961  integrals using the mesh position given by @a x. */
1962  void EnableNormalization(const GridFunction &x);
1963 #ifdef MFEM_USE_MPI
1964  void ParEnableNormalization(const ParGridFunction &x);
1965 #endif
1966 
1967  /** @brief Enables FD-based approximation and computes dx. */
1968  void EnableFiniteDifferences(const GridFunction &x);
1969 #ifdef MFEM_USE_MPI
1971 #endif
1972 
1973  void SetFDhScale(double dxscale_) { dxscale = dxscale_; }
1974  bool GetFDFlag() const { return fdflag; }
1975  double GetFDh() const { return dx; }
1976 
1977  /** @brief Flag to control if exact action of Integration is effected. */
1978  void SetExactActionFlag(bool flag_) { exact_action = flag_; }
1979 
1980  /// Update the surface fitting weight as surf_fit_coeff *= factor;
1981  void UpdateSurfaceFittingWeight(double factor);
1982 
1983  /// Get the surface fitting weight.
1984  double GetSurfaceFittingWeight();
1985 
1986  /// Computes quantiles needed for UntangleMetrics. Note that in parallel,
1987  /// the ParFiniteElementSpace must be passed as argument for consistency
1988  /// across MPI ranks.
1989  void ComputeUntangleMetricQuantiles(const Vector &x,
1990  const FiniteElementSpace &fes);
1991 };
1992 
1994 {
1995 protected:
1996  // Integrators in the combination. Owned.
1998 
1999 public:
2001 
2003  {
2004  for (int i = 0; i < tmopi.Size(); i++) { delete tmopi[i]; }
2005  }
2006 
2007  /// Adds a new TMOP_Integrator to the combination.
2008  void AddTMOPIntegrator(TMOP_Integrator *ti) { tmopi.Append(ti); }
2009 
2011 
2012  /// Adds the limiting term to the first integrator. Disables it for the rest.
2013  void EnableLimiting(const GridFunction &n0, const GridFunction &dist,
2014  Coefficient &w0, TMOP_LimiterFunction *lfunc = NULL);
2015 
2016  /** @brief Adds the limiting term to the first integrator. Disables it for
2017  the rest (@a dist in the general version of the method) equal to 1. */
2018  void EnableLimiting(const GridFunction &n0, Coefficient &w0,
2019  TMOP_LimiterFunction *lfunc = NULL);
2020 
2021  /// Update the original/reference nodes used for limiting.
2022  void SetLimitingNodes(const GridFunction &n0);
2023 
2024  virtual double GetElementEnergy(const FiniteElement &el,
2026  const Vector &elfun);
2027  virtual void AssembleElementVector(const FiniteElement &el,
2029  const Vector &elfun, Vector &elvect);
2030  virtual void AssembleElementGrad(const FiniteElement &el,
2032  const Vector &elfun, DenseMatrix &elmat);
2033 
2034  virtual double GetRefinementElementEnergy(const FiniteElement &el,
2036  const Vector &elfun,
2037  const IntegrationRule &irule);
2038 
2039  virtual double GetDerefinementElementEnergy(const FiniteElement &el,
2041  const Vector &elfun);
2042 
2043  /// Normalization factor that considers all integrators in the combination.
2044  void EnableNormalization(const GridFunction &x);
2045 #ifdef MFEM_USE_MPI
2046  void ParEnableNormalization(const ParGridFunction &x);
2047 #endif
2048 
2049  // PA extension
2051  virtual void AssemblePA(const FiniteElementSpace&);
2052  virtual void AssembleGradPA(const Vector&, const FiniteElementSpace&);
2053  virtual double GetLocalStateEnergyPA(const Vector&) const;
2054  virtual void AddMultPA(const Vector&, Vector&) const;
2055  virtual void AddMultGradPA(const Vector&, Vector&) const;
2056  virtual void AssembleGradDiagonalPA(Vector&) const;
2057 };
2058 
2059 /// Interpolates the @a metric's values at the nodes of @a metric_gf.
2060 /** Assumes that @a metric_gf's FiniteElementSpace is initialized. */
2061 void InterpolateTMOP_QualityMetric(TMOP_QualityMetric &metric,
2062  const TargetConstructor &tc,
2063  const Mesh &mesh, GridFunction &metric_gf);
2064 }
2065 
2066 #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:269
void ComputeUntangleMetricQuantiles(const Vector &x, const FiniteElementSpace &fes)
Definition: tmop.cpp:4125
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:478
virtual void AddMultGradPA(const Vector &, Vector &) const
Method for partially assembled gradient action.
Definition: tmop_pa.cpp:271
virtual double GetDerefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
Definition: tmop.cpp:3064
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:416
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:78
virtual void SetParDiscreteTargetSize(const ParGridFunction &tspec_)
Definition: tmop.cpp:1689
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:116
Shifted barrier form of metric 56 (area, ideal barrier metric), 2D.
Definition: tmop.hpp:565
TMOP_Metric_347(double gamma_)
Definition: tmop.hpp:904
3D shifted barrier form of metric 316 (not typed).
Definition: tmop.hpp:920
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:1277
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:386
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:671
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
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:768
void AssembleElemVecSurfFit(const FiniteElement &el_x, IsoparametricTransformation &Tpr, DenseMatrix &mat)
Definition: tmop.cpp:3491
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:239
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
VectorCoefficient * vector_tspec
Definition: tmop.hpp:1311
3D non-barrier Shape (S) metric.
Definition: tmop.hpp:939
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition: tmop.cpp:245
void EnableFiniteDifferences(const GridFunction &x)
Enables FD-based approximation and computes dx.
Definition: tmop.cpp:3993
const TargetType target_type
Definition: tmop.hpp:1193
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:1027
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:1024
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
DenseMatrix PMatO
Definition: tmop.hpp:1638
virtual void AddMultPA(const Vector &, Vector &) const
Method for partially assembled action.
Definition: tmop.cpp:4321
DenseTensor Jtr
Definition: tmop.hpp:1668
double & min_detT
Definition: tmop.hpp:366
3D barrier Shape+Size (VS) metric, well-posed (invex).
Definition: tmop.hpp:744
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:690
const DenseMatrix * Jtr
Definition: tmop.hpp:26
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:747
2D non-barrier Shape+Size (VS) metric.
Definition: tmop.hpp:452
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:875
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:776
double GetSurfaceFittingWeight()
Get the surface fitting weight.
Definition: tmop.cpp:3818
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:537
TMOP_Metric_066(double gamma_)
Definition: tmop.hpp:460
3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
Definition: tmop.hpp:872
Base class for vector Coefficients that optionally depend on time and space.
double GetLocalStateEnergyPA_3D(const Vector &) const
Definition: tmop_pa_w3.cpp:155
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:266
virtual ~TMOP_Metric_334()
Definition: tmop.hpp:892
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:991
DenseMatrix Jrt
Definition: tmop.hpp:1638
const GridFunction * lim_dist
Definition: tmop.hpp:1590
virtual void SetSerialDiscreteTargetOrientation(const GridFunction &tspec_)
Definition: tmop.cpp:1802
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
Compute the local energy.
Definition: tmop.cpp:4194
double epsilon
Definition: ex25.cpp:140
struct mfem::TMOP_Integrator::@23 PA
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field, int new_nodes_ordering=Ordering::byNODES)=0
virtual ~TMOP_Metric_333()
Definition: tmop.hpp:868
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:1097
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:489
double GetGamma() const
Definition: tmop.hpp:845
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:328
TMOP_AMetric_126(double gamma_)
Definition: tmop.hpp:1041
virtual double GetLocalStateEnergyPA(const Vector &) const
Compute the local (to the MPI rank) energy with partial assembly.
Definition: tmop_pa.cpp:296
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:401
const double eps
Definition: tmop.hpp:670
void EnableNormalization(const GridFunction &x)
Computes the normalization factors of the metric and limiting integrals using the mesh position given...
Definition: tmop.cpp:3829
TMOP_QualityMetric & GetAMRQualityMetric()
Definition: tmop.hpp:1937
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
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:4347
3D barrier Shape (S) metric, well-posed (polyconvex &amp; invex).
Definition: tmop.hpp:625
TMOP_QualityMetric * metric
Definition: tmop.hpp:1572
void SetAdaptivityEvaluator(AdaptivityEvaluator *ae)
Definition: tmop.hpp:1495
virtual ~TMOP_Metric_332()
Definition: tmop.hpp:847
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:418
virtual ~TMOP_LimiterFunction()
Virtual destructor.
Definition: tmop.hpp:1074
Base class for limiting functions to be used in class TMOP_Integrator.
Definition: tmop.hpp:1057
virtual double EvalW(const DenseMatrix &Jpt) const =0
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:1251
const GridFunction * GetNodes() const
Get the nodes to be used in the target-matrix construction.
Definition: tmop.hpp:1249
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:468
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:2164
3D non-barrier Skew metric.
Definition: tmop.hpp:230
2D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:493
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:1030
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:941
void ParUpdateAfterMeshTopologyChange()
Definition: tmop.cpp:1658
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:250
DenseMatrix Jpr
Definition: tmop.hpp:1638
IntegrationRules * IntegRules
Definition: tmop.hpp:1576
double GetGamma() const
Definition: tmop.hpp:914
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:498
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:1038
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:667
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:1309
double GetGamma() const
Definition: tmop.hpp:469
void SetFDhScale(double dxscale_)
Definition: tmop.hpp:1973
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:200
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:1532
virtual void ComputeElementTargetsGradient(const IntegrationRule &ir, const Vector &elfun, IsoparametricTransformation &Tpr, DenseTensor &dJtr) const
Definition: tmop.cpp:1586
Container class for integration rules.
Definition: intrules.hpp:311
ParFiniteElementSpace * pfes
Definition: tmop.hpp:1122
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:766
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
3D barrier Shape (S) metric, well-posed (polyconvex &amp; invex).
Definition: tmop.hpp:604
void GetSurfaceFittingErrors(double &err_avg, double &err_max)
Definition: tmop.cpp:2786
int Size() const
Returns the size of the vector.
Definition: vector.hpp:200
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:310
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:1066
void SetVolumeScale(double vol_scale)
Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
Definition: tmop.hpp:1252
2D barrier shape (S) metric (polyconvex).
Definition: tmop.hpp:275
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...
Definition: tmop_pa.cpp:137
void EnableAdaptiveLimiting(const GridFunction &z0, Coefficient &coeff, AdaptivityEvaluator &ae)
Restriction of the node positions to certain regions.
Definition: tmop.cpp:2715
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:901
3D non-barrier metric without a type.
Definition: tmop.hpp:707
2D barrier size (V) metric (polyconvex).
Definition: tmop.hpp:475
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:885
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition: tmop.hpp:1823
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:1018
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:211
void EnableNormalization(const GridFunction &x)
Normalization factor that considers all integrators in the combination.
Definition: tmop.cpp:4264
2D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:1033
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:679
bool IsSurfaceFittingEnabled()
Definition: tmop.hpp:1900
double & min_detT
Definition: tmop.hpp:689
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:251
3D barrier Shape+Size (VS) metric, well-posed (invex).
Definition: tmop.hpp:786
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:603
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:236
virtual void AssemblePA(const FiniteElementSpace &)
Method defining partial assembly.
Definition: tmop.cpp:4296
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:292
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:335
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:354
void SetSerialMetaInfo(const Mesh &m, const FiniteElementSpace &f)
Definition: tmop.cpp:2623
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
Definition: mesh.hpp:1837
void UpdateGradientTargetSpecification(const Vector &x, double dx, bool use_flag=false, int x_ordering=Ordering::byNODES)
Definition: tmop.cpp:2511
const IntegrationRule * ir
Definition: tmop.hpp:1676
void SetExactActionFlag(bool flag_)
Flag to control if exact action of Integration is effected.
Definition: tmop.hpp:1978
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:436
void AssembleElementVectorFD(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
Definition: tmop.cpp:3643
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:621
const GridFunction * nodes
Definition: tmop.hpp:1190
TMOP_WorstCaseUntangleOptimizer_Metric(TMOP_QualityMetric &tmop_metric_, int exponent_=1, double alpha_=1.5, double detT_ep_=0.0001, double muT_ep_=0.0001, BarrierType btype_=BarrierType::None, WorstCaseType wctype_=WorstCaseType::None)
Definition: tmop.hpp:155
virtual void SetSerialDiscreteTargetSkew(const GridFunction &tspec_)
Definition: tmop.cpp:1782
virtual double EvalWBarrier(const DenseMatrix &Jpt) const
Definition: tmop.cpp:88
void UpdateSurfaceFittingWeight(double factor)
Update the surface fitting weight as surf_fit_coeff *= factor;.
Definition: tmop.cpp:3806
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:721
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:525
void SetRefinementSubElement(int amr_el_)
Definition: tmop.hpp:1553
2D non-barrier Aspect ratio metric.
Definition: tmop.hpp:245
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:854
TMOP_Metric_080(double gamma_)
Definition: tmop.hpp:501
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:1019
2D barrier Shape+Size (VS) metric (not polyconvex).
Definition: tmop.hpp:314
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:856
Default limiter function in TMOP_Integrator.
Definition: tmop.hpp:1078
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:1949
void AssembleElementGradExact(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
Definition: tmop.cpp:3286
Array< Vector * > ElemDer
Definition: tmop.hpp:1625
3D barrier Shape+Size (VS) metric, well-posed (invex).
Definition: tmop.hpp:765
virtual ~TMOP_AMetric_126()
Definition: tmop.hpp:1050
const Vector & GetTspecPertMixH()
Definition: tmop.hpp:1503
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:4165
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:1153
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:810
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:504
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:457
double GetFDh() const
Definition: tmop.hpp:1975
2D Shifted barrier form of shape metric (mu_2).
Definition: tmop.hpp:363
TMOP_Metric_313(double &mindet)
Definition: tmop.hpp:693
TMOP_Metric_352(double &t0)
Definition: tmop.hpp:927
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:1552
Coefficient * scalar_tspec
Definition: tmop.hpp:1310
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:713
virtual void SetSerialDiscreteTargetAspectRatio(const GridFunction &tspec_)
Definition: tmop.cpp:1792
TMOPMatrixCoefficient(int dim)
Definition: tmop.hpp:1295
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:868
2D barrier shape (S) metric (not polyconvex).
Definition: tmop.hpp:433
const DofToQuad * maps
Definition: tmop.hpp:1672
virtual ~TMOP_QualityMetric()
Definition: tmop.hpp:35
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:150
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:822
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:520
double ComputeUntanglerMaxMuBarrier(const Vector &x, const FiniteElementSpace &fes)
Definition: tmop.cpp:4063
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
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:43
const GeometricFactors * geom
Definition: tmop.hpp:1674
3D barrier Shape (S) metric, well-posed (polyconvex &amp; invex).
Definition: tmop.hpp:646
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:177
virtual void AddMultGradPA(const Vector &, Vector &) const
Method for partially assembled gradient action.
Definition: tmop.cpp:4329
virtual void ComputeElementTargetsGradient(const IntegrationRule &ir, const Vector &elfun, IsoparametricTransformation &Tpr, DenseTensor &dJtr) const
Definition: tmop.cpp:1529
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
2D non-barrier shape (S) metric.
Definition: tmop.hpp:296
2D barrier Size (V) metric (polyconvex).
Definition: tmop.hpp:979
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:1812
Coefficient * adapt_lim_coeff
Definition: tmop.hpp:1602
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:848
virtual ~TMOPMatrixCoefficient()
Definition: tmop.hpp:1303
virtual ~AdaptivityEvaluator()
Definition: tmop.cpp:2654
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:710
virtual void AddQualityMetric(TMOP_QualityMetric *tq, double wt=1.0)
Definition: tmop.hpp:89
void SetLimitingNodes(const GridFunction &n0)
Update the original/reference nodes used for limiting.
Definition: tmop.cpp:4186
virtual ~TMOP_Metric_080()
Definition: tmop.hpp:512
void SetDiscreteTargetBase(const GridFunction &tspec_)
Definition: tmop.cpp:1735
virtual ~TMOP_Metric_066()
Definition: tmop.hpp:471
AdaptivityEvaluator * adapt_lim_eval
Definition: tmop.hpp:1603
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_pa.cpp:23
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition: tmop.cpp:757
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:976
DiscreteAdaptTC * discr_tc
Definition: tmop.hpp:1612
void EnableLimiting(const GridFunction &n0, const GridFunction &dist, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
Limiting of the mesh displacements (general version).
Definition: tmop.cpp:2690
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:637
void AddMultGradPA_C0_2D(const Vector &, Vector &) const
virtual WorstCaseType GetWorstCaseType()
Definition: tmop.hpp:193
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:964
void SetCoefficient(Coefficient &w1)
Sets a scaling Coefficient for the quality metric term of the integrator.
Definition: tmop.hpp:1835
2D non-barrier size (V) metric (not polyconvex).
Definition: tmop.hpp:398
void UpdateTargetSpecificationAtNode(const FiniteElement &el, ElementTransformation &T, int nodenum, int idir, const Vector &IntData)
Definition: tmop.cpp:1884
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:877
void AssembleDiagonalPA_C0_2D(Vector &) const
void SetTspecFromIntRule(int e_id, const IntegrationRule &intrule)
Definition: tmop.cpp:1915
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:462
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:1218
TargetType GetTargetType() const
Definition: tmop.hpp:1254
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:56
void UpdateAfterMeshPositionChange(const Vector &new_x, int new_x_ordering=Ordering::byNODES)
Definition: tmop.cpp:3957
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:491
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:336
3D barrier metric without a type.
Definition: tmop.hpp:725
DiscreteAdaptTC(TargetType ttype)
Definition: tmop.hpp:1405
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition: tmop.cpp:697
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:901
void ResetRefinementTspecData()
Definition: tmop.hpp:1538
void AssembleDiagonalPA_3D(Vector &) const
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:550
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:1229
void SetNodes(const GridFunction &n)
Set the nodes to be used in the target-matrix construction.
Definition: tmop.hpp:1246
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:807
void AssembleElementVectorExact(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
Definition: tmop.cpp:3147
void AssembleDiagonalPA_2D(Vector &) const
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:812
virtual void AssembleGradDiagonalPA(Vector &) const
Method for computing the diagonal of the gradient with partial assembly.
Definition: tmop.cpp:4313
const TargetConstructor * targetC
Definition: tmop.hpp:1573
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:841
Array< Vector * > ElemPertEnergy
Definition: tmop.hpp:1626
const GridFunction * lim_nodes0
Definition: tmop.hpp:1587
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:646
const ParGridFunction * adapt_lim_pgf0
Definition: tmop.hpp:1599
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:856
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition: tmop.cpp:800
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:832
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:457
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
Definition: nonlininteg.cpp:25
3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
Definition: tmop.hpp:896
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
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:561
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:844
Abstract class used to define combination of metrics with constant coefficients.
Definition: tmop.hpp:82
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:367
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:4304
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:190
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:789
virtual void SetParDiscreteTargetAspectRatio(const ParGridFunction &tspec_)
Definition: tmop.cpp:1709
virtual void AssemblePA(const FiniteElementSpace &)
Method defining partial assembly.
Definition: tmop_pa.cpp:178
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:1137
virtual ~TargetConstructor()
Definition: tmop.hpp:1232
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:574
void ParEnableNormalization(const ParGridFunction &x)
Definition: tmop.cpp:4280
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:569
double GetLocalStateEnergyPA_2D(const Vector &) const
Definition: tmop_pa_w2.cpp:145
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
Definition: tmop.cpp:3118
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:299
virtual double GetLocalStateEnergyPA(const Vector &) const
Compute the local (to the MPI rank) energy with partial assembly.
Definition: tmop.cpp:4337
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:531
GridFunction * GetTSpecData()
Get the entire tspec.
Definition: tmop.hpp:1455
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:925
FiniteElementSpace * tspec_fesv
Definition: tmop.hpp: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:311
Coefficient * metric_coeff
Definition: tmop.hpp:1580
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:110
void AssembleElemGradSurfFit(const FiniteElement &el_x, IsoparametricTransformation &Tpr, DenseMatrix &mat)
Definition: tmop.cpp:3545
virtual bool ContainsVolumeInfo() const
Checks if the target matrices contain non-trivial size specification.
Definition: tmop.cpp:1444
TMOP_QualityMetric * h_metric
Definition: tmop.hpp:1571
Array< TMOP_QualityMetric * > tmop_q_arr
Definition: tmop.hpp:85
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition: tmop.hpp:47
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:782
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:988
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition: tmop.cpp:998
virtual BarrierType GetBarrierType()
Definition: tmop.hpp:191
void RestoreTargetSpecificationAtNode(ElementTransformation &T, int nodenum)
Definition: tmop.cpp:1901
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:913
void AssemblePA_Limiting()
Definition: tmop_pa.cpp:51
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:303
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition: tmop.cpp:534
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:498
2D barrier Shape+Orientation (OS) metric (polyconvex).
Definition: tmop.hpp:1015
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:522
void SetTspecAtIndex(int idx, const GridFunction &tspec_)
Definition: tmop.cpp:1759
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:2978
void ComputeAvgVolume() const
Definition: tmop.cpp:1330
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:1008
Coefficient * lim_coeff
Definition: tmop.hpp:1588
3D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:807
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:33
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition: tmop.cpp:1130
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:496
TMOP_Metric_328(double gamma_)
Definition: tmop.hpp:815
virtual double Eval(const Vector &x, const Vector &x0, double dist) const
Returns the limiting function, f(x, x0, d).
Definition: tmop.hpp:1081
AdaptivityEvaluator * surf_fit_eval
Definition: tmop.hpp:1609
virtual ~TMOP_QuadraticLimiter()
Definition: tmop.hpp:1105
FiniteElementSpace * GetTSpecFESpace()
Get the FESpace associated with tspec.
Definition: tmop.hpp:1453
2D barrier (not a shape) metric (polyconvex).
Definition: tmop.hpp:382
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:1038
virtual void SetTargetJacobian(const DenseMatrix &Jtr_)
Specify the reference-element -&gt; target-element Jacobian matrix for the point of interest.
Definition: tmop.hpp:43
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:899
TMOP_Metric_211(double epsilon=1e-4)
Definition: tmop.hpp:553
void AssembleElemGradAdaptLim(const FiniteElement &el, IsoparametricTransformation &Tpr, const IntegrationRule &ir, const Vector &weights, DenseMatrix &m)
Definition: tmop.cpp:3430
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
Definition: tmop.cpp:3132
const DofToQuad * maps_lim
Definition: tmop.hpp:1673
2D untangling metric.
Definition: tmop.hpp:546
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:221
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
3D Size (V) untangling metric.
Definition: tmop.hpp:667
const double eps
Definition: tmop.hpp:549
virtual void SetTargetJacobian(const DenseMatrix &Jtr_)
Specify the reference-element -&gt; target-element Jacobian matrix for the point of interest.
Definition: tmop.hpp:95
2D barrier size (V) metric (polyconvex).
Definition: tmop.hpp:415
void SetLimitingNodes(const GridFunction &n0)
Update the original/reference nodes used for limiting.
Definition: tmop.hpp:1903
DenseMatrix tspec_refine
Definition: tmop.hpp:1364
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:581
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:924
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition: vector.cpp:460
TMOPMatrixCoefficient * matrix_tspec
Definition: tmop.hpp:1312
TMOP_Metric_333(double gamma_)
Definition: tmop.hpp:859
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:1458
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:455
const IntegrationRule & ActionIntegrationRule(const FiniteElement &el) const
Definition: tmop.hpp:1746
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:856
Base class for Matrix Coefficients that optionally depend on time and space.
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:385
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:126
void ReleasePADeviceMemory(bool copy_to_host=true)
Definition: tmop.cpp:2664
bool UsesPhysicalCoordinates() const
Return true if the methods ComputeElementTargets(), ComputeAllElementTargets(), and ComputeElementTar...
Definition: tmop.hpp:1259
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:1077
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:262
void Clear()
Delete the matrix data array (if owned) and reset the matrix state.
Definition: densemat.hpp:98
void DisableLimiting()
Definition: tmop.hpp:1731
double GetLocalStateEnergyPA_C0_2D(const Vector &) const
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:956
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
Definition: tmop.cpp:4222
double DistanceSquaredTo(const double *p) const
Compute the square of the Euclidean distance to another vector.
Definition: vector.hpp:649
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:443
DiscreteAdaptTC * GetDiscreteAdaptTC() const
Definition: tmop.hpp:1958
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
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:1222
TMOP_Metric_252(double &t0)
Note that t0 is stored by reference.
Definition: tmop.hpp:573
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:475
DenseTensor Jtrcomp
Definition: tmop.hpp:1371
void SetParMetaInfo(const ParMesh &m, const ParFiniteElementSpace &f)
Parallel version of SetSerialMetaInfo.
Definition: tmop.cpp:2634
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:659
const IntegrationRule * IntRule
Definition: nonlininteg.hpp:30
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:509
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:1144
void ResetUpdateFlags()
Used in combination with the Update methods to avoid extra computations.
Definition: tmop.hpp:1447
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:761
AnalyticAdaptTC(TargetType ttype)
Definition: tmop.hpp:1315
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:1000
void AddTMOPIntegrator(TMOP_Integrator *ti)
Adds a new TMOP_Integrator to the combination.
Definition: tmop.hpp:2008
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:545
GridFunction * tspec_gf
Definition: tmop.hpp:1377
2D non-barrier Skew metric.
Definition: tmop.hpp:215
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:900
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:370
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:296
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:353
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:949
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:553
void UpdateAfterMeshTopologyChange()
Definition: tmop.cpp:2815
void AddMultPA_3D(const Vector &, Vector &) const
Definition: tmop_pa_p3.cpp:192
virtual void SetMaxMuT(double max_muT_)
Definition: tmop.hpp:189
DenseMatrix PMatI
Definition: tmop.hpp:1638
virtual void SetMinDetT(double min_detT_)
Definition: tmop.hpp:187
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:728
2D barrier Shape+Size (VS) metric (not polyconvex).
Definition: tmop.hpp:332
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:344
void AddMultPA_2D(const Vector &, Vector &) const
Definition: tmop_pa_p2.cpp:167
void SetTransformation(ElementTransformation &)
The method HyperelasticModel::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:3977
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition: tmop.cpp:1215
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:3852
TMOP_Metric_334(double gamma_)
Definition: tmop.hpp:880
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:1608
Class for integration point with weight.
Definition: intrules.hpp:25
void UpdateHessianTargetSpecification(const Vector &x, double dx, bool use_flag=false, int x_ordering=Ordering::byNODES)
Definition: tmop.cpp:2546
virtual ~TMOP_Metric_347()
Definition: tmop.hpp:916
AdaptivityEvaluator * adapt_eval
Definition: tmop.hpp:1394
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition: tmop.cpp:962
ParFiniteElementSpace * ptspec_fesv
Definition: tmop.hpp:1380
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:588
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:540
void AssembleDiagonalPA_C0_3D(Vector &) const
const FiniteElementSpace * fes
Definition: tmop.hpp:1675
Array< TMOP_Integrator * > tmopi
Definition: tmop.hpp:1997
A standard isoparametric element transformation.
Definition: eltrans.hpp:361
void ComputeMinJac(const Vector &x, const FiniteElementSpace &fes)
Definition: tmop.cpp:3923
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:784
TMOP_Integrator(TMOP_QualityMetric *m, TargetConstructor *tc, TMOP_QualityMetric *hm)
Definition: tmop.hpp:1797
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition: tmop.cpp:834
void Diag(double c, int n)
Creates n x n diagonal matrix with diagonal elements c.
Definition: densemat.cpp:1351
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:1009
const IntegrationRule & EnergyIntegrationRule(const FiniteElement &el) const
Definition: tmop.hpp:1737
const IntegrationRule & GradientIntegrationRule(const FiniteElement &el) const
Definition: tmop.hpp:1751
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:278
double GetFDDerivative(const FiniteElement &el, ElementTransformation &T, Vector &elfun, const int nodenum, const int idir, const double baseenergy, bool update_stored)
Definition: tmop.cpp:3621
double GetGamma() const
Definition: tmop.hpp:890
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.hpp:1006
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:997
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:973
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:512
double surf_fit_normal
Definition: tmop.hpp:1610
3D Shape (S) metric, untangling version of 303.
Definition: tmop.hpp:686
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:2841
void AssembleElemVecAdaptLim(const FiniteElement &el, IsoparametricTransformation &Tpr, const IntegrationRule &ir, const Vector &weights, DenseMatrix &mat)
Definition: tmop.cpp:3394
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:104
void ParUpdateAfterMeshTopologyChange()
Definition: tmop.cpp:2828
virtual double GetDerefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
Definition: tmop.cpp:4251
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:889
void AssembleGradPA_3D(const Vector &) const
void Destroy()
Destroy a vector.
Definition: vector.hpp:590
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:1372
3D non-barrier Aspect ratio metric.
Definition: tmop.hpp:260
const Vector & GetTspecPert2H()
Definition: tmop.hpp:1502
const Vector & GetTspecPert1H()
Definition: tmop.hpp:1501
void AssembleElementGradFD(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat)
Definition: tmop.cpp:3709
ParFiniteElementSpace * GetTSpecParFESpace()
Definition: tmop.hpp:1460
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:877
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:968
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:663
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:278
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:722
virtual void SetParDiscreteTargetSpec(const ParGridFunction &tspec_)
Definition: tmop.cpp:1729
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:317
TMOP_QualityMetric * sh_metric
Definition: tmop.hpp:812
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:803
Abstract class for hyperelastic models.
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:210
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...
virtual void AssembleGradDiagonalPA(Vector &) const
Method for computing the diagonal of the gradient with partial assembly.
Definition: tmop_pa.cpp:224
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:406
GridFunction * adapt_lim_gf
Definition: tmop.hpp:1601
TargetConstructor(TargetType ttype, MPI_Comm mpicomm)
Constructor for use in parallel.
Definition: tmop.hpp:1228
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:984
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:272
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:357
GridFunction * surf_fit_gf
Definition: tmop.hpp:1606
void ComputeAllElementTargets(const Vector &xe=Vector()) const
Definition: tmop_pa.cpp:166
TMOP_Integrator(TMOP_QualityMetric *m, TargetConstructor *tc)
Definition: tmop.hpp:1812
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:706
void FinalizeParDiscreteTargetSpec(const ParGridFunction &tspec_)
Definition: tmop.cpp:1633
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:703
ParGridFunction * tspec_pgf
Definition: tmop.hpp:1382
bool GetFDFlag() const
Definition: tmop.hpp:1974
Vector data type.
Definition: vector.hpp:60
void GetDiscreteTargetSpec(GridFunction &tspec_, int idx)
Get one of the discrete fields from tspec.
Definition: tmop.cpp:1833
DenseMatrix DSh
Definition: tmop.hpp:1638
const Array< bool > * surf_fit_marker
Definition: tmop.hpp:1607
DenseMatrix DS
Definition: tmop.hpp:1638
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition: tmop.cpp:1056
TMOP_Metric_022(double &t0)
Definition: tmop.hpp:370
3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
Definition: tmop.hpp:851
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:628
void SetTspecDataForDerefinement(FiniteElementSpace *fes)
Computes target specification data with respect to the coarse FE space.
Definition: tmop.cpp:1941
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:71
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:621
TargetType
Target-matrix construction algorithms supported by this class.
Definition: tmop.hpp:1170
void UpdateTargetSpecification(const Vector &new_x, bool use_flag=false, int new_x_ordering=Ordering::byNODES)
Definition: tmop.cpp:1864
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:483
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:2750
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:942
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:651
virtual double EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition: tmop.cpp:23
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:877
TMOP_Metric_311(double epsilon=1e-4)
Definition: tmop.hpp:674
virtual void SetSerialDiscreteTargetSpec(const GridFunction &tspec_)
Definition: tmop.cpp:1858
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:982
virtual void SetParDiscreteTargetSkew(const ParGridFunction &tspec_)
Definition: tmop.cpp:1699
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:454
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:1088
virtual double GetRefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, const IntegrationRule &irule)
Definition: tmop.cpp:4238
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:254
InvariantsEvaluator2D< double > ie
Definition: tmop.hpp:1036
virtual int Id() const
Return the metric ID.
Definition: tmop.hpp:642
virtual void SetParDiscreteTargetOrientation(const ParGridFunction &tspec_)
Definition: tmop.cpp:1719
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:1166
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:1293
void AddMultPA_C0_3D(const Vector &, Vector &) const
virtual ~DiscreteAdaptTC()
Definition: tmop.cpp:2613
TMOP_Metric_332(double gamma_)
Definition: tmop.hpp:835
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:256
virtual void SetAnalyticTargetSpec(Coefficient *sspec, VectorCoefficient *vspec, TMOPMatrixCoefficient *mspec)
Definition: tmop.cpp:1543
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:284
3D barrier Shape+Size (VS) metric (polyconvex).
Definition: tmop.hpp:828
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:953
Class for parallel meshes.
Definition: pmesh.hpp:32
virtual void SetSerialDiscreteTargetSize(const GridFunction &tspec_)
Definition: tmop.cpp:1772
TMOP_QualityMetric * sz_metric
Definition: tmop.hpp:832
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:1592
const Array< TMOP_Integrator * > & GetTMOPIntegrators() const
Definition: tmop.hpp:2010
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:933
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:814
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:607
const GridFunction * adapt_lim_gf0
Definition: tmop.hpp:1597
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
Definition: tmop.cpp:4206
double ComputeMinDetT(const Vector &x, const FiniteElementSpace &fes)
Definition: tmop.cpp:4021
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:1093
DenseMatrix Jpt
Definition: tmop.hpp:1638
double GetGamma() const
Definition: tmop.hpp:510
2D barrier Shape+Orientation (OS) metric (polyconvex).
Definition: tmop.hpp:516
Array< double > wt_arr
Definition: tmop.hpp:86
void ParEnableNormalization(const ParGridFunction &x)
Definition: tmop.cpp:3839
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:1177
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:224
virtual void AddMultPA(const Vector &, Vector &) const
Method for partially assembled action.
Definition: tmop_pa.cpp:249
DenseMatrix P
Definition: tmop.hpp:1638
MPI_Comm GetComm() const
Definition: tmop.hpp:1236
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:1189
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:1238
void UpdateAfterMeshTopologyChange()
Update all discrete fields based on tspec and update for AMR.
Definition: tmop.cpp:1847
void ResetDerefinementTspecData()
Definition: tmop.hpp:1545
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition: tmop.cpp:920
bool Parallel() const
Definition: tmop.hpp:1235
3D barrier Shape (S) metric, well-posed (polyconvex &amp; invex).
Definition: tmop.hpp:585
void AssembleGradPA_C0_3D(const Vector &) const
2D non-barrier metric without a type.
Definition: tmop.hpp:197
virtual double EvalW(const DenseMatrix &Jpt) const
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants...
Definition: tmop.cpp:1169
FiniteElementSpace * fes
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:590
bool ComputeAllElementTargets(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
Definition: tmop_pa.cpp:123
InvariantsEvaluator3D< double > ie
Definition: tmop.hpp:649
virtual ~TMOP_Metric_328()
Definition: tmop.hpp:824
FiniteElementSpace * coarse_tspec_fesv
Definition: tmop.hpp:1376
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition: tmop.hpp:1565
2D non-barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition: tmop.hpp:348