MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
tmop.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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
16#include "nonlininteg.hpp"
17#include "../linalg/dual.hpp"
18
19namespace mfem
20{
21
22using AD1Type = future::dual<real_t, real_t>;
23using AD2Type = future::dual<AD1Type, AD1Type>;
24
25/** @brief Abstract class for local mesh quality metrics in the target-matrix
26 optimization paradigm (TMOP) by P. Knupp et al. */
28{
29protected:
30 const DenseMatrix *Jtr; /**< Jacobian of the reference-element to
31 target-element transformation. */
32
33 /** @brief The method HyperelasticModel::SetTransformation() is hidden
34 for TMOP_QualityMetric%s, because it is not used. */
36
37 /** @brief See AssembleH(). This is a default implementation for the case
38 when the 2nd derivatives of the metric are pre-computed and stored into
39 @a H. This function is used in combination with AD-based computations. */
40 void DefaultAssembleH(const DenseTensor &H, const DenseMatrix &DS,
41 const real_t weight, DenseMatrix &A) const;
42
43public:
44 TMOP_QualityMetric() : Jtr(NULL) { }
45 virtual ~TMOP_QualityMetric() { }
46
47 /** @brief Specify the reference-element -> target-element Jacobian matrix
48 for the point of interest.
49
50 The specified Jacobian matrix, #Jtr, can be used by metrics that cannot
51 be written just as a function of the target->physical Jacobian matrix,
52 Jpt. */
53 virtual void SetTargetJacobian(const DenseMatrix &Jtr_) { Jtr = &Jtr_; }
54
55 /** @brief Evaluates the metric in matrix form (opposed to invariant form).
56 Used for validating the invariant evaluations. */
57 virtual real_t EvalWMatrixForm(const DenseMatrix &Jpt) const
58 { return -1.0; /* not implemented -> checks would fail. */ }
59
60 /** @brief Evaluate the strain energy density function, W = W(Jpt), by using
61 the 2D or 3D matrix invariants, see linalg/invariants.hpp.
62 @param[in] Jpt Represents the target->physical transformation
63 Jacobian matrix. */
64 virtual real_t EvalW(const DenseMatrix &Jpt) const = 0;
65
66 /** @brief Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
67 @param[in] Jpt Represents the target->physical transformation
68 Jacobian matrix.
69 @param[out] P The evaluated 1st Piola-Kirchhoff stress tensor. */
70 virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const = 0;
71
72 /** Compute dmu/dW */
73 virtual void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const
74 { PW = 0.0;}
75
76 /// @brief First-derivative hook for AD-based computations.
77 /// @warning Not for public use. Internal use for AD-based computations.
78 virtual AD1Type EvalW_AD1(const std::vector<AD1Type> &T,
79 const std::vector<AD1Type> &W) const
80 {
81 MFEM_ABORT("EvalW_AD1 not implemented for this metric");
82 }
83
84 /// @brief Second-derivative hook for AD-based computations.
85 /// @warning Not for public use. Internal use for AD-based computations.
86 virtual AD2Type EvalW_AD2(const std::vector<AD2Type> &T,
87 const std::vector<AD2Type> &W) const
88 {
89 MFEM_ABORT("EvalW_AD2 not implemented for this metric");
90 }
91
92 /** @brief Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor
93 and assemble its contribution to the local gradient matrix 'A'.
94 @param[in] Jpt Represents the target->physical transformation
95 Jacobian matrix.
96 @param[in] DS Gradient of the basis matrix (dof x dim).
97 @param[in] weight Quadrature weight coefficient for the point.
98 @param[in,out] A Local gradient matrix where the contribution from this
99 point will be added.
100
101 Computes weight * d(dW_dxi)_d(xj) at the current point, for all i and j,
102 where x1 ... xn are the FE dofs. This function is usually defined using
103 the matrix invariants and their derivatives. */
104 virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
105 const real_t weight, DenseMatrix &A) const = 0;
106
107 /** @brief Return the metric ID. */
108 virtual int Id() const { return 0; }
109};
110
111class TargetConstructor;
112
113/// Abstract class used to define explicit combination of metrics with constant
114/// coefficients.
116{
117protected:
118 Array<TMOP_QualityMetric *> tmop_q_arr; //the metrics are not owned
120
121public:
122 virtual void AddQualityMetric(TMOP_QualityMetric *tq, real_t wt = 1.0)
123 {
124 tmop_q_arr.Append(tq);
125 wt_arr.Append(wt);
126 }
127
128 void SetTargetJacobian(const DenseMatrix &Jtr_) override
129 {
130 for (int i = 0; i < tmop_q_arr.Size(); i++)
131 {
132 tmop_q_arr[i]->SetTargetJacobian(Jtr_);
133 }
134 }
135
136 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
137
138 real_t EvalW(const DenseMatrix &Jpt) const override;
139
140 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
141
142 void EvalPW(const DenseMatrix &Jpt, DenseMatrix &P) const override;
143
144 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
145 const real_t weight, DenseMatrix &A) const override;
146
147 AD1Type EvalW_AD1(const std::vector<AD1Type> &T,
148 const std::vector<AD1Type> &W) const override;
149
150 AD2Type EvalW_AD2(const std::vector<AD2Type> &T,
151 const std::vector<AD2Type> &W) const override;
152
153 /// Computes the averages of all metrics (integral of metric / volume).
154 /// Works in parallel when called with a ParGridFunction.
156 const TargetConstructor &tc, Vector &averages,
157 bool use_pa = false,
158 const IntegrationRule *IntRule = nullptr) const;
159
160 /// Computes weights so that the averages of all metrics are equal, and the
161 /// weights sum to one. Works in parallel when called with a ParGridFunction.
163 const TargetConstructor &tc, Vector &weights,
164 bool use_pa = false,
165 const IntegrationRule *IntRule = nullptr) const;
166
168 const TargetConstructor &tc,
169 int m_index, real_t &energy, real_t &vol,
170 const IntegrationRule &ir) const;
172 const TargetConstructor &tc,
173 int m_index, real_t &energy, real_t &vol,
174 const IntegrationRule &ir) const;
175
176 void GetWeights(Array<real_t> &weights) const { weights = wt_arr; }
177
178 /// Changes the weights of the metrics in the combination.
179 void SetWeights(const Vector &weights)
180 {
181 MFEM_VERIFY(tmop_q_arr.Size() == weights.Size(), "Incorrect #weights");
182 for (int i = 0; i < tmop_q_arr.Size(); i++) { wt_arr[i] = weights(i); }
183 }
184};
185
186/// Simultaneous Untangler + Worst Case Improvement Metric
187/// Uses a base metric mu and is defined as:
188/// mu_tilde = mu_hat, when WorstCaseType = None,
189/// = mu_hat/(beta - mu_hat), when WorstCaseType = Beta,
190/// = mu_hat^p, when WorstCaseType = PMean,
191/// where beta = max(mu_hat) + muT_ep,
192/// and mu_hat = (mu/2phi(tau,ep)) where
193/// 2phi(tau,ep) = 1, when when BarrierType = None,
194/// = 2*(tau - min(alpha*min(tau)-detT_ep,0)), when BarrierType = Shifted
195/// = tau^2 + sqrt(tau^2 + ep^2), when BarrierType = Pseudo
196/// where tau = det(T), and max(mu_hat) and min(tau) are computed over the
197/// entire mesh.
198/// Ultimately, this metric can be used for mesh untangling with the BarrierType
199/// option and for worst case quality improvement with the WorstCaseType option.
201{
202public:
203 enum class BarrierType
204 {
205 None,
206 Shifted,
207 Pseudo
208 };
209 enum class WorstCaseType
210 {
211 None,
212 Beta,
213 PMean
214 };
215
216protected:
217 TMOP_QualityMetric &tmop_metric; // non-barrier metric to use
218 real_t min_detT; // minimum Jacobian in the mesh
219 real_t max_muT; // max mu_k/phi(tau,ep) in the mesh
220 int exponent; // used for p-mean metrics
221 real_t alpha; // scaling factor for min(det(T))
222 real_t detT_ep; // small constant subtracted from min(detT)
223 real_t muT_ep; // small constant added to muT term
226
227public:
229 int exponent_ = 1,
230 real_t alpha_ = 1.5,
231 real_t detT_ep_ = 0.0001,
232 real_t muT_ep_ = 0.0001,
235 tmop_metric(tmop_metric_), exponent(exponent_), alpha(alpha_),
236 detT_ep(detT_ep_), muT_ep(muT_ep_), btype(btype_), wctype(wctype_)
237 {
238 MFEM_VERIFY(wctype == WorstCaseType::None,
239 "Worst-case optimization has not been fully developed!");
241 {
242 const int m_id = tmop_metric.Id();
243 MFEM_VERIFY(m_id == 4 || m_id == 14 || m_id == 66,
244 "Incorrect input barrier metric -- must be 4 / 14 / 66");
245 }
246 }
247
248 real_t EvalW(const DenseMatrix &Jpt) const override;
249
250 AD1Type EvalW_AD1(const std::vector<AD1Type> &T,
251 const std::vector<AD1Type> &W) const override;
252
253 AD2Type EvalW_AD2(const std::vector<AD2Type> &T,
254 const std::vector<AD2Type> &W) const override;
255
256 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
257
258 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
259 const real_t weight, DenseMatrix &A) const override;
260
261 // Compute mu_hat.
262 real_t EvalWBarrier(const DenseMatrix &Jpt) const;
263
264 void SetMinDetT(real_t min_detT_) { min_detT = min_detT_; }
265
266 void SetMaxMuT(real_t max_muT_) { max_muT = max_muT_; }
267
269
271};
272
273/// 0 metric
275{
276public:
277 // W = 0.
278 real_t EvalW(const DenseMatrix &Jpt) const override {return 0.0;}
279
280 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override {P=0.0;}
281
282 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
283 const real_t weight, DenseMatrix &A) const override {A=0.0;}
284
285 int Id() const override { return 0; }
286};
287
288/// 2D non-barrier metric without a type.
290{
291protected:
293
294public:
295 // W = |J|^2.
296 real_t EvalW(const DenseMatrix &Jpt) const override;
297
298 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
299
300 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
301 const real_t weight, DenseMatrix &A) const override;
302
303 int Id() const override { return 1; }
304};
305
306// TODO: Remove in MFEM 5.0
307/// (DEPRECATED) 2D non-barrier Skew metric.
308/// We recommend TMOP_AMetric_050 metric for controlling skewness only. In
309/// general, the Size+Skew metric TMOP_AMetric_051 may produce better
310/// meshes than a purely Skew metric.
311class MFEM_DEPRECATED TMOP_Metric_skew2D : public TMOP_QualityMetric
312{
313public:
314 // W = 0.5 (1 - cos(angle_Jpr - angle_Jtr)).
315 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
316 real_t EvalW(const DenseMatrix &Jpt) const override
317 { return EvalWMatrixForm(Jpt); }
318
319 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
320
321 void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const override;
322
323 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
324 const real_t weight, DenseMatrix &A) const override;
325};
326
327/// 3D non-barrier Skew metric.
329{
330public:
331 // W = 1/6 (3 - sum_i cos(angle_Jpr_i - angle_Jtr_i)), i = 1..3.
332 real_t EvalW(const DenseMatrix &Jpt) const override;
333
334 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
335 { MFEM_ABORT("Not implemented"); }
336
337 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
338 const real_t weight, DenseMatrix &A) const override
339 { MFEM_ABORT("Not implemented"); }
340};
341
342/// 2D non-barrier Aspect ratio metric.
344{
345public:
346 // W = 0.5 (ar_Jpr/ar_Jtr + ar_Jtr/ar_Jpr) - 1.
347 real_t EvalW(const DenseMatrix &Jpt) const override;
348
349 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
350 { MFEM_ABORT("Not implemented"); }
351
352 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
353 const real_t weight, DenseMatrix &A) const override
354 { MFEM_ABORT("Not implemented"); }
355};
356
357/// 3D non-barrier Aspect ratio metric.
359{
360public:
361 // W = 1/3 sum [0.5 (ar_Jpr_i/ar_Jtr_i + ar_Jtr_i/ar_Jpr_i) - 1], i = 1..3.
362 real_t EvalW(const DenseMatrix &Jpt) const override;
363
364 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
365 { MFEM_ABORT("Not implemented"); }
366
367 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
368 const real_t weight, DenseMatrix &A) const override
369 { MFEM_ABORT("Not implemented"); }
370};
371
372/// 2D barrier shape (S) metric (polyconvex).
373/// Grade - A.
375{
376protected:
378
379public:
380 // W = 0.5 |J|^2 / det(J) - 1.
381 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
382
383 // W = 0.5 I1b - 1.
384 real_t EvalW(const DenseMatrix &Jpt) const override;
385
386 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
387
388 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
389 const real_t weight, DenseMatrix &A) const override;
390
391 int Id() const override { return 2; }
392};
393
394/// 2D non-barrier shape (S) metric.
395/// Grade - F.
397{
398protected:
400
401 template<typename type>
402 type EvalW_AD_impl(const std::vector<type> &T,
403 const std::vector<type> &W) const;
404
405public:
406 // W = |J|^2 - 2*det(J)
407 real_t EvalW(const DenseMatrix &Jpt) const override;
408
409 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
410
411 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
412 const real_t weight, DenseMatrix &A) const override;
413
414 AD1Type EvalW_AD1(const std::vector<AD1Type> &T,
415 const std::vector<AD1Type> &W) const override;
416
417 AD2Type EvalW_AD2(const std::vector<AD2Type> &T,
418 const std::vector<AD2Type> &W) const override;
419
420 int Id() const override { return 4; }
421};
422
423/// 2D barrier Shape+Size (VS) metric (not polyconvex).
425{
426protected:
428
429public:
430 // W = |J - J^-t|^2.
431 real_t EvalW(const DenseMatrix &Jpt) const override;
432
433 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
434
435 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
436 const real_t weight, DenseMatrix &A) const override;
437
438 int Id() const override { return 7; }
439};
440
441/// 2D barrier Shape+Size (VS) metric (not polyconvex).
443{
444protected:
446
447public:
448 // W = det(J) * |J - J^-t|^2.
449 real_t EvalW(const DenseMatrix &Jpt) const override;
450
451 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
452
453 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
454 const real_t weight, DenseMatrix &A) const override;
455
456 int Id() const override { return 9; }
457};
458
459/// 2D non-barrier Shape+Size+Orientation (VOS) metric (polyconvex).
461{
462protected:
464
465 template <typename type>
466 type EvalW_AD_impl(const std::vector<type> &T,
467 const std::vector<type> &W) const;
468
469public:
470 // W = |J - I|^2.
471 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
472
473 // W = I1[J-I].
474 real_t EvalW(const DenseMatrix &Jpt) const override;
475
476 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
477
478 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
479 const real_t weight, DenseMatrix &A) const override;
480
481 AD1Type EvalW_AD1(const std::vector<AD1Type> &T,
482 const std::vector<AD1Type> &W) const override;
483
484 AD2Type EvalW_AD2(const std::vector<AD2Type> &T,
485 const std::vector<AD2Type> &W) const override;
486
487 int Id() const override { return 14; }
488};
489
490/// 2D Shifted barrier form of shape metric (mu_2).
492{
493protected:
496
497public:
499
500 // W = 0.5(|J|^2 - 2det(J)) / (det(J) - tau0).
501 real_t EvalW(const DenseMatrix &Jpt) const override;
502
503 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
504
505 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
506 const real_t weight, DenseMatrix &A) const override;
507
508 int Id() const override { return 22; }
509};
510
511/// 2D barrier shape metric (polyconvex).
512/// Grade - A.
514{
515protected:
517
518public:
519 // W = 0.5 |J^t J|^2 / det(J)^2 - 1.
520 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
521
522 // W = 0.5 I1b^2 - 2.
523 real_t EvalW(const DenseMatrix &Jpt) const override;
524
525 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
526
527 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
528 const real_t weight, DenseMatrix &A) const override;
529
530 int Id() const override { return 50; }
531};
532
533/// 2D non-barrier size (V) metric (not polyconvex).
534/// Grade - F.
536{
537protected:
539
540 template<typename type>
541 type EvalW_AD_impl(const std::vector<type> &T,
542 const std::vector<type> &W) const;
543
544public:
545 // W = (det(J) - 1)^2.
546 real_t EvalW(const DenseMatrix &Jpt) const override;
547
548 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
549
550 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
551 const real_t weight, DenseMatrix &A) const override;
552
553 AD1Type EvalW_AD1(const std::vector<AD1Type> &T,
554 const std::vector<AD1Type> &W) const override;
555
556 AD2Type EvalW_AD2(const std::vector<AD2Type> &T,
557 const std::vector<AD2Type> &W) const override;
558
559 int Id() const override { return 55; }
560
561};
562
563/// 2D barrier size (V) metric (polyconvex).
564/// Grade - C.
566{
567protected:
569
570public:
571 // W = 0.5 (det(J) + 1 / det(J)) - 1.
572 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
573
574 // W = 0.5 (I2b + 1/I2b) - 1.
575 real_t EvalW(const DenseMatrix &Jpt) const override;
576
577 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
578
579 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
580 const real_t weight, DenseMatrix &A) const override;
581
582 int Id() const override { return 56; }
583};
584
585/// 2D barrier shape (S) metric (not polyconvex).
587{
588protected:
590
591public:
592 // W = |J^t J|^2 / det(J)^2 - 2|J|^2 / det(J) + 2
593 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
594
595 // W = I1b (I1b - 2).
596 real_t EvalW(const DenseMatrix &Jpt) const override;
597
598 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
599
600 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
601 const real_t weight, DenseMatrix &A) const override;
602
603 int Id() const override { return 58; }
604};
605
606/// 2D non-barrier Shape+Size (VS) metric.
607/// Grade - F.
609{
610protected:
613
614public:
617 {
618 // (1-gamma) mu_4 + gamma mu_55
619 AddQualityMetric(sh_metric, 1.-gamma);
621 }
622 int Id() const override { return 66; }
623 real_t GetGamma() const { return wt_arr[1]; }
624
625 virtual ~TMOP_Metric_066() { delete sh_metric; delete sz_metric; }
626};
627
628/// 2D barrier size (V) metric (polyconvex).
629/// Grade - C.
631{
632protected:
634
635public:
636 // W = 0.5 (det(J) - 1 / det(J))^2.
637 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
638
639 // W = 0.5 (I2 + 1 / I2) - 1.0.
640 real_t EvalW(const DenseMatrix &Jpt) const override;
641
642 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
643
644 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
645 const real_t weight, DenseMatrix &A) const override;
646
647 int Id() const override { return 77; }
648};
649
650/// 2D barrier Shape+Size (VS) metric (polyconvex).
651/// Grade - A.
653{
654protected:
657
658public:
661 {
662 // (1-gamma) mu_2 + gamma mu_77
663 AddQualityMetric(sh_metric, 1.0 - gamma);
665 }
666
667 int Id() const override { return 80; }
668 real_t GetGamma() const { return wt_arr[1]; }
669
670 virtual ~TMOP_Metric_080() { delete sh_metric; delete sz_metric; }
671};
672
673/// 2D barrier Shape+Orientation (OS) metric (polyconvex).
675{
676public:
677 // W = |T-T'|^2, where T'= |T|*I/sqrt(2).
678 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
679 real_t EvalW(const DenseMatrix &Jpt) const override
680 { return EvalWMatrixForm(Jpt); }
681
682 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
683
684 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
685 const real_t weight, DenseMatrix &A) const override;
686
687 int Id() const override { return 85; }
688};
689
690/// 2D compound barrier Shape+Size (VS) metric (balanced).
692{
693protected:
696
697public:
700 {
701 // mu_50 + lambda mu_77.
702 // 1 <= lambda <= 4 should produce best asymptotic balance.
705 }
706
707 int Id() const override { return 90; }
708 virtual ~TMOP_Metric_090() { delete sh_metric; delete sz_metric; }
709};
710
711/// 2D compound barrier Shape+Size (VS) metric (balanced).
713{
714protected:
717
718public:
721 {
722 // mu_2 + lambda mu_56.
723 // 1 <= lambda <= 2 should produce best asymptotic balance.
726 }
727
728 int Id() const override { return 94; }
729 virtual ~TMOP_Metric_094() { delete sh_metric; delete sz_metric; }
730};
731
732/// 2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
734{
735public:
736 // W = 1/tau |T-I|^2.
737 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
738 real_t EvalW(const DenseMatrix &Jpt) const override
739 { return EvalWMatrixForm(Jpt); }
740
741 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
742
743 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
744 const real_t weight, DenseMatrix &A) const override;
745
746 int Id() const override { return 98; }
747};
748
749/// 2D untangling metric.
751{
752protected:
753 const real_t eps;
755
756public:
758
759 // W = (det(J) - 1)^2 - det(J) + sqrt(det(J)^2 + eps).
760 real_t EvalW(const DenseMatrix &Jpt) const override;
761
762 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
763
764 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
765 const real_t weight, DenseMatrix &A) const override;
766
767 int Id() const override { return 211; }
768};
769
770/// Shifted barrier form of metric 56 (area, ideal barrier metric), 2D
772{
773protected:
776
777public:
778 /// Note that @a t0 is stored by reference
780
781 // W = 0.5(det(J) - 1)^2 / (det(J) - tau0).
782 real_t EvalW(const DenseMatrix &Jpt) const override;
783
784 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
785
786 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
787 const real_t weight, DenseMatrix &A) const override;
788
789 int Id() const override { return 252; }
790};
791
792/// 3D barrier Shape (S) metric, well-posed (polyconvex & invex).
794{
795protected:
797
798public:
799 // W = 1/3 |J| |J^-1| - 1.
800 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
801
802 // W = 1/3 sqrt(I1b * I2b) - 1
803 real_t EvalW(const DenseMatrix &Jpt) const override;
804
805 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
806
807 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
808 const real_t weight, DenseMatrix &A) const override;
809
810 int Id() const override { return 301; }
811};
812
813/// 3D barrier Shape (S) metric, well-posed (polyconvex & invex).
815{
816protected:
818
819public:
820 // W = |J|^2 |J^{-1}|^2 / 9 - 1.
821 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
822
823 // W = I1b * I2b / 9 - 1.
824 real_t EvalW(const DenseMatrix &Jpt) const override;
825
826 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
827
828 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
829 const real_t weight, DenseMatrix &A) const override;
830
831 int Id() const override { return 302; }
832};
833
834/// 3D barrier Shape (S) metric, well-posed (polyconvex & invex).
836{
837protected:
839
840public:
841 // W = |J|^2 / 3 / det(J)^(2/3) - 1.
842 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
843
844 // W = I1b / 3 - 1.
845 real_t EvalW(const DenseMatrix &Jpt) const override;
846
847 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
848
849 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
850 const real_t weight, DenseMatrix &A) const override;
851
852 int Id() const override { return 303; }
853};
854
855/// 3D barrier Shape (S) metric, well-posed (polyconvex & invex).
857{
858protected:
860
861public:
862 // W = |J|^3 / 3^(3/2) / det(J) - 1.
863 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
864
865 // W = (I1b/3)^3/2 - 1.
866 real_t EvalW(const DenseMatrix &Jpt) const override;
867
868 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
869
870 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
871 const real_t weight, DenseMatrix &A) const override;
872
873 int Id() const override { return 304; }
874};
875
876/// 3D Size (V) untangling metric.
878{
879protected:
880 const real_t eps;
882
883public:
885
886 // W = (det(J) - 1)^2 - det(J) + (det(J)^2 + eps)^(1/2).
887 real_t EvalW(const DenseMatrix &Jpt) const override;
888
889 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
890
891 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
892 const real_t weight, DenseMatrix &A) const override;
893
894 int Id() const override { return 311; }
895};
896
897/// 3D Shape (S) metric, untangling version of 303.
899{
900protected:
903
904public:
905 TMOP_Metric_313(real_t &mindet) : min_detT(mindet) { }
906
907 // W = 1/3 |J|^2 / [det(J)-tau0]^(-2/3).
908 real_t EvalW(const DenseMatrix &Jpt) const override;
909
910 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
911
912 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
913 const real_t weight, DenseMatrix &A) const override;
914
915 int Id() const override { return 313; }
916};
917
918/// 3D Size (V) metric.
920{
921protected:
923
924public:
925 // W = (det(J) - 1)^2.
926 real_t EvalW(const DenseMatrix &Jpt) const override;
927
928 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
929
930 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
931 const real_t weight, DenseMatrix &A) const override;
932
933 int Id() const override { return 315; }
934};
935
936/// 3D Size (V) metric.
938{
939protected:
941
942public:
943 // W = 0.5 (det(J) + 1/det(J)) - 1.
944 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
945
946 // W = 0.5 (I3b + 1/I3b) - 1.
947 real_t EvalW(const DenseMatrix &Jpt) const override;
948
949 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
950
951 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
952 const real_t weight, DenseMatrix &A) const override;
953
954 int Id() const override { return 316; }
955};
956
957/// 3D Size (V) metric.
959{
960protected:
962
963public:
964 // W = 0.5 (det(J)^2 + 1/det(J)^2) - 1.
965 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
966
967 // W = 0.5 (I3 + 1/I3) - 1.
968 real_t EvalW(const DenseMatrix &Jpt) const override;
969
970 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
971
972 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
973 const real_t weight, DenseMatrix &A) const override;
974
975 int Id() const override { return 318; }
976};
977
978/// 3D barrier Shape+Size (VS) metric, well-posed (invex).
980{
981protected:
983
984public:
985 // W = |J - J^-t|^2.
986 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
987
988 // W = I1 + I2/I3 - 6.
989 real_t EvalW(const DenseMatrix &Jpt) const override;
990
991 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
992
993 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
994 const real_t weight, DenseMatrix &A) const override;
995
996 int Id() const override { return 321; }
997};
998
999/// 3D barrier Shape+Size (VS) metric, well-posed (invex).
1001{
1002protected:
1004
1005public:
1006 // W = |J - adjJ^-t|^2.
1007 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
1008
1009 // W = I1b / (I3b^-1/3) / 6 + I2b (I3b^1/3) / 6 - 1
1010 real_t EvalW(const DenseMatrix &Jpt) const override;
1011
1012 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
1013
1014 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
1015 const real_t weight, DenseMatrix &A) const override;
1016
1017 int Id() const override { return 322; }
1018};
1019
1020/// 3D barrier Shape+Size (VS) metric, well-posed (invex).
1022{
1023protected:
1025
1026public:
1027 // W = |J|^3 - 3 sqrt(3) ln(det(J)) - 3 sqrt(3).
1028 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
1029
1030 // W = I1^3/2 - 3 sqrt(3) ln(I3b) - 3 sqrt(3).
1031 real_t EvalW(const DenseMatrix &Jpt) const override;
1032
1033 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
1034
1035 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
1036 const real_t weight, DenseMatrix &A) const override;
1037
1038 int Id() const override { return 323; }
1039};
1040
1041/// 3D compound barrier Shape+Size (VS) metric (polyconvex, balanced).
1043{
1044protected:
1047
1048public:
1051 {
1052 // lambda mu_301 + mu_316.
1053 // 3/8 <= lambda <= 9/8 should produce best asymptotic balance.
1056 }
1057
1058 int Id() const override { return 328; }
1059 virtual ~TMOP_Metric_328() { delete sh_metric; delete sz_metric; }
1060};
1061
1062/// 3D compound barrier Shape+Size (VS) metric (polyconvex).
1064{
1065protected:
1067
1068public:
1071 {
1072 // (1-gamma) mu_302 + gamma mu_315
1073 AddQualityMetric(sh_metric, 1.-gamma);
1075 }
1076
1077 int Id() const override { return 332; }
1078 real_t GetGamma() const { return wt_arr[1]; }
1079
1080 virtual ~TMOP_Metric_332() { delete sh_metric; delete sz_metric; }
1081};
1082
1083/// 3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
1085{
1086protected:
1089
1090public:
1093 {
1094 // (1-gamma) mu_302 + gamma mu_316
1095 AddQualityMetric(sh_metric, 1.-gamma);
1097 }
1098
1099 int Id() const override { return 333; }
1100 virtual ~TMOP_Metric_333() { delete sh_metric; delete sz_metric; }
1101};
1102
1103/// 3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
1105{
1106protected:
1109
1110public:
1113 {
1114 // (1-gamma) mu_303 + gamma mu_316
1115 AddQualityMetric(sh_metric, 1.-gamma);
1117 }
1118
1119 int Id() const override { return 334; }
1120 real_t GetGamma() const { return wt_arr[1]; }
1121
1122 virtual ~TMOP_Metric_334() { delete sh_metric; delete sz_metric; }
1123};
1124
1125/// 3D compound barrier Shape+Size (VS) metric (polyconvex, balanced).
1127{
1128protected:
1131
1132public:
1135 {
1136 // mu_302 + lambda mu_318.
1137 // 4/9 <= lambda <= 3 should produce best asymptotic balance.
1139 AddQualityMetric(sz_metric, 0.5 * (4.0/9.0 + 3.0));
1140 }
1141
1142 int Id() const override { return 338; }
1143 virtual ~TMOP_Metric_338() { delete sh_metric; delete sz_metric; }
1144};
1145
1146/// 3D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
1148{
1149protected:
1151
1152public:
1153 // W = 1/(tau^0.5) |T-I|^2.
1154 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
1155 real_t EvalW(const DenseMatrix &Jpt) const override
1156 { return EvalWMatrixForm(Jpt); }
1157
1158 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
1159
1160 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
1161 const real_t weight, DenseMatrix &A) const override;
1162
1163 int Id() const override { return 342; }
1164};
1165
1166/// 3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
1168{
1169protected:
1172
1173public:
1176 {
1177 // (1-gamma) mu_304 + gamma mu_316
1178 AddQualityMetric(sh_metric, 1.-gamma);
1180 }
1181
1182 int Id() const override { return 347; }
1183 real_t GetGamma() const { return wt_arr[1]; }
1184
1185 virtual ~TMOP_Metric_347() { delete sh_metric; delete sz_metric; }
1186};
1187
1188/// 3D shifted barrier form of metric 316 (not typed).
1190{
1191protected:
1194
1195public:
1197
1198 // W = 0.5(det(J) - 1)^2 / (det(J) - tau0).
1199 real_t EvalW(const DenseMatrix &Jpt) const override;
1200
1201 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
1202
1203 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
1204 const real_t weight, DenseMatrix &A) const override;
1205
1206 int Id() const override { return 352; }
1207};
1208
1209/// 3D non-barrier Shape (S) metric.
1211{
1212protected:
1214
1215public:
1216 // W = |J|^3 / 3^(3/2) - det(J).
1217 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
1218
1219 // W = (I1b/3)^3/2 - 1.
1220 real_t EvalW(const DenseMatrix &Jpt) const override;
1221
1222 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
1223
1224 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
1225 const real_t weight, DenseMatrix &A) const override;
1226
1227 int Id() const override { return 360; }
1228};
1229
1230/// A-metrics
1231/// 2D barrier Shape (S) metric (polyconvex).
1233{
1234protected:
1236
1237public:
1238 // (1/4 alpha) | A - (adj A)^t W^t W / omega |^2
1239 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
1240 real_t EvalW(const DenseMatrix &Jpt) const override
1241 { return EvalWMatrixForm(Jpt); }
1242
1243 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
1244
1245 void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const override;
1246
1247 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
1248 const real_t weight, DenseMatrix &A) const override;
1249};
1250
1251/// 2D barrier Size (V) metric (polyconvex).
1253{
1254protected:
1256
1257public:
1258 // 0.5 * ( sqrt(alpha/omega) - sqrt(omega/alpha) )^2
1259 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
1260 real_t EvalW(const DenseMatrix &Jpt) const override
1261 { return EvalWMatrixForm(Jpt); }
1262
1263 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
1264
1265 void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const override;
1266
1267 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
1268 const real_t weight, DenseMatrix &A) const override;
1269};
1270
1271/// 2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
1273{
1274protected:
1276
1277public:
1278 // (1/alpha) | A - W |^2
1279 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
1280 real_t EvalW(const DenseMatrix &Jpt) const override
1281 { return EvalWMatrixForm(Jpt); }
1282
1283 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
1284
1285 void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const override;
1286
1287 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
1288 const real_t weight, DenseMatrix &A) const override;
1289};
1290
1291/// 2D barrier Skew (Q) metric.
1293{
1294protected:
1296
1297public:
1298 // [ 1.0 - cos( phi_A - phi_W ) ] / (sin phi_A * sin phi_W)
1299 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
1300 real_t EvalW(const DenseMatrix &Jpt) const override
1301 { return EvalWMatrixForm(Jpt); }
1302
1303 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
1304
1305 void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const override;
1306
1307 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
1308 const real_t weight, DenseMatrix &A) const override;
1309};
1310
1311/// 2D barrier Size+Skew (VQ) metric.
1313{
1314protected:
1316
1317public:
1318 // [ 0.5 * (ups_A / ups_W + ups_W / ups_A) - cos(phi_A - phi_W) ] /
1319 // (sin phi_A * sin phi_W)
1320 // where ups = l_1 l_2 sin(phi)
1321 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
1322 real_t EvalW(const DenseMatrix &Jpt) const override
1323 { return EvalWMatrixForm(Jpt); }
1324
1325 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
1326
1327 void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const override;
1328
1329 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
1330 const real_t weight, DenseMatrix &A) const override;
1331};
1332
1333/// 2D barrier Shape+Orientation (OS) metric (polyconvex).
1335{
1336protected:
1338
1339public:
1340 // (1/2 alpha) | A - (|A|/|W|) W |^2
1341 real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override;
1342 real_t EvalW(const DenseMatrix &Jpt) const override
1343 { return EvalWMatrixForm(Jpt); }
1344
1345 void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override;
1346
1347 void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const override;
1348
1349 void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
1350 const real_t weight, DenseMatrix &A) const override;
1351};
1352
1353/// 2D barrier Shape+Size (VS) metric (polyconvex).
1355{
1356protected:
1359
1360public:
1363 {
1364 // (1-gamma) nu_11 + gamma nu_14
1365 AddQualityMetric(sh_metric, 1.-gamma);
1367 }
1368
1369 virtual ~TMOP_AMetric_126() { delete sh_metric; delete sz_metric; }
1370};
1371
1372/// 2D barrier Shape+Skew (SQ) metric.
1373/// gamma is recommended to be in (0, 0.9) as a pure skew metric has poor
1374/// convergence properties.
1376{
1377protected:
1380
1381public:
1384 {
1385 // (1-gamma) mu_2 + gamma nu_50
1386 AddQualityMetric(sh_metric, 1.0 - gamma);
1388 }
1389
1390 int Id() const override { return 49; }
1391 real_t GetGamma() const { return wt_arr[1]; }
1392
1393 virtual ~TMOP_AMetric_049() { delete sh_metric; delete sk_metric; }
1394};
1395
1396/// Base class for limiting functions to be used in class TMOP_Integrator.
1397/** This class represents a scalar function f(x, x0, d), where x and x0 are
1398 positions in physical space, and d is a reference physical distance
1399 associated with the point x0. */
1401{
1402public:
1403 /// Returns the limiting function, f(x, x0, d).
1404 virtual real_t Eval(const Vector &x, const Vector &x0, real_t d) const = 0;
1405
1406 /** @brief Returns the gradient of the limiting function f(x, x0, d) with
1407 respect to x. */
1408 virtual void Eval_d1(const Vector &x, const Vector &x0, real_t dist,
1409 Vector &d1) const = 0;
1410
1411 /** @brief Returns the Hessian of the limiting function f(x, x0, d) with
1412 respect to x. */
1413 virtual void Eval_d2(const Vector &x, const Vector &x0, real_t dist,
1414 DenseMatrix &d2) const = 0;
1415
1416 /// Virtual destructor.
1418};
1419
1420/// Default limiter function in TMOP_Integrator.
1422{
1423public:
1424 real_t Eval(const Vector &x, const Vector &x0, real_t dist) const override
1425 {
1426 MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
1427
1428 return 0.5 * x.DistanceSquaredTo(x0) / (dist * dist);
1429 }
1430
1431 void Eval_d1(const Vector &x, const Vector &x0, real_t dist,
1432 Vector &d1) const override
1433 {
1434 MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
1435
1436 d1.SetSize(x.Size());
1437 subtract(1.0 / (dist * dist), x, x0, d1);
1438 }
1439
1440 void Eval_d2(const Vector &x, const Vector &x0, real_t dist,
1441 DenseMatrix &d2) const override
1442 {
1443 MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
1444
1445 d2.Diag(1.0 / (dist * dist), x.Size());
1446 }
1447
1449};
1450
1451/// Exponential limiter function in TMOP_Integrator.
1453{
1454public:
1455 real_t Eval(const Vector &x, const Vector &x0, real_t dist) const override
1456 {
1457 MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
1458
1459 return exp(10.0*((x.DistanceSquaredTo(x0) / (dist * dist))-1.0));
1460 }
1461
1462 void Eval_d1(const Vector &x, const Vector &x0, real_t dist,
1463 Vector &d1) const override
1464 {
1465 MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
1466
1467 d1.SetSize(x.Size());
1468 real_t dist_squared = dist*dist;
1469 subtract(20.0*exp(10.0*((x.DistanceSquaredTo(x0) / dist_squared) - 1.0)) /
1470 dist_squared, x, x0, d1);
1471 }
1472
1473 void Eval_d2(const Vector &x, const Vector &x0, real_t dist,
1474 DenseMatrix &d2) const override
1475 {
1476 MFEM_ASSERT(x.Size() == x0.Size(), "Bad input.");
1477 Vector tmp;
1478 tmp.SetSize(x.Size());
1479 real_t dist_squared = dist*dist;
1480 real_t dist_squared_squared = dist_squared*dist_squared;
1481 real_t f = exp(10.0*((x.DistanceSquaredTo(x0) / dist_squared)-1.0));
1482
1483 subtract(x,x0,tmp);
1484 d2.SetSize(x.Size());
1485 d2(0,0) = ((400.0*tmp(0)*tmp(0)*f)/dist_squared_squared)+(20.0*f/dist_squared);
1486 d2(1,1) = ((400.0*tmp(1)*tmp(1)*f)/dist_squared_squared)+(20.0*f/dist_squared);
1487 d2(0,1) = (400.0*tmp(0)*tmp(1)*f)/dist_squared_squared;
1488 d2(1,0) = d2(0,1);
1489
1490 if (x.Size() == 3)
1491 {
1492 d2(0,2) = (400.0*tmp(0)*tmp(2)*f)/dist_squared_squared;
1493 d2(1,2) = (400.0*tmp(1)*tmp(2)*f)/dist_squared_squared;
1494 d2(2,0) = d2(0,2);
1495 d2(2,1) = d2(1,2);
1496 d2(2,2) = ((400.0*tmp(2)*tmp(2)*f)/dist_squared_squared)+(20.0*f/dist_squared);
1497 }
1498
1499 }
1500
1502};
1503
1504class FiniteElementCollection;
1505class FiniteElementSpace;
1506class ParFiniteElementSpace;
1507
1509{
1510protected:
1511 // Owned.
1514
1515#ifdef MFEM_USE_MPI
1516 // Owned.
1519#endif
1520
1521public:
1523 {
1524#ifdef MFEM_USE_MPI
1525 pmesh = NULL;
1526 pfes = NULL;
1527#endif
1528 }
1529 virtual ~AdaptivityEvaluator();
1530
1531 /// Specifies the Mesh and FiniteElementSpace of the solution that will
1532 /// be evaluated. The given mesh will be copied into the internal object.
1533 void SetSerialMetaInfo(const Mesh &m,
1534 const FiniteElementSpace &f);
1535
1536#ifdef MFEM_USE_MPI
1537 /// Parallel version of SetSerialMetaInfo.
1538 void SetParMetaInfo(const ParMesh &m,
1539 const ParFiniteElementSpace &f);
1540#endif
1541
1542 // TODO use GridFunctions to make clear it's on the ldofs? Then do we
1543 // need the SetMetaInfo at all -- the space and mesh can be extracted?
1544 virtual void SetInitialField(const Vector &init_nodes,
1545 const Vector &init_field) = 0;
1546
1547 /// Called when the FE space of the final field is different than
1548 /// the FE space of the initial field.
1550
1551 /** @brief Perform field transfer between the original and a new mesh. The
1552 source mesh and field are given by SetInitialField().
1553
1554 @param[in] new_mesh_nodes Mesh node positions of the new mesh (ldofs).
1555 It is assumed that this is the field's mesh.
1556 @param[out] new_field Result of the transfer (ldofs).
1557 @param[in] nodes_ordering Ordering of new_mesh_nodes. */
1558 virtual void ComputeAtNewPosition(const Vector &new_mesh_nodes,
1559 Vector &new_field,
1560 int nodes_ordering = Ordering::byNODES) = 0;
1561
1562 /** @brief Using the source mesh and field given by SetInitialField(),
1563 compute corresponding values at specified physical positions.
1564
1565 @param[in] positions Physical positions to compute values.
1566 @param[out] values Computed field values.
1567 @param[in] p_ordering Ordering of the positions Vector. */
1568 virtual void ComputeAtGivenPositions(const Vector &positions,
1569 Vector &values,
1570 int p_ordering = Ordering::byNODES) = 0;
1571
1572 void ClearGeometricFactors();
1573};
1574
1575/** @brief Base class representing target-matrix construction algorithms for
1576 mesh optimization via the target-matrix optimization paradigm (TMOP). */
1577/** This class is used by class TMOP_Integrator to construct the target Jacobian
1578 matrices (reference-element to target-element) at quadrature points. It
1579 supports a set of algorithms chosen by the #TargetType enumeration.
1580
1581 New target-matrix construction algorithms can be defined by deriving new
1582 classes and overriding the methods ComputeElementTargets() and
1583 ContainsVolumeInfo(). */
1585{
1586public:
1587 /// Target-matrix construction algorithms supported by this class.
1589 {
1591 Ideal shape, unit size; the nodes are not used. */
1593 Ideal shape, equal size/volume; the given nodes define the total target
1594 volume; for each mesh element, the target volume is the average volume
1595 multiplied by the volume scale, set with SetVolumeScale(). */
1597 Ideal shape, given size/volume; the given nodes define the target
1598 volume at all quadrature points. */
1600 Given shape, given size/volume; the given nodes define the exact target
1601 Jacobian matrix at all quadrature points. */
1602 GIVEN_FULL /**<
1603 Full target tensor is specified at every quadrature point. */
1605
1606protected:
1607 // Nodes that are used in ComputeElementTargets(), depending on target_type.
1608 const GridFunction *nodes; // not owned
1612 bool uses_phys_coords; // see UsesPhysicalCoordinates()
1613
1614 /// Cached copy of GeomToPerfGeomJac used on device.
1616 /// Geometry type of current W matrix (used for cache invalidation).
1618
1619#ifdef MFEM_USE_MPI
1620 MPI_Comm comm;
1621#endif
1622
1623 // should be called only if avg_volume == 0.0, i.e. avg_volume is not
1624 // computed yet
1625 void ComputeAvgVolume() const;
1626
1627 template<int DIM>
1629 const IntegrationRule &ir,
1630 const Vector &xe,
1631 DenseTensor &Jtr) const;
1632
1633 // CPU fallback that uses ComputeElementTargets()
1635 const IntegrationRule &ir,
1636 const Vector &xe,
1637 DenseTensor &Jtr) const;
1638
1639public:
1640 /// Constructor for use in serial
1642 : nodes(NULL), avg_volume(), volume_scale(1.0), target_type(ttype),
1643 uses_phys_coords(false)
1644 {
1645#ifdef MFEM_USE_MPI
1646 comm = MPI_COMM_NULL;
1647#endif
1648 }
1649#ifdef MFEM_USE_MPI
1650 /// Constructor for use in parallel
1651 TargetConstructor(TargetType ttype, MPI_Comm mpicomm)
1652 : nodes(NULL), avg_volume(), volume_scale(1.0), target_type(ttype),
1653 uses_phys_coords(false), comm(mpicomm) { }
1654#endif
1656
1657#ifdef MFEM_USE_MPI
1658 bool Parallel() const { return (comm != MPI_COMM_NULL); }
1659 MPI_Comm GetComm() const { return comm; }
1660#else
1661 bool Parallel() const { return false; }
1662#endif
1663
1664 /** @brief Set the nodes to be used in the target-matrix construction.
1665
1666 This method should be called every time the target nodes are updated
1667 externally and recomputation of the target average volume is needed. The
1668 nodes are used by all target types except IDEAL_SHAPE_UNIT_SIZE. */
1669 void SetNodes(const GridFunction &n) { nodes = &n; avg_volume = 0.0; }
1670
1671 /** @brief Get the nodes to be used in the target-matrix construction. */
1672 const GridFunction *GetNodes() const { return nodes; }
1673
1674 /// Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
1675 void SetVolumeScale(real_t vol_scale) { volume_scale = vol_scale; }
1676
1678
1679 /** @brief Return true if the methods ComputeElementTargets(),
1680 ComputeAllElementTargets(), and ComputeElementTargetsGradient() use the
1681 physical node coordinates provided by the parameters 'elfun', or 'xe'. */
1683
1684 /// Checks if the target matrices contain non-trivial size specification.
1685 virtual bool ContainsVolumeInfo() const;
1686
1687 /** @brief Given an element and quadrature rule, computes ref->target
1688 transformation Jacobians for each quadrature point in the element.
1689 The physical positions of the element's nodes are given by @a elfun. */
1690 virtual void ComputeElementTargets(int e_id, const FiniteElement &fe,
1691 const IntegrationRule &ir,
1692 const Vector &elfun,
1693 DenseTensor &Jtr) const;
1694
1695 /** @brief Computes reference-to-target transformation Jacobians for all
1696 quadrature points in all elements.
1697
1698 @param[in] fes The nodal FE space
1699 @param[in] ir The quadrature rule to use for all elements
1700 @param[in] xe E-vector with the current physical coordinates/positions;
1701 this parameter is used only when needed by the target
1702 constructor, see UsesPhysicalCoordinates()
1703 @param[out] Jtr The computed ref->target Jacobian matrices. */
1704 virtual void ComputeAllElementTargets(const FiniteElementSpace &fes,
1705 const IntegrationRule &ir,
1706 const Vector &xe,
1707 DenseTensor &Jtr) const;
1708
1709 virtual void ComputeElementTargetsGradient(const IntegrationRule &ir,
1710 const Vector &elfun,
1712 DenseTensor &dJtr) const;
1713};
1714
1716{
1717public:
1719
1720 /** @brief Evaluate the derivative of the matrix coefficient with respect to
1721 @a comp in the element described by @a T at the point @a ip, storing the
1722 result in @a K. */
1724 const IntegrationPoint &ip, int comp) = 0;
1725
1727};
1728
1730{
1731protected:
1732 // Analytic target specification.
1736
1737public:
1739 : TargetConstructor(ttype),
1740 scalar_tspec(NULL), vector_tspec(NULL), matrix_tspec(NULL)
1741 { uses_phys_coords = true; }
1742
1743 virtual void SetAnalyticTargetSpec(Coefficient *sspec,
1744 VectorCoefficient *vspec,
1745 TMOPMatrixCoefficient *mspec);
1746
1747 /** @brief Given an element and quadrature rule, computes ref->target
1748 transformation Jacobians for each quadrature point in the element.
1749 The physical positions of the element's nodes are given by @a elfun. */
1750 void ComputeElementTargets(int e_id, const FiniteElement &fe,
1751 const IntegrationRule &ir,
1752 const Vector &elfun,
1753 DenseTensor &Jtr) const override;
1754
1756 const IntegrationRule &ir,
1757 const Vector &xe,
1758 DenseTensor &Jtr) const override;
1759
1761 const Vector &elfun,
1763 DenseTensor &dJtr) const override;
1764};
1765
1766#ifdef MFEM_USE_MPI
1767class ParGridFunction;
1768#endif
1769
1771{
1772protected:
1773 // Discrete target specification.
1774 // Data is owned, updated by UpdateTargetSpecification.
1776 Vector tspec; //eta(x) - we enforce Ordering::byNODES
1779 Vector tspec_pert2h; //eta(x+2*h)
1780 Vector tspec_pertmix; //eta(x+h,y+h)
1781 // The order inside these perturbation vectors (e.g. in 2D) is
1782 // eta1(x+h,y), eta2(x+h,y) ... etan(x+h,y), eta1(x,y+h), eta2(x,y+h) ...
1783 // same for tspec_pert2h and tspec_pertmix.
1784
1785 // DenseMatrix to hold target_spec values for the (children of the)
1786 // element being refined to consider for h-refinement.
1788 // Vector to hold the target_spec values for the coarse version of the
1789 // current mesh. Used for derefinement decision with hr-adaptivity.
1791
1792 // Components of Target Jacobian at each quadrature point of an element. This
1793 // is required for computation of the derivative using chain rule.
1795
1796 // Note: do not use the Nodes of this space as they may not be on the
1797 // positions corresponding to the values of tspec.
1799 FiniteElementSpace *coarse_tspec_fesv; //not owned, derefinement FESpace
1800 GridFunction *tspec_gf; //owned, uses tspec and tspec_fes
1801 // discrete adaptivity
1802#ifdef MFEM_USE_MPI
1803 ParFiniteElementSpace *ptspec_fesv; //owned, needed for derefinement to
1804 // get update operator.
1805 ParGridFunction *tspec_pgf; // similar to tspec_gf
1806#endif
1807
1810
1811 // These flags can be used by outside functions to avoid recomputing the
1812 // tspec and tspec_perth fields again on the same mesh.
1814
1815 // Evaluation of the discrete target specification on different meshes.
1816 // Owned.
1818
1819 void SetDiscreteTargetBase(const GridFunction &tspec_);
1820 void SetTspecAtIndex(int idx, const GridFunction &tspec_);
1822#ifdef MFEM_USE_MPI
1823 void SetTspecAtIndex(int idx, const ParGridFunction &tspec_);
1825#endif
1826
1827public:
1829 : TargetConstructor(ttype),
1830 ncomp(0),
1831 sizeidx(-1), skewidx(-1), aspectratioidx(-1), orientationidx(-1),
1834 tspec_fesv(NULL), coarse_tspec_fesv(NULL), tspec_gf(NULL),
1835#ifdef MFEM_USE_MPI
1836 ptspec_fesv(NULL), tspec_pgf(NULL),
1837#endif
1838 amr_el(-1), lim_min_size(-0.1),
1839 good_tspec(false), good_tspec_grad(false), good_tspec_hess(false),
1840 adapt_eval(NULL) { }
1841
1842 virtual ~DiscreteAdaptTC();
1843
1844 /** @name Target specification methods.
1845 The following methods are used to specify geometric parameters of the
1846 targets when these parameters are given by discrete FE functions.
1847 Note that every GridFunction given to the Set methods must use a
1848 H1_FECollection of the same order. The number of components must
1849 correspond to the type of geometric parameter and dimension.
1850
1851 @param[in] tspec_ Input values of a geometric parameter. Note that
1852 the methods in this class support only functions that
1853 use H1_FECollection collection of the same order. */
1854 ///@{
1855 virtual void SetSerialDiscreteTargetSpec(const GridFunction &tspec_);
1856 virtual void SetSerialDiscreteTargetSize(const GridFunction &tspec_);
1857 virtual void SetSerialDiscreteTargetSkew(const GridFunction &tspec_);
1858 virtual void SetSerialDiscreteTargetAspectRatio(const GridFunction &tspec_);
1859 virtual void SetSerialDiscreteTargetOrientation(const GridFunction &tspec_);
1860#ifdef MFEM_USE_MPI
1861 virtual void SetParDiscreteTargetSpec(const ParGridFunction &tspec_);
1862 virtual void SetParDiscreteTargetSize(const ParGridFunction &tspec_);
1863 virtual void SetParDiscreteTargetSkew(const ParGridFunction &tspec_);
1864 virtual void SetParDiscreteTargetAspectRatio(const ParGridFunction &tspec_);
1865 virtual void SetParDiscreteTargetOrientation(const ParGridFunction &tspec_);
1866#endif
1867 ///@}
1868
1869 /// Used in combination with the Update methods to avoid extra computations.
1872
1873 /// Get one of the discrete fields from tspec.
1874 void GetDiscreteTargetSpec(GridFunction &tspec_, int idx);
1875 /// Get the FESpace associated with tspec.
1877 /// Get the entire tspec.
1879 /// Update all discrete fields based on tspec and update for AMR
1881
1882#ifdef MFEM_USE_MPI
1885#endif
1886
1887 /** Used to update the target specification after the mesh has changed. The
1888 new mesh positions are given by new_x. If @a reuse_flag is true,
1889 repeated calls won't do anything until ResetUpdateFlags() is called. */
1890 void UpdateTargetSpecification(const Vector &new_x, bool reuse_flag = false,
1891 int new_x_ordering=Ordering::byNODES);
1892
1893 void UpdateTargetSpecification(Vector &new_x, Vector &IntData,
1894 int new_x_ordering=Ordering::byNODES);
1895
1898 int nodenum, int idir,
1899 const Vector &IntData);
1900
1902
1903 /** Used for finite-difference based computations. Computes the target
1904 specifications after a mesh perturbation in x or y direction.
1905 If @a reuse_flag is true, repeated calls won't do anything until
1906 ResetUpdateFlags() is called. */
1908 bool reuse_flag = false,
1909 int x_ordering = Ordering::byNODES);
1910 /** Used for finite-difference based computations. Computes the target
1911 specifications after two mesh perturbations in x and/or y direction.
1912 If @a reuse_flag is true, repeated calls won't do anything until
1913 ResetUpdateFlags() is called. */
1915 bool reuse_flag = false,
1916 int x_ordering = Ordering::byNODES);
1917
1919 {
1920 if (adapt_eval) { delete adapt_eval; }
1921 adapt_eval = ae;
1922 }
1923
1925 {
1926 return adapt_eval;
1927 }
1928
1932
1933 /** @brief Given an element and quadrature rule, computes ref->target
1934 transformation Jacobians for each quadrature point in the element.
1935 The physical positions of the element's nodes are given by @a elfun.
1936 Note that this function assumes that UpdateTargetSpecification() has
1937 been called with the position vector corresponding to @a elfun. */
1938 void ComputeElementTargets(int e_id, const FiniteElement &fe,
1939 const IntegrationRule &ir,
1940 const Vector &elfun,
1941 DenseTensor &Jtr) const override;
1942
1944 const IntegrationRule &ir,
1945 const Vector &xe,
1946 DenseTensor &Jtr) const override;
1947
1949 const Vector &elfun,
1951 DenseTensor &dJtr) const override;
1952
1953 // Generates tspec_vals for target construction using intrule
1954 // Used for the refinement component in hr-adaptivity.
1955 void SetTspecFromIntRule(int e_id, const IntegrationRule &intrule);
1956
1957 // Targets based on discrete functions can result in invalid (negative)
1958 // size at the quadrature points. This method can be used to set a
1959 // minimum target size.
1960 void SetMinSizeForTargets(real_t min_size_) { lim_min_size = min_size_; }
1961
1962 /// Computes target specification data with respect to the coarse FE space.
1964
1965 // Reset refinement data associated with h-adaptivity component.
1967 {
1969 amr_el = -1;
1970 }
1971
1972 // Reset derefinement data associated with h-adaptivity component.
1974 {
1976 coarse_tspec_fesv = NULL;
1977 }
1978
1979 // Used to specify the fine element for determining energy of children of a
1980 // parent element.
1981 void SetRefinementSubElement(int amr_el_) { amr_el = amr_el_; }
1982};
1983
1984class TMOPNewtonSolver;
1985
1986/** @brief A TMOP integrator class based on any given TMOP_QualityMetric and
1987 TargetConstructor.
1988
1989 Represents $ \int W(Jpt) dx $ over a target zone, where W is the
1990 metric's strain energy density function, and Jpt is the Jacobian of the
1991 target->physical coordinates transformation. The virtual target zone is
1992 defined by the TargetConstructor. */
1994{
1995protected:
1996 friend class TMOPNewtonSolver;
1998 friend class TMOPEnergyPA2D;
1999 friend class TMOPEnergyPA3D;
2002 friend class TMOPAddMultPA2D;
2003 friend class TMOPAddMultPA3D;
2004
2005 // Initial positions of the mesh nodes. Not owned. The pointer is set at the
2006 // start of the solve by TMOPNewtonSolver::Mult(), and unset at the end.
2007 // When x_0 == nullptr, the integrator works on the mesh positions.
2008 // When x_0 != nullptr, the integrator works on the displacements.
2009 // TODO in MFEM-5.0 make it always work with displacements.
2011 bool periodic = false;
2012
2015 const TargetConstructor *targetC; // not owned
2016
2017 // Custom integration rules.
2021
2022 // Weight Coefficient multiplying the quality metric term.
2023 Coefficient *metric_coeff; // not owned, if NULL -> metric_coeff is 1.
2024 // Normalization factor for the metric term.
2026
2027 // Nodes and weight Coefficient used for "limiting" the TMOP_Integrator.
2028 // These are both NULL when there is no limiting.
2029 // The class doesn't own lim_nodes0 and lim_coeff.
2032 // Limiting reference distance. Not owned.
2034 // Limiting function. Owned.
2036 // Normalization factor for the limiting term.
2038
2039 // Adaptive limiting.
2040 const GridFunction *adapt_lim_gf0; // Not owned.
2041#ifdef MFEM_USE_MPI
2043#endif
2044 GridFunction *adapt_lim_gf; // Owned. Updated by adapt_lim_eval.
2047
2048 // Surface fitting.
2049 const Array<bool> *surf_fit_marker; // Not owned. Nodes to fit.
2050 Coefficient *surf_fit_coeff; // Not owned. Fitting term scaling.
2051 // Fitting to a discrete level set.
2052 GridFunction *surf_fit_gf; // Owned. Updated by surf_fit_eval.
2054 // Fitting to given physical positions.
2055 TMOP_QuadraticLimiter *surf_fit_limiter; // Owned. Created internally.
2056 const GridFunction *surf_fit_pos; // Not owned. Positions to fit.
2057 real_t surf_fit_normal; // Normalization factor.
2058 GridFunction *surf_fit_grad, *surf_fit_hess; // Owned. Created internally.
2060 Array<int> surf_fit_dof_count; // Number of dofs per node.
2061 Array<int> surf_fit_marker_dof_index; // Indices of nodes to fit.
2062
2064
2065 // Parameters for FD-based Gradient & Hessian calculation.
2069 // Specifies that ComputeElementTargets is being called by a FD function.
2070 // It's used to skip terms that have exact derivative calculations.
2072 // Compute the exact action of the Integrator (includes derivative of the
2073 // target with respect to spatial position)
2075
2076 Array <Vector *> ElemDer; //f'(x)
2077 Array <Vector *> ElemPertEnergy; //f(x+h)
2078
2079 // Jrt: the inverse of the ref->target Jacobian, Jrt = Jtr^{-1}.
2080 // Jpr: the ref->physical transformation Jacobian, Jpr = PMatI^t DS.
2081 // Jpt: the target->physical transformation Jacobian, Jpt = Jpr Jrt.
2082 // P: represents dW_d(Jtp) (dim x dim).
2083 // DSh: gradients of reference shape functions (dof x dim).
2084 // DS: gradients of the shape functions in the target configuration,
2085 // DS = DSh Jrt (dof x dim).
2086 // PMatI: current coordinates of the nodes (dof x dim).
2087 // PMat0: reshaped view into the local element contribution to the operator
2088 // output - the result of AssembleElementVector() (dof x dim).
2090
2091 // PA extension
2092 // ------------
2093 // Jtr: all ref->target Jacobians, (dim x dim) Q-Vector as DenseTensor.
2094 // updated when needed, based on Jtr_needs_update.
2095 //
2096 // E: Q-vector for TMOP-energy
2097 // Used as temporary storage when the total energy is computed.
2098 // O: Q-Vector of 1.0, used to compute sums using the dot product kernel.
2099 // X0: E-vector for initial nodal coordinates.
2100 // Does not change during the TMOP iteration.
2101 // XL: E-vector for nodal coordinates used for limiting.
2102 // Does not change during the TMOP iteration.
2103 // H: Q-Vector for Hessian associated with the metric term.
2104 // Updated by every call to PANonlinearFormExtension::GetGradient().
2105 // C0: Q-Vector for spatial weight used for the limiting term.
2106 // Updated when the mesh nodes change.
2107 // LD: E-Vector constructed using limiting distance grid function (delta).
2108 // Does not change during the TMOP iteration.
2109 // H0: Q-Vector for Hessian associated with the limiting term.
2110 // Updated by every call to PANonlinearFormExtension::GetGradient().
2111 // MC: Q-Vector for the metric Coefficient.
2112 // Updated when the mesh nodes change.
2113 //
2114 // maps: Dof2Quad map for fes associated with the nodal coordinates.
2115 // maps_lim: Dof2Quad map for fes associated with the limiting dist GridFunc.
2116 //
2117 // Jtr_debug_grad
2118 // We keep track if Jtr was set by AssembleGradPA() in Jtr_debug_grad: it
2119 // is set to true by AssembleGradPA(); any other call to
2120 // ComputeAllElementTargets() will set the flag to false. This flag will
2121 // be used to check that Jtr is the one set by AssembleGradPA() when
2122 // performing operations with the gradient like AddMultGradPA() and
2123 // AssembleGradDiagonalPA().
2124 //
2125 // TODO:
2126 // * Merge LD, C0, H0 into one scalar Q-vector
2127 struct
2128 {
2130 int dim, ne, nq;
2132 mutable bool Jtr_needs_update;
2133 mutable bool Jtr_debug_grad;
2134 mutable Vector E, O, X0, XL, H, C0, LD, H0, MC;
2136 const DofToQuad *maps_lim = nullptr;
2141
2143 real_t &metric_energy, real_t &lim_energy);
2144
2147 const Vector &d_el, Vector &elvect);
2148
2151 const Vector &d_el, DenseMatrix &elmat);
2152
2155 const Vector &d_el, Vector &elvect);
2156
2157 // Assumes that AssembleElementVectorFD has been called.
2158 void AssembleElementGradFD(const FiniteElement &el,
2160 const Vector &d_el, DenseMatrix &elmat);
2161
2164 const IntegrationRule &ir,
2165 const Vector &weights, DenseMatrix &mat);
2168 const IntegrationRule &ir,
2169 const Vector &weights, DenseMatrix &m);
2170
2171 // First derivative of the surface fitting term.
2172 void AssembleElemVecSurfFit(const FiniteElement &el_x,
2174 DenseMatrix &mat);
2175
2176 // Second derivative of the surface fitting term.
2177 void AssembleElemGradSurfFit(const FiniteElement &el_x,
2179 DenseMatrix &mat);
2180
2183 Vector &d_el, const int nodenum, const int idir,
2184 const real_t baseenergy, bool update_stored);
2185
2186 /** @brief Determines the perturbation, h, for FD-based approximation. */
2187 void ComputeFDh(const Vector &d, const FiniteElementSpace &fes);
2188 void ComputeMinJac(const Vector &x, const FiniteElementSpace &fes);
2189
2191 const FiniteElementSpace &d_fes);
2192
2194 {
2195 lim_nodes0 = NULL; lim_coeff = NULL; lim_dist = NULL;
2196 delete lim_func; lim_func = NULL;
2197 }
2198
2200 {
2201 if (IntegRules)
2202 {
2203 return IntegRules->Get(el.GetGeomType(), integ_order);
2204 }
2205 return (IntRule) ? *IntRule
2206 /* */ : IntRules.Get(el.GetGeomType(), 2*el.GetOrder() + 3);
2207 }
2209 {
2210 // TODO the energy most likely needs less integration points.
2211 return EnergyIntegrationRule(el);
2212 }
2214 {
2215 // TODO the action and energy most likely need less integration points.
2216 return EnergyIntegrationRule(el);
2217 }
2218
2219 // Auxiliary PA methods
2220 void AssembleGradPA_2D(const Vector&) const;
2221 void AssembleGradPA_3D(const Vector&) const;
2222 void AssembleGradPA_C0_2D(const Vector&) const;
2223 void AssembleGradPA_C0_3D(const Vector&) const;
2224
2225 void GetLocalStateEnergyPA_2D(const Vector &x, real_t &energy) const;
2226 void GetLocalStateEnergyPA_3D(const Vector&, real_t &energy) const;
2230 real_t &met_energy,
2231 real_t &lim_energy) const;
2233 real_t &met_energy,
2234 real_t &lim_energy) const;
2235
2236 void AddMultPA_2D(const Vector&, Vector&) const;
2237 void AddMultPA_3D(const Vector&, Vector&) const;
2238 void AddMultPA_C0_2D(const Vector&, Vector&) const;
2239 void AddMultPA_C0_3D(const Vector&, Vector&) const;
2240
2241 void AddMultGradPA_2D(const Vector&, Vector&) const;
2242 void AddMultGradPA_3D(const Vector&, Vector&) const;
2243 void AddMultGradPA_C0_2D(const Vector&, Vector&) const;
2244 void AddMultGradPA_C0_3D(const Vector&, Vector&) const;
2245
2250
2251 void AssemblePA_Limiting();
2252 void ComputeAllElementTargets(const Vector &xe = Vector()) const;
2253 // Updates the Q-vectors for the metric_coeff and lim_coeff, based on the
2254 // new physical positions of the quadrature points.
2255 void UpdateCoefficientsPA(const Vector &d_loc);
2256
2257 // Compute Min(Det(Jpt)) in the mesh, does not reduce over MPI.
2259 // Compute Max(mu_hat) for the TMOP_WorstCaseUntangleOptimizer_Metric,
2260 // does not reduce over MPI.
2262 const FiniteElementSpace &fes);
2263
2264 // Remaps the internal surface fitting grid function object at provided
2265 // locations.
2267 int new_x_ordering);
2268public:
2269 /** @param[in] m TMOP_QualityMetric for r-adaptivity (not owned).
2270 @param[in] tc Target-matrix construction algorithm to use (not owned).
2271 @param[in] hm TMOP_QualityMetric for h-adaptivity (not owned). */
2274 : x_0(nullptr), h_metric(hm), metric(m), targetC(tc), IntegRules(NULL),
2275 integ_order(-1), metric_coeff(NULL), metric_normal(1.0),
2276 lim_nodes0(NULL), lim_coeff(NULL),
2277 lim_dist(NULL), lim_func(NULL), lim_normal(1.0),
2278 adapt_lim_gf0(NULL), adapt_lim_gf(NULL), adapt_lim_coeff(NULL),
2279 adapt_lim_eval(NULL),
2280 surf_fit_marker(NULL), surf_fit_coeff(NULL),
2281 surf_fit_gf(NULL), surf_fit_eval(NULL),
2282 surf_fit_limiter(NULL), surf_fit_pos(NULL),
2283 surf_fit_normal(1.0), surf_fit_grad(NULL), surf_fit_hess(NULL),
2285 discr_tc(dynamic_cast<DiscreteAdaptTC *>(tc)),
2286 fdflag(false), fd_h_scale(1.0e3), fd_call_flag(false), exact_action(false)
2287 { PA.enabled = false; }
2288
2291
2293
2294 /// Release the device memory of large PA allocations. This will copy device
2295 /// memory back to the host before releasing.
2296 void ReleasePADeviceMemory(bool copy_to_host = true);
2297
2298 /// Prescribe a set of integration rules; relevant for mixed meshes.
2299 /** This function has priority over SetIntRule(), if both are called. */
2300 void SetIntegrationRules(IntegrationRules &irules, int order)
2301 {
2302 IntegRules = &irules;
2303 integ_order = order;
2304 }
2305
2306 /// As the integrator operates on mesh displacements, this function is needed
2307 /// to set the initial mesh positions. For periodic meshes, the function is
2308 /// called with L2 positions.
2309 /// Called with nullptr to unset the x_0 after the problem is solved.
2310 // TODO should not be public in mfem 5.0.
2311 void SetInitialMeshPos(const GridFunction *x0);
2312
2313 /// The TMOP integrals can be computed over the reference element or the
2314 /// target elements. This function is used to switch between the two options.
2315 /// By default integration is performed over the target elements.
2316 void IntegrateOverTarget(bool integ_over_target_)
2317 {
2318 MFEM_VERIFY(metric_normal == 1.0 && lim_normal == 1.0,
2319 "This function must be called before EnableNormalization, as "
2320 "the normalization computations must know how to integrate.");
2321
2322 integ_over_target = integ_over_target_;
2323 }
2324
2325 /// Sets a scaling Coefficient for the quality metric term of the integrator.
2326 /** With this addition, the integrator becomes
2327 $ \int w1 W(Jpt) dx $.
2328
2329 Note that the Coefficient is evaluated in the physical configuration and
2330 not in the target configuration which may be undefined. */
2332
2333 /** @brief Limiting of the mesh displacements (general version).
2334
2335 Adds the term $ \int w_0 f(x, x_0, d) dx $, where f is a measure of
2336 the displacement between x and x_0, given the max allowed displacement d.
2337
2338 @param[in] n0 Original mesh node coordinates (x0 above).
2339 @param[in] dist Allowed displacement in physical space (d above).
2340 @param[in] w0 Coefficient scaling the limiting integral.
2341 @param[in] lfunc TMOP_LimiterFunction defining the function f. If
2342 NULL, a TMOP_QuadraticLimiter will be used. The
2343 TMOP_Integrator assumes ownership of this pointer. */
2344 void EnableLimiting(const GridFunction &n0, const GridFunction &dist,
2345 Coefficient &w0, TMOP_LimiterFunction *lfunc = NULL);
2346
2347 /** @brief Adds a limiting term to the integrator with limiting distance
2348 function (@a dist in the general version of the method) equal to 1. */
2349 void EnableLimiting(const GridFunction &n0, Coefficient &w0,
2350 TMOP_LimiterFunction *lfunc = NULL);
2351
2352 /** @brief Restriction of the node positions to certain regions.
2353
2354 Adds the term $ \int c (z(x) - z_0(x_0))^2 $, where z0(x0) is a given
2355 function on the starting mesh, and z(x) is its image on the new mesh.
2356 Minimizing this term means that a node at x0 is allowed to move to a
2357 position x(x0) only if z(x) ~ z0(x0).
2358 Such term can be used for tangential mesh relaxation.
2359
2360 @param[in] z0 Function z0 that controls the adaptive limiting.
2361 @param[in] coeff Coefficient c for the above integral.
2362 @param[in] ae AdaptivityEvaluator to compute z(x) from z0(x0). */
2363 void EnableAdaptiveLimiting(const GridFunction &z0, Coefficient &coeff,
2365#ifdef MFEM_USE_MPI
2366 /// Parallel support for adaptive limiting.
2367 void EnableAdaptiveLimiting(const ParGridFunction &z0, Coefficient &coeff,
2369#endif
2370
2371 /** @brief Fitting of certain DOFs to the zero level set of a function.
2372
2373 Having a level set function s0(x0) on the starting mesh, and a set of
2374 marked nodes (or DOFs), we move these nodes to the zero level set of s0.
2375 If s(x) is the image of s0(x0) on the current mesh, this function adds to
2376 the TMOP functional the term $ \int c \bar{s}(x))^2 $, where
2377 $\bar{s}(x)$ is the restriction of s(x) on the aligned DOFs.
2378 Minimizing this term means that a marked node at x0 is allowed to move to
2379 a position x(x0) only if s(x) ~ 0.
2380 Such term can be used for surface fitting and tangential relaxation.
2381
2382 @param[in] s0 The level set function on the initial mesh.
2383 @param[in] smarker Indicates which DOFs will be aligned.
2384 @param[in] coeff Coefficient c for the above integral.
2385 @param[in] ae AdaptivityEvaluator to compute s(x) from s0(x0). */
2386 void EnableSurfaceFitting(const GridFunction &s0,
2387 const Array<bool> &smarker, Coefficient &coeff,
2389
2390#ifdef MFEM_USE_MPI
2391 /// Parallel support for surface fitting to the zero level set of a function.
2392 /// Here, we add two optional inputs: @a aegrad and @a aehess. When provided,
2393 /// the first and second derivative of the input level set are computed on
2394 /// the initial mesh, and @a aegrad and @a aehess are used to remap grad_s(x)
2395 /// from grad_s0(x0) and hess_s(x) from hess_s0(x0), respectively.
2397 const Array<bool> &smarker, Coefficient &coeff,
2399 AdaptivityEvaluator *aegrad = NULL,
2400 AdaptivityEvaluator *aehess = NULL);
2401
2402 /** @brief Fitting of certain DOFs in the current mesh to the zero level set
2403 of a function defined on another (finer) source mesh.
2404
2405 Having a level set function s_bg(x_bg) on a source/background mesh,
2406 a set of marked nodes (or DOFs) in the current mesh, we move the marked
2407 nodes to the zero level set of s_bg. This functionality is used for
2408 surface fitting and tangential relaxation.
2409
2410 @param[in] s_bg The level set function on the background mesh.
2411 @param[in] s0 The level set function (automatically) interpolated
2412 on the initial mesh.
2413 @param[in] smarker Marker for aligned DOFs in the current mesh.
2414 @param[in] coeff Coefficient c for the fitting penalty term.
2415 @param[in] ae Interpolates s(x) from s_bg(x_bg).
2416 @param[in] s_bg_grad Gradient of s_bg on the background mesh.
2417 @param[in] s0_grad Gradient of s0 on the initial mesh.
2418 @param[in] age Interpolates s_grad(x) from s_bg_grad(x_bg).
2419 @param[in] s_bg_hess Hessian of s(x) on the background mesh.
2420 @param[in] s0_hess Hessian of s0 on the initial mesh.
2421 @param[in] ahe Interpolates s_hess(x) from s_bg_hess(x_bg).
2422 See the pmesh-fitting miniapp for details on usage. */
2424 ParGridFunction &s0,
2425 const Array<bool> &smarker,
2426 Coefficient &coeff,
2428 const ParGridFunction &s_bg_grad,
2429 ParGridFunction &s0_grad,
2431 const ParGridFunction &s_bg_hess,
2432 ParGridFunction &s0_hess,
2433 AdaptivityEvaluator &ahe);
2434#endif
2435 /** @brief Fitting of certain DOFs to given positions in physical space.
2436
2437 Having a set S of marked nodes (or DOFs) and their target positions in
2438 physical space x_t, we move these nodes to the target positions during
2439 the optimization process.
2440 This function adds to the TMOP functional the term
2441 $ \sum_{i \in S} c \frac{1}{2} (x_i - x_{t,i})^2 $,
2442 where $c$ corresponds to @a coeff below and is evaluated at the
2443 DOF locations.
2444
2445 @param[in] pos The desired positions for the mesh nodes.
2446 @param[in] smarker Indicates which DOFs will be aligned.
2447 @param[in] coeff Coefficient c for the above integral. */
2448 void EnableSurfaceFitting(const GridFunction &pos,
2449 const Array<bool> &smarker, Coefficient &coeff);
2450 void GetSurfaceFittingErrors(const Vector &d_loc,
2451 real_t &err_avg, real_t &err_max);
2453 {
2454 return surf_fit_gf != NULL || surf_fit_pos != NULL;
2455 }
2456
2457 /// Update the original/reference nodes used for limiting.
2458 void SetLimitingNodes(const GridFunction &n0) { lim_nodes0 = &n0; }
2459
2460 /** @brief Computes the integral of W(Jacobian(Trt)) over a target zone.
2461 @param[in] el Type of FiniteElement (TMOP assumes H1 elements).
2462 @param[in] T Mesh element transformation.
2463 @param[in] d_el Physical displacement of the zone w.r.t. x_0. */
2466 const Vector &d_el) override;
2467
2468 /** @brief Computes the mean of the energies of the given element's children.
2469
2470 In addition to the inputs for GetElementEnergy, this function requires an
2471 IntegrationRule to be specified that will give the decomposition of the
2472 given element based on the refinement type being considered. */
2475 const Vector &elfun,
2476 const IntegrationRule &irule);
2477
2478 /// This function is similar to GetElementEnergy, but ignores components
2479 /// such as limiting etc. to compute the element energy.
2482 const Vector &elfun);
2483
2484 /// First defivative of GetElementEnergy() w.r.t. each local H1 DOF.
2485 void AssembleElementVector(const FiniteElement &el,
2487 const Vector &d_el, Vector &elvect) override;
2488
2489 /// Second derivative of GetElementEnergy() w.r.t. each local H1 DOF.
2490 void AssembleElementGrad(const FiniteElement &el,
2492 const Vector &d_el, DenseMatrix &elmat) override;
2493
2495
2497#ifdef MFEM_USE_MPI
2499#endif
2500
2501 // PA extension
2503 void AssemblePA(const FiniteElementSpace&) override;
2504
2505 void AssembleGradPA(const Vector &, const FiniteElementSpace &) override;
2506
2507 real_t GetLocalStateEnergyPA(const Vector&) const override;
2508
2509 void AddMultPA(const Vector&, Vector&) const override;
2510
2511 void AddMultGradPA(const Vector&, Vector&) const override;
2512
2513 void AssembleGradDiagonalPA(Vector&) const override;
2514
2516
2517 /** @brief Computes the normalization factors of the metric and limiting
2518 integrals using the mesh position given by @a x. */
2519 void EnableNormalization(const GridFunction &x);
2520#ifdef MFEM_USE_MPI
2522#endif
2523
2524 /** @brief Get the normalization factors of the metric */
2528
2529 /** @brief Enables FD-based approximation and computes dx. */
2531#ifdef MFEM_USE_MPI
2533#endif
2534
2535 void SetFDhScale(real_t scale) { fd_h_scale = scale; }
2536 bool GetFDFlag() const { return fdflag; }
2537 real_t GetFDh() const { return fd_h; }
2538
2539 /** @brief Flag to control if exact action of Integration is effected. */
2540 void SetExactActionFlag(bool flag_) { exact_action = flag_; }
2541
2542 /// Update the surface fitting weight as surf_fit_coeff *= factor;
2544
2545 /// Get the surface fitting weight.
2547
2548 /// Computes quantiles needed for UntangleMetrics. Note that in parallel,
2549 /// the ParFiniteElementSpace must be passed as argument for consistency
2550 /// across MPI ranks.
2552 const FiniteElementSpace &fes);
2553};
2554
2556{
2557protected:
2558 friend class TMOPNewtonSolver;
2559
2560 // Integrators in the combination. Owned.
2562
2564 {
2565 for (int i = 0; i < tmopi.Size(); i++)
2566 { tmopi[i]->SetInitialMeshPos(x0); }
2567 }
2568
2569public:
2571
2573 {
2574 for (int i = 0; i < tmopi.Size(); i++) { delete tmopi[i]; }
2575 }
2576
2577 /// Adds a new TMOP_Integrator to the combination.
2578 void AddTMOPIntegrator(TMOP_Integrator *ti) { tmopi.Append(ti); }
2579
2581
2582 /// Adds the limiting term to the first integrator. Disables it for the rest.
2583 void EnableLimiting(const GridFunction &n0, const GridFunction &dist,
2584 Coefficient &w0, TMOP_LimiterFunction *lfunc = NULL);
2585
2586 /** @brief Adds the limiting term to the first integrator. Disables it for
2587 the rest (@a dist in the general version of the method) equal to 1. */
2588 void EnableLimiting(const GridFunction &n0, Coefficient &w0,
2589 TMOP_LimiterFunction *lfunc = NULL);
2590
2591 /// Update the original/reference nodes used for limiting.
2592 void SetLimitingNodes(const GridFunction &n0);
2593
2596 const Vector &elfun) override;
2597 void AssembleElementVector(const FiniteElement &el,
2599 const Vector &elfun, Vector &elvect) override;
2600 void AssembleElementGrad(const FiniteElement &el,
2602 const Vector &elfun, DenseMatrix &elmat) override;
2603
2606 const Vector &elfun,
2607 const IntegrationRule &irule);
2608
2611 const Vector &elfun);
2612
2613 /// Normalization factor that considers all integrators in the combination.
2614 void EnableNormalization(const GridFunction &x);
2615#ifdef MFEM_USE_MPI
2617#endif
2618
2619 // PA extension
2621 void AssemblePA(const FiniteElementSpace&) override;
2622 void AssembleGradPA(const Vector&, const FiniteElementSpace&) override;
2623 real_t GetLocalStateEnergyPA(const Vector&) const override;
2624 void AddMultPA(const Vector&, Vector&) const override;
2625 void AddMultGradPA(const Vector&, Vector&) const override;
2626 void AssembleGradDiagonalPA(Vector&) const override;
2627};
2628
2629/// Interpolates the @a metric's values at the nodes of @a metric_gf.
2630/** Assumes that @a metric_gf's FiniteElementSpace is initialized. */
2631void InterpolateTMOP_QualityMetric(TMOP_QualityMetric &metric,
2632 const TargetConstructor &tc,
2633 const Mesh &mesh, GridFunction &metric_gf);
2634}
2635
2636#endif
virtual ~AdaptivityEvaluator()
Definition tmop.cpp:3717
virtual void SetNewFieldFESpace(const FiniteElementSpace &fes)=0
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)=0
virtual void ComputeAtGivenPositions(const Vector &positions, Vector &values, int p_ordering=Ordering::byNODES)=0
Using the source mesh and field given by SetInitialField(), compute corresponding values at specified...
void SetSerialMetaInfo(const Mesh &m, const FiniteElementSpace &f)
Definition tmop.cpp:3686
FiniteElementSpace * fes
Definition tmop.hpp:1513
void SetParMetaInfo(const ParMesh &m, const ParFiniteElementSpace &f)
Parallel version of SetSerialMetaInfo.
Definition tmop.cpp:3697
virtual void ComputeAtNewPosition(const Vector &new_mesh_nodes, Vector &new_field, int nodes_ordering=Ordering::byNODES)=0
Perform field transfer between the original and a new mesh. The source mesh and field are given by Se...
ParFiniteElementSpace * pfes
Definition tmop.hpp:1518
AnalyticAdaptTC(TargetType ttype)
Definition tmop.hpp:1738
Coefficient * scalar_tspec
Definition tmop.hpp:1733
virtual void SetAnalyticTargetSpec(Coefficient *sspec, VectorCoefficient *vspec, TMOPMatrixCoefficient *mspec)
Definition tmop.cpp:2603
TMOPMatrixCoefficient * matrix_tspec
Definition tmop.hpp:1735
void ComputeAllElementTargets(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const override
Computes reference-to-target transformation Jacobians for all quadrature points in all elements.
Definition pa.cpp:142
VectorCoefficient * vector_tspec
Definition tmop.hpp:1734
void ComputeElementTargets(int e_id, const FiniteElement &fe, const IntegrationRule &ir, const Vector &elfun, DenseTensor &Jtr) const override
Given an element and quadrature rule, computes ref->target transformation Jacobians for each quadratu...
Definition tmop.cpp:2612
void ComputeElementTargetsGradient(const IntegrationRule &ir, const Vector &elfun, IsoparametricTransformation &Tpr, DenseTensor &dJtr) const override
Definition tmop.cpp:2646
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:912
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:116
void Diag(real_t c, int n)
Creates n x n diagonal matrix with diagonal elements c.
void Clear()
Delete the matrix data array (if owned) and reset the matrix state.
Definition densemat.hpp:106
Rank 3 tensor (array of matrices)
const Vector & GetTspecPert1H()
Definition tmop.hpp:1929
void UpdateTargetSpecification(const Vector &new_x, bool reuse_flag=false, int new_x_ordering=Ordering::byNODES)
Definition tmop.cpp:2923
virtual void SetSerialDiscreteTargetSkew(const GridFunction &tspec_)
Definition tmop.cpp:2841
virtual void SetParDiscreteTargetAspectRatio(const ParGridFunction &tspec_)
Definition tmop.cpp:2769
ParGridFunction * tspec_pgf
Definition tmop.hpp:1805
DenseTensor Jtrcomp
Definition tmop.hpp:1794
FiniteElementSpace * GetTSpecFESpace()
Get the FESpace associated with tspec.
Definition tmop.hpp:1876
virtual void SetSerialDiscreteTargetAspectRatio(const GridFunction &tspec_)
Definition tmop.cpp:2851
const AdaptivityEvaluator * GetAdaptivityEvaluator() const
Definition tmop.hpp:1924
void SetTspecDataForDerefinement(FiniteElementSpace *fes)
Computes target specification data with respect to the coarse FE space.
Definition tmop.cpp:3000
void ResetDerefinementTspecData()
Definition tmop.hpp:1973
const Vector & GetTspecPertMixH()
Definition tmop.hpp:1931
void ComputeElementTargetsGradient(const IntegrationRule &ir, const Vector &elfun, IsoparametricTransformation &Tpr, DenseTensor &dJtr) const override
Definition tmop.cpp:3223
void SetMinSizeForTargets(real_t min_size_)
Definition tmop.hpp:1960
virtual void SetParDiscreteTargetSpec(const ParGridFunction &tspec_)
Definition tmop.cpp:2789
AdaptivityEvaluator * adapt_eval
Definition tmop.hpp:1817
void FinalizeSerialDiscreteTargetSpec(const GridFunction &tspec_)
Definition tmop.cpp:2871
void FinalizeParDiscreteTargetSpec(const ParGridFunction &tspec_)
Definition tmop.cpp:2693
ParFiniteElementSpace * ptspec_fesv
Definition tmop.hpp:1803
GridFunction * tspec_gf
Definition tmop.hpp:1800
FiniteElementSpace * coarse_tspec_fesv
Definition tmop.hpp:1799
virtual ~DiscreteAdaptTC()
Definition tmop.cpp:3676
GridFunction * GetTSpecData()
Get the entire tspec.
Definition tmop.hpp:1878
virtual void SetSerialDiscreteTargetSize(const GridFunction &tspec_)
Definition tmop.cpp:2831
virtual void SetParDiscreteTargetSize(const ParGridFunction &tspec_)
Definition tmop.cpp:2749
const Vector & GetTspecPert2H()
Definition tmop.hpp:1930
void ComputeElementTargets(int e_id, const FiniteElement &fe, const IntegrationRule &ir, const Vector &elfun, DenseTensor &Jtr) const override
Given an element and quadrature rule, computes ref->target transformation Jacobians for each quadratu...
Definition tmop.cpp:3008
void UpdateAfterMeshTopologyChange()
Update all discrete fields based on tspec and update for AMR.
Definition tmop.cpp:2906
void SetAdaptivityEvaluator(AdaptivityEvaluator *ae)
Definition tmop.hpp:1918
virtual void SetSerialDiscreteTargetSpec(const GridFunction &tspec_)
Definition tmop.cpp:2917
virtual void SetParDiscreteTargetOrientation(const ParGridFunction &tspec_)
Definition tmop.cpp:2779
void ResetRefinementTspecData()
Definition tmop.hpp:1966
virtual void SetParDiscreteTargetSkew(const ParGridFunction &tspec_)
Definition tmop.cpp:2759
void ParUpdateAfterMeshTopologyChange()
Definition tmop.cpp:2718
void UpdateTargetSpecificationAtNode(const FiniteElement &el, ElementTransformation &T, int nodenum, int idir, const Vector &IntData)
Definition tmop.cpp:2943
void ResetUpdateFlags()
Used in combination with the Update methods to avoid extra computations.
Definition tmop.hpp:1870
void RestoreTargetSpecificationAtNode(ElementTransformation &T, int nodenum)
Definition tmop.cpp:2960
FiniteElementSpace * tspec_fesv
Definition tmop.hpp:1798
virtual void SetSerialDiscreteTargetOrientation(const GridFunction &tspec_)
Definition tmop.cpp:2861
void SetTspecAtIndex(int idx, const GridFunction &tspec_)
Definition tmop.cpp:2818
ParFiniteElementSpace * GetTSpecParFESpace()
Definition tmop.hpp:1883
void SetRefinementSubElement(int amr_el_)
Definition tmop.hpp:1981
void GetDiscreteTargetSpec(GridFunction &tspec_, int idx)
Get one of the discrete fields from tspec.
Definition tmop.cpp:2892
void UpdateHessianTargetSpecification(const Vector &x, real_t dx, bool reuse_flag=false, int x_ordering=Ordering::byNODES)
Definition tmop.cpp:3608
void ComputeAllElementTargets(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const override
Computes reference-to-target transformation Jacobians for all quadrature points in all elements.
void SetDiscreteTargetBase(const GridFunction &tspec_)
Definition tmop.cpp:2795
DenseMatrix tspec_refine
Definition tmop.hpp:1787
DiscreteAdaptTC(TargetType ttype)
Definition tmop.hpp:1828
void UpdateGradientTargetSpecification(const Vector &x, real_t dx, bool reuse_flag=false, int x_ordering=Ordering::byNODES)
Definition tmop.cpp:3571
void SetTspecFromIntRule(int e_id, const IntegrationRule &intrule)
Definition tmop.cpp:2974
Structure representing the matrices/tensors needed to evaluate (in reference space) the values,...
Definition fe_base.hpp:141
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:208
Abstract class for all finite elements.
Definition fe_base.hpp:247
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition fe_base.hpp:341
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition fe_base.hpp:334
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
Definition mesh.hpp:3076
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:32
Abstract class for hyperelastic models.
Class for integration point with weight.
Definition intrules.hpp:35
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Container class for integration rules.
Definition intrules.hpp:422
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
const IntegrationRule * IntRule
Auxiliary class for evaluating the 2x2 matrix invariants and their first and second derivatives.
Auxiliary class for evaluating the 3x3 matrix invariants and their first and second derivatives.
A standard isoparametric element transformation.
Definition eltrans.hpp:629
Base class for Matrix Coefficients that optionally depend on time and space.
Mesh data type.
Definition mesh.hpp:65
This class is used to express the local action of a general nonlinear finite element operator....
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
Abstract parallel finite element space.
Definition pfespace.hpp:31
Class for parallel grid function.
Definition pgridfunc.hpp:50
Class for parallel meshes.
Definition pmesh.hpp:34
virtual real_t GetDerefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
Definition tmop.cpp:5996
void AssembleGradPA(const Vector &, const FiniteElementSpace &) override
Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at th...
Definition tmop.cpp:6049
void AssembleGradDiagonalPA(Vector &) const override
Method for computing the diagonal of the gradient with partial assembly.
Definition tmop.cpp:6058
real_t GetLocalStateEnergyPA(const Vector &) const override
Compute the local (to the MPI rank) energy with partial assembly.
Definition tmop.cpp:6082
void AddMultPA(const Vector &, Vector &) const override
Method for partially assembled action.
Definition tmop.cpp:6066
void ParEnableNormalization(const ParGridFunction &x)
Definition tmop.cpp:6025
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:5910
void AssembleElementGrad(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, DenseMatrix &elmat) override
Assemble the local gradient matrix.
Definition tmop.cpp:5967
real_t GetElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun) override
Compute the local energy.
Definition tmop.cpp:5939
void AddMultGradPA(const Vector &, Vector &) const override
Method for partially assembled gradient action.
Definition tmop.cpp:6074
virtual real_t GetRefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, const IntegrationRule &irule)
Definition tmop.cpp:5983
void AssemblePA(const FiniteElementSpace &) override
Method defining partial assembly.
Definition tmop.cpp:6041
Array< TMOP_Integrator * > tmopi
Definition tmop.hpp:2561
void SetInitialMeshPos(const GridFunction *x0)
Definition tmop.hpp:2563
void AssembleElementVector(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, Vector &elvect) override
Perform the local action of the NonlinearFormIntegrator.
Definition tmop.cpp:5951
void EnableNormalization(const GridFunction &x)
Normalization factor that considers all integrators in the combination.
Definition tmop.cpp:6009
void AddTMOPIntegrator(TMOP_Integrator *ti)
Adds a new TMOP_Integrator to the combination.
Definition tmop.hpp:2578
const Array< TMOP_Integrator * > & GetTMOPIntegrators() const
Definition tmop.hpp:2580
void SetLimitingNodes(const GridFunction &n0)
Update the original/reference nodes used for limiting.
Definition tmop.cpp:5931
virtual ~TMOPMatrixCoefficient()
Definition tmop.hpp:1726
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 ...
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:2167
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:1235
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.hpp:1240
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:2193
void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const override
Definition tmop.cpp:2187
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:2181
2D barrier Size (V) metric (polyconvex).
Definition tmop.hpp:1253
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:2218
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.hpp:1260
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:1255
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:2230
void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const override
Definition tmop.cpp:2224
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:2204
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition tmop.hpp:1273
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:1275
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:2241
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:2267
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:2255
void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const override
Definition tmop.cpp:2261
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.hpp:1280
real_t GetGamma() const
Definition tmop.hpp:1391
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:1378
TMOP_QualityMetric * sh_metric
Definition tmop.hpp:1379
int Id() const override
Return the metric ID.
Definition tmop.hpp:1390
TMOP_AMetric_049(real_t gamma)
Definition tmop.hpp:1382
TMOP_QualityMetric * sk_metric
Definition tmop.hpp:1379
virtual ~TMOP_AMetric_049()
Definition tmop.hpp:1393
2D barrier Skew (Q) metric.
Definition tmop.hpp:1293
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:2292
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:2304
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:2278
void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const override
Definition tmop.cpp:2298
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.hpp:1300
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:1295
2D barrier Size+Skew (VQ) metric.
Definition tmop.hpp:1313
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:2329
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:2315
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.hpp:1322
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:1315
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:2341
void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const override
Definition tmop.cpp:2335
2D barrier Shape+Orientation (OS) metric (polyconvex).
Definition tmop.hpp:1335
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:2352
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:2366
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:2378
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:1337
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.hpp:1342
void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const override
Definition tmop.cpp:2372
2D barrier Shape+Size (VS) metric (polyconvex).
Definition tmop.hpp:1355
virtual ~TMOP_AMetric_126()
Definition tmop.hpp:1369
TMOP_QualityMetric * sh_metric
Definition tmop.hpp:1358
TMOP_AMetric_126(real_t gamma)
Definition tmop.hpp:1361
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:1357
TMOP_QualityMetric * sz_metric
Definition tmop.hpp:1358
void SetTargetJacobian(const DenseMatrix &Jtr_) override
Specify the reference-element -> target-element Jacobian matrix for the point of interest.
Definition tmop.hpp:128
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:595
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:537
void GetLocalEnergyPA_3D(const GridFunction &nodes, const TargetConstructor &tc, int m_index, real_t &energy, real_t &vol, const IntegrationRule &ir) const
void GetLocalEnergyPA_2D(const GridFunction &nodes, const TargetConstructor &tc, int m_index, real_t &energy, real_t &vol, const IntegrationRule &ir) const
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:547
void SetWeights(const Vector &weights)
Changes the weights of the metrics in the combination.
Definition tmop.hpp:179
void ComputeAvgMetrics(const GridFunction &nodes, const TargetConstructor &tc, Vector &averages, bool use_pa=false, const IntegrationRule *IntRule=nullptr) const
Definition tmop.cpp:639
AD2Type EvalW_AD2(const std::vector< AD2Type > &T, const std::vector< AD2Type > &W) const override
Second-derivative hook for AD-based computations.
Definition tmop.cpp:583
void ComputeBalancedWeights(const GridFunction &nodes, const TargetConstructor &tc, Vector &weights, bool use_pa=false, const IntegrationRule *IntRule=nullptr) const
Definition tmop.cpp:609
void GetWeights(Array< real_t > &weights) const
Definition tmop.hpp:176
AD1Type EvalW_AD1(const std::vector< AD1Type > &T, const std::vector< AD1Type > &W) const override
First-derivative hook for AD-based computations.
Definition tmop.cpp:571
Array< TMOP_QualityMetric * > tmop_q_arr
Definition tmop.hpp:118
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:527
virtual void AddQualityMetric(TMOP_QualityMetric *tq, real_t wt=1.0)
Definition tmop.hpp:122
void EvalPW(const DenseMatrix &Jpt, DenseMatrix &P) const override
Definition tmop.cpp:559
Exponential limiter function in TMOP_Integrator.
Definition tmop.hpp:1453
void Eval_d1(const Vector &x, const Vector &x0, real_t dist, Vector &d1) const override
Returns the gradient of the limiting function f(x, x0, d) with respect to x.
Definition tmop.hpp:1462
real_t Eval(const Vector &x, const Vector &x0, real_t dist) const override
Returns the limiting function, f(x, x0, d).
Definition tmop.hpp:1455
void Eval_d2(const Vector &x, const Vector &x0, real_t dist, DenseMatrix &d2) const override
Returns the Hessian of the limiting function f(x, x0, d) with respect to x.
Definition tmop.hpp:1473
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition tmop.hpp:1994
void AddMultGradPA(const Vector &, Vector &) const override
Method for partially assembled gradient action.
Definition pa.cpp:385
void SetExactActionFlag(bool flag_)
Flag to control if exact action of Integration is effected.
Definition tmop.hpp:2540
real_t ComputeUntanglerMaxMuBarrier(const Vector &x, const FiniteElementSpace &fes)
Definition tmop.cpp:5800
void UpdateAfterMeshPositionChange(const Vector &d, const FiniteElementSpace &d_fes)
Definition tmop.cpp:5633
DenseMatrix Jpr
Definition tmop.hpp:2089
void AssembleElemVecAdaptLim(const FiniteElement &el, IsoparametricTransformation &Tpr, const IntegrationRule &ir, const Vector &weights, DenseMatrix &mat)
Definition tmop.cpp:4824
void AssembleElemGradAdaptLim(const FiniteElement &el, IsoparametricTransformation &Tpr, const IntegrationRule &ir, const Vector &weights, DenseMatrix &m)
Definition tmop.cpp:4860
const DofToQuad * maps_lim
Definition tmop.hpp:2136
void GetNormalizationFactors(real_t &metric_normal, real_t &lim_normal, real_t &surf_fit_normal)
Get the normalization factors of the metric.
Definition tmop.cpp:5371
friend class TMOPAssembleGradPA2D
Definition tmop.hpp:2000
real_t GetLocalStateEnergyPA_C0_2D(const Vector &) const
TMOP_LimiterFunction * lim_func
Definition tmop.hpp:2035
void ComputeAllElementTargets(const Vector &xe=Vector()) const
Definition pa.cpp:170
void ParUpdateAfterMeshTopologyChange()
Definition tmop.cpp:4158
friend class TMOPEnergyPA2D
Definition tmop.hpp:1998
DenseMatrix PMatO
Definition tmop.hpp:2089
TMOP_QualityMetric * metric
Definition tmop.hpp:2014
const DofToQuad * maps
Definition tmop.hpp:2135
const ParGridFunction * adapt_lim_pgf0
Definition tmop.hpp:2042
void SetFDhScale(real_t scale)
Definition tmop.hpp:2535
void AddMultPA_C0_3D(const Vector &, Vector &) const
void GetSurfaceFittingErrors(const Vector &d_loc, real_t &err_avg, real_t &err_max)
Definition tmop.cpp:4074
real_t ComputeMinDetT(const Vector &x, const FiniteElementSpace &fes)
Definition tmop.cpp:5757
AdaptivityEvaluator * surf_fit_eval
Definition tmop.hpp:2053
void EnableNormalization(const GridFunction &x)
Computes the normalization factors of the metric and limiting integrals using the mesh position given...
Definition tmop.cpp:5349
void ComputeMinJac(const Vector &x, const FiniteElementSpace &fes)
Definition tmop.cpp:5470
void AssembleElementVectorExact(const FiniteElement &el, ElementTransformation &T, const Vector &d_el, Vector &elvect)
Definition tmop.cpp:4519
friend class TMOPAddMultPA2D
Definition tmop.hpp:2002
TMOP_QualityMetric * h_metric
Definition tmop.hpp:2013
void AssembleElementGradFD(const FiniteElement &el, ElementTransformation &T, const Vector &d_el, DenseMatrix &elmat)
Definition tmop.cpp:5217
DiscreteAdaptTC * discr_tc
Definition tmop.hpp:2063
void AssembleDiagonalPA_2D(Vector &) const
friend class TMOPEnergyPA3D
Definition tmop.hpp:1999
Coefficient * lim_coeff
Definition tmop.hpp:2031
GridFunction * surf_fit_gf
Definition tmop.hpp:2052
const GridFunction * lim_dist
Definition tmop.hpp:2033
void AssembleDiagonalPA_C0_3D(Vector &) const
Array< Vector * > ElemPertEnergy
Definition tmop.hpp:2077
const GridFunction * adapt_lim_gf0
Definition tmop.hpp:2040
void SetCoefficient(Coefficient &w1)
Sets a scaling Coefficient for the quality metric term of the integrator.
Definition tmop.hpp:2331
void AddMultGradPA_2D(const Vector &, Vector &) const
friend class TMOPAddMultPA3D
Definition tmop.hpp:2003
void ComputeUntangleMetricQuantiles(const Vector &d, const FiniteElementSpace &fes)
Definition tmop.cpp:5862
void AssembleGradPA_C0_2D(const Vector &) const
void EnableSurfaceFittingFromSource(const ParGridFunction &s_bg, ParGridFunction &s0, const Array< bool > &smarker, Coefficient &coeff, AdaptivityEvaluator &ae, const ParGridFunction &s_bg_grad, ParGridFunction &s0_grad, AdaptivityEvaluator &age, const ParGridFunction &s_bg_hess, ParGridFunction &s0_hess, AdaptivityEvaluator &ahe)
Fitting of certain DOFs in the current mesh to the zero level set of a function defined on another (f...
Definition tmop.cpp:3999
const TargetConstructor * targetC
Definition tmop.hpp:2015
AdaptivityEvaluator * surf_fit_eval_grad
Definition tmop.hpp:2059
void GetLocalStateEnergyPA_2D(const Vector &x, real_t &energy) const
const GridFunction * lim_nodes0
Definition tmop.hpp:2030
const IntegrationRule * ir
Definition tmop.hpp:2139
TMOP_QuadraticLimiter * surf_fit_limiter
Definition tmop.hpp:2055
void EnableFiniteDifferences(const GridFunction &x)
Enables FD-based approximation and computes dx.
Definition tmop.cpp:5703
virtual real_t GetRefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun, const IntegrationRule &irule)
Computes the mean of the energies of the given element's children.
Definition tmop.cpp:4350
real_t GetSurfaceFittingWeight()
Get the surface fitting weight.
Definition tmop.cpp:5338
void SetInitialMeshPos(const GridFunction *x0)
Definition tmop.cpp:3741
void AssemblePA(const FiniteElementSpace &) override
Method defining partial assembly.
Definition pa.cpp:234
DenseMatrix DSh
Definition tmop.hpp:2089
void AssembleGradDiagonalPA(Vector &) const override
Method for computing the diagonal of the gradient with partial assembly.
Definition pa.cpp:327
void AssembleDiagonalPA_C0_2D(Vector &) const
bool IsSurfaceFittingEnabled()
Definition tmop.hpp:2452
const GridFunction * surf_fit_pos
Definition tmop.hpp:2056
Array< int > surf_fit_marker_dof_index
Definition tmop.hpp:2061
const GeometricFactors * geom
Definition tmop.hpp:2137
real_t GetLocalStateEnergyPA_C0_3D(const Vector &) const
DenseTensor Jtr
Definition tmop.hpp:2131
void ParEnableNormalization(const ParGridFunction &x)
Definition tmop.cpp:5358
void AddMultGradPA_3D(const Vector &, Vector &) const
real_t GetElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &d_el) override
Computes the integral of W(Jacobian(Trt)) over a target zone.
Definition tmop.cpp:4171
Coefficient * adapt_lim_coeff
Definition tmop.hpp:2045
real_t GetFDh() const
Definition tmop.hpp:2537
const FiniteElementSpace * fes
Definition tmop.hpp:2138
void AssemblePA_Limiting()
Definition pa.cpp:58
const IntegrationRule & EnergyIntegrationRule(const FiniteElement &el) const
Definition tmop.hpp:2199
GridFunction * surf_fit_grad
Definition tmop.hpp:2058
GridFunction * surf_fit_hess
Definition tmop.hpp:2058
friend class TMOPAssembleGradPA3D
Definition tmop.hpp:2001
void ReleasePADeviceMemory(bool copy_to_host=true)
Definition tmop.cpp:3727
void AddMultPA_2D(const Vector &, Vector &) const
bool GetFDFlag() const
Definition tmop.hpp:2536
virtual real_t GetDerefinementElementEnergy(const FiniteElement &el, ElementTransformation &T, const Vector &elfun)
Definition tmop.cpp:4436
TMOP_QualityMetric & GetAMRQualityMetric()
Definition tmop.hpp:2494
void AssembleElementGrad(const FiniteElement &el, ElementTransformation &T, const Vector &d_el, DenseMatrix &elmat) override
Second derivative of GetElementEnergy() w.r.t. each local H1 DOF.
Definition tmop.cpp:4504
void AssembleElementVectorFD(const FiniteElement &el, ElementTransformation &T, const Vector &d_el, Vector &elvect)
Definition tmop.cpp:5139
void AddMultPA_3D(const Vector &, Vector &) const
DenseMatrix PMatI
Definition tmop.hpp:2089
const IntegrationRule & GradientIntegrationRule(const FiniteElement &el) const
Definition tmop.hpp:2213
void EnableAdaptiveLimiting(const GridFunction &z0, Coefficient &coeff, AdaptivityEvaluator &ae)
Restriction of the node positions to certain regions.
Definition tmop.cpp:3798
void AssembleGradPA_C0_3D(const Vector &) const
void AssembleGradPA_2D(const Vector &) const
struct mfem::TMOP_Integrator::@26 PA
AdaptivityEvaluator * surf_fit_eval_hess
Definition tmop.hpp:2059
void SetLimitingNodes(const GridFunction &n0)
Update the original/reference nodes used for limiting.
Definition tmop.hpp:2458
Coefficient * metric_coeff
Definition tmop.hpp:2023
void SetIntegrationRules(IntegrationRules &irules, int order)
Prescribe a set of integration rules; relevant for mixed meshes.
Definition tmop.hpp:2300
void GetLocalNormalizationEnergiesPA_3D(const Vector &x, real_t &met_energy, real_t &lim_energy) const
const GridFunction * x_0
Definition tmop.hpp:2010
void AddMultPA(const Vector &, Vector &) const override
Method for partially assembled action.
Definition pa.cpp:353
GridFunction * adapt_lim_gf
Definition tmop.hpp:2044
real_t GetFDDerivative(const FiniteElement &el, ElementTransformation &T, Vector &d_el, const int nodenum, const int idir, const real_t baseenergy, bool update_stored)
Definition tmop.cpp:5117
DiscreteAdaptTC * GetDiscreteAdaptTC() const
Definition tmop.hpp:2515
void UpdateAfterMeshTopologyChange()
Definition tmop.cpp:4145
void ComputeNormalizationEnergies(const GridFunction &x, real_t &metric_energy, real_t &lim_energy)
Definition tmop.cpp:5380
void AddMultGradPA_C0_2D(const Vector &, Vector &) const
IntegrationRules * IntegRules
Definition tmop.hpp:2018
void RemapSurfaceFittingLevelSetAtNodes(const Vector &new_x, int new_x_ordering)
Definition tmop.cpp:5504
void AssembleDiagonalPA_3D(Vector &) const
DenseMatrix Jrt
Definition tmop.hpp:2089
void UpdateCoefficientsPA(const Vector &d_loc)
Definition pa.cpp:187
void AssembleGradPA_3D(const Vector &) const
Coefficient * surf_fit_coeff
Definition tmop.hpp:2050
Array< int > surf_fit_dof_count
Definition tmop.hpp:2060
DenseMatrix Jpt
Definition tmop.hpp:2089
void GetLocalStateEnergyPA_3D(const Vector &, real_t &energy) const
void AssembleElementGradExact(const FiniteElement &el, ElementTransformation &T, const Vector &d_el, DenseMatrix &elmat)
Definition tmop.cpp:4695
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:3833
void AssembleGradPA(const Vector &, const FiniteElementSpace &) override
Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at th...
Definition pa.cpp:20
void AssembleElemGradSurfFit(const FiniteElement &el_x, IsoparametricTransformation &Tpr, DenseMatrix &mat)
Definition tmop.cpp:5002
void AssembleElementVector(const FiniteElement &el, ElementTransformation &T, const Vector &d_el, Vector &elvect) override
First defivative of GetElementEnergy() w.r.t. each local H1 DOF.
Definition tmop.cpp:4490
TMOP_Integrator(TMOP_QualityMetric *m, TargetConstructor *tc)
Definition tmop.hpp:2289
void AddMultPA_C0_2D(const Vector &, Vector &) const
void UpdateSurfaceFittingWeight(real_t factor)
Update the surface fitting weight as surf_fit_coeff *= factor;.
Definition tmop.cpp:5326
Array< Vector * > ElemDer
Definition tmop.hpp:2076
void IntegrateOverTarget(bool integ_over_target_)
Definition tmop.hpp:2316
TMOP_Integrator(TMOP_QualityMetric *m, TargetConstructor *tc, TMOP_QualityMetric *hm)
Definition tmop.hpp:2272
const Array< bool > * surf_fit_marker
Definition tmop.hpp:2049
real_t GetLocalStateEnergyPA(const Vector &) const override
Compute the local (to the MPI rank) energy with partial assembly.
Definition pa.cpp:411
const IntegrationRule & ActionIntegrationRule(const FiniteElement &el) const
Definition tmop.hpp:2208
void AddMultGradPA_C0_3D(const Vector &, Vector &) const
AdaptivityEvaluator * adapt_lim_eval
Definition tmop.hpp:2046
void ComputeFDh(const Vector &d, const FiniteElementSpace &fes)
Determines the perturbation, h, for FD-based approximation.
Definition tmop.cpp:5681
void AssembleElemVecSurfFit(const FiniteElement &el_x, IsoparametricTransformation &Tpr, DenseMatrix &mat)
Definition tmop.cpp:4921
void EnableLimiting(const GridFunction &n0, const GridFunction &dist, Coefficient &w0, TMOP_LimiterFunction *lfunc=NULL)
Limiting of the mesh displacements (general version).
Definition tmop.cpp:3773
void GetLocalNormalizationEnergiesPA_2D(const Vector &x, real_t &met_energy, real_t &lim_energy) const
Base class for limiting functions to be used in class TMOP_Integrator.
Definition tmop.hpp:1401
virtual real_t Eval(const Vector &x, const Vector &x0, real_t d) const =0
Returns the limiting function, f(x, x0, d).
virtual void Eval_d1(const Vector &x, const Vector &x0, real_t dist, Vector &d1) const =0
Returns the gradient of the limiting function f(x, x0, d) with respect to x.
virtual ~TMOP_LimiterFunction()
Virtual destructor.
Definition tmop.hpp:1417
virtual void Eval_d2(const Vector &x, const Vector &x0, real_t dist, DenseMatrix &d2) const =0
Returns the Hessian of the limiting function f(x, x0, d) with respect to x.
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.hpp:278
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.hpp:280
int Id() const override
Return the metric ID.
Definition tmop.hpp:285
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.hpp:282
2D non-barrier metric without a type.
Definition tmop.hpp:290
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:829
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:292
int Id() const override
Return the metric ID.
Definition tmop.hpp:303
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:823
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:817
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:988
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:377
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:976
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:971
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:982
int Id() const override
Return the metric ID.
Definition tmop.hpp:391
AD1Type EvalW_AD1(const std::vector< AD1Type > &T, const std::vector< AD1Type > &W) const override
First-derivative hook for AD-based computations.
Definition tmop.cpp:1029
AD2Type EvalW_AD2(const std::vector< AD2Type > &T, const std::vector< AD2Type > &W) const override
Second-derivative hook for AD-based computations.
Definition tmop.cpp:1035
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:998
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1004
int Id() const override
Return the metric ID.
Definition tmop.hpp:420
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:399
type EvalW_AD_impl(const std::vector< type > &T, const std::vector< type > &W) const
Definition tmop.cpp:1023
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1010
2D barrier Shape+Size (VS) metric (not polyconvex).
Definition tmop.hpp:425
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1056
int Id() const override
Return the metric ID.
Definition tmop.hpp:438
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1041
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1048
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:427
2D barrier Shape+Size (VS) metric (not polyconvex).
Definition tmop.hpp:443
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1098
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1081
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:445
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1089
int Id() const override
Return the metric ID.
Definition tmop.hpp:456
2D non-barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition tmop.hpp:461
AD2Type EvalW_AD2(const std::vector< AD2Type > &T, const std::vector< AD2Type > &W) const override
Second-derivative hook for AD-based computations.
Definition tmop.cpp:1175
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1123
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:463
int Id() const override
Return the metric ID.
Definition tmop.hpp:487
AD1Type EvalW_AD1(const std::vector< AD1Type > &T, const std::vector< AD1Type > &W) const override
First-derivative hook for AD-based computations.
Definition tmop.cpp:1169
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1114
type EvalW_AD_impl(const std::vector< type > &T, const std::vector< type > &W) const
Definition tmop.cpp:1163
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1146
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1134
2D Shifted barrier form of shape metric (mu_2).
Definition tmop.hpp:492
TMOP_Metric_022(real_t &t0)
Definition tmop.hpp:498
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1201
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:495
int Id() const override
Return the metric ID.
Definition tmop.hpp:508
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1211
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1181
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:516
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1260
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1249
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1238
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1268
int Id() const override
Return the metric ID.
Definition tmop.hpp:530
AD2Type EvalW_AD2(const std::vector< AD2Type > &T, const std::vector< AD2Type > &W) const override
Second-derivative hook for AD-based computations.
Definition tmop.cpp:1323
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1281
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:538
type EvalW_AD_impl(const std::vector< type > &T, const std::vector< type > &W) const
Definition tmop.cpp:1311
AD1Type EvalW_AD1(const std::vector< AD1Type > &T, const std::vector< AD1Type > &W) const override
First-derivative hook for AD-based computations.
Definition tmop.cpp:1317
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1297
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1289
int Id() const override
Return the metric ID.
Definition tmop.hpp:559
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1352
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1329
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:568
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1344
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1336
int Id() const override
Return the metric ID.
Definition tmop.hpp:582
2D barrier shape (S) metric (not polyconvex).
Definition tmop.hpp:587
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1377
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1385
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:589
int Id() const override
Return the metric ID.
Definition tmop.hpp:603
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1393
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1366
TMOP_QualityMetric * sz_metric
Definition tmop.hpp:612
real_t GetGamma() const
Definition tmop.hpp:623
TMOP_QualityMetric * sh_metric
Definition tmop.hpp:612
virtual ~TMOP_Metric_066()
Definition tmop.hpp:625
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:611
TMOP_Metric_066(real_t gamma)
Definition tmop.hpp:615
int Id() const override
Return the metric ID.
Definition tmop.hpp:622
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:633
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1429
int Id() const override
Return the metric ID.
Definition tmop.hpp:647
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1406
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1420
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1412
TMOP_QualityMetric * sh_metric
Definition tmop.hpp:656
real_t GetGamma() const
Definition tmop.hpp:668
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:655
virtual ~TMOP_Metric_080()
Definition tmop.hpp:670
int Id() const override
Return the metric ID.
Definition tmop.hpp:667
TMOP_Metric_080(real_t gamma)
Definition tmop.hpp:659
TMOP_QualityMetric * sz_metric
Definition tmop.hpp:656
2D barrier Shape+Orientation (OS) metric (polyconvex).
Definition tmop.hpp:675
int Id() const override
Return the metric ID.
Definition tmop.hpp:687
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1442
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1456
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1450
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.hpp:679
2D compound barrier Shape+Size (VS) metric (balanced).
Definition tmop.hpp:692
int Id() const override
Return the metric ID.
Definition tmop.hpp:707
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:694
TMOP_QualityMetric * sh_metric
Definition tmop.hpp:695
virtual ~TMOP_Metric_090()
Definition tmop.hpp:708
TMOP_QualityMetric * sz_metric
Definition tmop.hpp:695
2D compound barrier Shape+Size (VS) metric (balanced).
Definition tmop.hpp:713
TMOP_QualityMetric * sz_metric
Definition tmop.hpp:716
TMOP_QualityMetric * sh_metric
Definition tmop.hpp:716
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:715
int Id() const override
Return the metric ID.
Definition tmop.hpp:728
virtual ~TMOP_Metric_094()
Definition tmop.hpp:729
2D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition tmop.hpp:734
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.hpp:738
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1482
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1476
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1468
int Id() const override
Return the metric ID.
Definition tmop.hpp:746
2D untangling metric.
Definition tmop.hpp:751
const real_t eps
Definition tmop.hpp:753
TMOP_Metric_211(real_t epsilon=1e-4)
Definition tmop.hpp:757
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1507
int Id() const override
Return the metric ID.
Definition tmop.hpp:767
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:754
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1502
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1493
Shifted barrier form of metric 56 (area, ideal barrier metric), 2D.
Definition tmop.hpp:772
int Id() const override
Return the metric ID.
Definition tmop.hpp:789
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1523
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:775
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1535
TMOP_Metric_252(real_t &t0)
Note that t0 is stored by reference.
Definition tmop.hpp:779
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1515
3D barrier Shape (S) metric, well-posed (polyconvex & invex).
Definition tmop.hpp:794
int Id() const override
Return the metric ID.
Definition tmop.hpp:810
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:796
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1578
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1569
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1562
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1553
3D barrier Shape (S) metric, well-posed (polyconvex & invex).
Definition tmop.hpp:815
int Id() const override
Return the metric ID.
Definition tmop.hpp:831
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1613
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1640
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:817
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1622
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1632
3D barrier Shape (S) metric, well-posed (polyconvex & invex).
Definition tmop.hpp:836
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1678
int Id() const override
Return the metric ID.
Definition tmop.hpp:852
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:838
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1656
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1670
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1663
3D barrier Shape (S) metric, well-posed (polyconvex & invex).
Definition tmop.hpp:857
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1690
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1697
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:859
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1704
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1712
int Id() const override
Return the metric ID.
Definition tmop.hpp:873
3D Size (V) untangling metric.
Definition tmop.hpp:878
TMOP_Metric_311(real_t epsilon=1e-4)
Definition tmop.hpp:884
const real_t eps
Definition tmop.hpp:880
int Id() const override
Return the metric ID.
Definition tmop.hpp:894
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1733
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1724
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1741
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:881
3D Shape (S) metric, untangling version of 303.
Definition tmop.hpp:899
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1776
TMOP_Metric_313(real_t &mindet)
Definition tmop.hpp:905
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:902
int Id() const override
Return the metric ID.
Definition tmop.hpp:915
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1781
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1756
3D Size (V) metric.
Definition tmop.hpp:920
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1789
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1805
int Id() const override
Return the metric ID.
Definition tmop.hpp:933
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:922
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1797
3D Size (V) metric.
Definition tmop.hpp:938
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1824
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1832
int Id() const override
Return the metric ID.
Definition tmop.hpp:954
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1818
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:940
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1840
3D Size (V) metric.
Definition tmop.hpp:959
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1877
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1869
int Id() const override
Return the metric ID.
Definition tmop.hpp:975
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:961
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1861
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1854
3D barrier Shape+Size (VS) metric, well-posed (invex).
Definition tmop.hpp:980
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1901
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1891
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:982
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1912
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1923
int Id() const override
Return the metric ID.
Definition tmop.hpp:996
3D barrier Shape+Size (VS) metric, well-posed (invex).
Definition tmop.hpp:1001
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:1986
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:1970
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:1003
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:1949
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:1959
int Id() const override
Return the metric ID.
Definition tmop.hpp:1017
3D barrier Shape+Size (VS) metric, well-posed (invex).
Definition tmop.hpp:1022
int Id() const override
Return the metric ID.
Definition tmop.hpp:1038
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:2023
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:2037
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:2030
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:1024
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:2046
3D compound barrier Shape+Size (VS) metric (polyconvex, balanced).
Definition tmop.hpp:1043
int Id() const override
Return the metric ID.
Definition tmop.hpp:1058
TMOP_QualityMetric * sz_metric
Definition tmop.hpp:1046
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:1045
TMOP_QualityMetric * sh_metric
Definition tmop.hpp:1046
virtual ~TMOP_Metric_328()
Definition tmop.hpp:1059
3D compound barrier Shape+Size (VS) metric (polyconvex).
Definition tmop.hpp:1064
TMOP_QualityMetric * sh_metric
Definition tmop.hpp:1066
int Id() const override
Return the metric ID.
Definition tmop.hpp:1077
real_t GetGamma() const
Definition tmop.hpp:1078
TMOP_QualityMetric * sz_metric
Definition tmop.hpp:1066
TMOP_Metric_332(real_t gamma)
Definition tmop.hpp:1069
virtual ~TMOP_Metric_332()
Definition tmop.hpp:1080
3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
Definition tmop.hpp:1085
int Id() const override
Return the metric ID.
Definition tmop.hpp:1099
virtual ~TMOP_Metric_333()
Definition tmop.hpp:1100
TMOP_QualityMetric * sh_metric
Definition tmop.hpp:1088
TMOP_Metric_333(real_t gamma)
Definition tmop.hpp:1091
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:1087
TMOP_QualityMetric * sz_metric
Definition tmop.hpp:1088
3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
Definition tmop.hpp:1105
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:1107
TMOP_QualityMetric * sz_metric
Definition tmop.hpp:1108
real_t GetGamma() const
Definition tmop.hpp:1120
virtual ~TMOP_Metric_334()
Definition tmop.hpp:1122
int Id() const override
Return the metric ID.
Definition tmop.hpp:1119
TMOP_Metric_334(real_t gamma)
Definition tmop.hpp:1111
TMOP_QualityMetric * sh_metric
Definition tmop.hpp:1108
3D compound barrier Shape+Size (VS) metric (polyconvex, balanced).
Definition tmop.hpp:1127
TMOP_QualityMetric * sh_metric
Definition tmop.hpp:1130
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:1129
TMOP_QualityMetric * sz_metric
Definition tmop.hpp:1130
virtual ~TMOP_Metric_338()
Definition tmop.hpp:1143
int Id() const override
Return the metric ID.
Definition tmop.hpp:1142
3D barrier Shape+Size+Orientation (VOS) metric (polyconvex).
Definition tmop.hpp:1148
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:1150
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.hpp:1155
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:2063
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:2071
int Id() const override
Return the metric ID.
Definition tmop.hpp:1163
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:2077
3D barrier Shape+Size (VS) metric, well-posed (polyconvex).
Definition tmop.hpp:1168
int Id() const override
Return the metric ID.
Definition tmop.hpp:1182
TMOP_QualityMetric * sz_metric
Definition tmop.hpp:1171
TMOP_QualityMetric * sh_metric
Definition tmop.hpp:1171
real_t GetGamma() const
Definition tmop.hpp:1183
virtual ~TMOP_Metric_347()
Definition tmop.hpp:1185
TMOP_Metric_347(real_t gamma)
Definition tmop.hpp:1174
InvariantsEvaluator2D< real_t > ie
Definition tmop.hpp:1170
3D shifted barrier form of metric 316 (not typed).
Definition tmop.hpp:1190
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:2108
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:2088
int Id() const override
Return the metric ID.
Definition tmop.hpp:1206
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:2096
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:1193
TMOP_Metric_352(real_t &t0)
Definition tmop.hpp:1196
3D non-barrier Shape (S) metric.
Definition tmop.hpp:1211
real_t EvalWMatrixForm(const DenseMatrix &Jpt) const override
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.cpp:2131
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:2145
InvariantsEvaluator3D< real_t > ie
Definition tmop.hpp:1213
int Id() const override
Return the metric ID.
Definition tmop.hpp:1227
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:2138
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:2154
2D non-barrier Aspect ratio metric.
Definition tmop.hpp:344
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.hpp:349
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.hpp:352
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:916
3D non-barrier Aspect ratio metric.
Definition tmop.hpp:359
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.hpp:364
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:936
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.hpp:367
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.hpp:316
3D non-barrier Skew metric.
Definition tmop.hpp:329
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.hpp:334
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:876
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.hpp:337
Default limiter function in TMOP_Integrator.
Definition tmop.hpp:1422
void Eval_d1(const Vector &x, const Vector &x0, real_t dist, Vector &d1) const override
Returns the gradient of the limiting function f(x, x0, d) with respect to x.
Definition tmop.hpp:1431
void Eval_d2(const Vector &x, const Vector &x0, real_t dist, DenseMatrix &d2) const override
Returns the Hessian of the limiting function f(x, x0, d) with respect to x.
Definition tmop.hpp:1440
real_t Eval(const Vector &x, const Vector &x0, real_t dist) const override
Returns the limiting function, f(x, x0, d).
Definition tmop.hpp:1424
virtual ~TMOP_QuadraticLimiter()
Definition tmop.hpp:1448
Abstract class for local mesh quality metrics in the target-matrix optimization paradigm (TMOP) by P....
Definition tmop.hpp:28
void DefaultAssembleH(const DenseTensor &H, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const
See AssembleH(). This is a default implementation for the case when the 2nd derivatives of the metric...
Definition tmop.cpp:492
const DenseMatrix * Jtr
Definition tmop.hpp:30
virtual real_t EvalWMatrixForm(const DenseMatrix &Jpt) const
Evaluates the metric in matrix form (opposed to invariant form). Used for validating the invariant ev...
Definition tmop.hpp:57
virtual void SetTargetJacobian(const DenseMatrix &Jtr_)
Specify the reference-element -> target-element Jacobian matrix for the point of interest.
Definition tmop.hpp:53
virtual AD1Type EvalW_AD1(const std::vector< AD1Type > &T, const std::vector< AD1Type > &W) const
First-derivative hook for AD-based computations.
Definition tmop.hpp:78
void SetTransformation(ElementTransformation &)
The method HyperelasticModel::SetTransformation() is hidden for TMOP_QualityMetrics,...
Definition tmop.hpp:35
virtual AD2Type EvalW_AD2(const std::vector< AD2Type > &T, const std::vector< AD2Type > &W) const
Second-derivative hook for AD-based computations.
Definition tmop.hpp:86
virtual ~TMOP_QualityMetric()
Definition tmop.hpp:45
virtual int Id() const
Return the metric ID.
Definition tmop.hpp:108
virtual void EvalPW(const DenseMatrix &Jpt, DenseMatrix &PW) const
Definition tmop.hpp:73
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const =0
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const =0
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual real_t EvalW(const DenseMatrix &Jpt) const =0
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
real_t EvalW(const DenseMatrix &Jpt) const override
Evaluate the strain energy density function, W = W(Jpt), by using the 2D or 3D matrix invariants,...
Definition tmop.cpp:726
AD1Type EvalW_AD1(const std::vector< AD1Type > &T, const std::vector< AD1Type > &W) const override
First-derivative hook for AD-based computations.
Definition tmop.cpp:759
void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const override
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
Definition tmop.cpp:775
TMOP_WorstCaseUntangleOptimizer_Metric(TMOP_QualityMetric &tmop_metric_, int exponent_=1, real_t alpha_=1.5, real_t detT_ep_=0.0001, real_t muT_ep_=0.0001, BarrierType btype_=BarrierType::None, WorstCaseType wctype_=WorstCaseType::None)
Definition tmop.hpp:228
void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
Definition tmop.cpp:793
real_t EvalWBarrier(const DenseMatrix &Jpt) const
Definition tmop.cpp:743
AD2Type EvalW_AD2(const std::vector< AD2Type > &T, const std::vector< AD2Type > &W) const override
Second-derivative hook for AD-based computations.
Definition tmop.cpp:767
Base class representing target-matrix construction algorithms for mesh optimization via the target-ma...
Definition tmop.hpp:1585
void SetVolumeScale(real_t vol_scale)
Used by target type IDEAL_SHAPE_EQUAL_SIZE. The default volume scale is 1.
Definition tmop.hpp:1675
const TargetType target_type
Definition tmop.hpp:1611
void SetNodes(const GridFunction &n)
Set the nodes to be used in the target-matrix construction.
Definition tmop.hpp:1669
Geometry::Type current_W_type
Geometry type of current W matrix (used for cache invalidation).
Definition tmop.hpp:1617
bool Parallel() const
Definition tmop.hpp:1658
void ComputeAvgVolume() const
Definition tmop.cpp:2389
virtual void ComputeElementTargetsGradient(const IntegrationRule &ir, const Vector &elfun, IsoparametricTransformation &Tpr, DenseTensor &dJtr) const
Definition tmop.cpp:2589
virtual void ComputeElementTargets(int e_id, const FiniteElement &fe, const IntegrationRule &ir, const Vector &elfun, DenseTensor &Jtr) const
Given an element and quadrature rule, computes ref->target transformation Jacobians for each quadratu...
Definition tmop.cpp:2518
virtual ~TargetConstructor()
Definition tmop.hpp:1655
bool UsesPhysicalCoordinates() const
Return true if the methods ComputeElementTargets(), ComputeAllElementTargets(), and ComputeElementTar...
Definition tmop.hpp:1682
TargetType GetTargetType() const
Definition tmop.hpp:1677
MPI_Comm GetComm() const
Definition tmop.hpp:1659
const GridFunction * GetNodes() const
Get the nodes to be used in the target-matrix construction.
Definition tmop.hpp:1672
const GridFunction * nodes
Definition tmop.hpp:1608
virtual bool ContainsVolumeInfo() const
Checks if the target matrices contain non-trivial size specification.
Definition tmop.cpp:2504
TargetType
Target-matrix construction algorithms supported by this class.
Definition tmop.hpp:1589
TargetConstructor(TargetType ttype)
Constructor for use in serial.
Definition tmop.hpp:1641
DenseMatrix current_W
Cached copy of GeomToPerfGeomJac used on device.
Definition tmop.hpp:1615
void ComputeAllElementTargets_Fallback(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
Definition tmop.cpp:2432
TargetConstructor(TargetType ttype, MPI_Comm mpicomm)
Constructor for use in parallel.
Definition tmop.hpp:1651
bool ComputeAllElementTargets(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
Base class for vector Coefficients that optionally depend on time and space.
Vector data type.
Definition vector.hpp:82
real_t DistanceSquaredTo(const real_t *p) const
Compute the square of the Euclidean distance to another vector.
Definition vector.hpp:737
void Destroy()
Destroy a vector.
Definition vector.hpp:673
int Size() const
Returns the size of the vector.
Definition vector.hpp:234
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:584
This file contains the declaration of a dual number class.
int dim
Definition ex24.cpp:53
real_t epsilon
Definition ex25.cpp:141
future::dual< AD1Type, AD1Type > AD2Type
MFEM native AD-type for second derivatives.
Definition tmop.cpp:27
future::dual< real_t, real_t > AD1Type
MFEM native AD-type for first derivatives.
Definition tmop.cpp:25
void InterpolateTMOP_QualityMetric(TMOP_QualityMetric &metric, const TargetConstructor &tc, const Mesh &mesh, GridFunction &metric_gf)
Interpolates the metric's values at the nodes of metric_gf.
Definition tmop.cpp:6092
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition vector.cpp:570
float real_t
Definition config.hpp:46
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:492
MFEM_HOST_DEVICE Complex exp(const Complex &q)
Dual number struct (value plus gradient)
Definition dual.hpp:36
std::array< int, NCMesh::MaxFaceNodes > nodes