MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
tmop_tools.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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#include "tmop_tools.hpp"
13#include "nonlinearform.hpp"
14#include "pnonlinearform.hpp"
16
17namespace mfem
18{
19
20using namespace mfem;
21
22void AdvectorCG::SetInitialField(const Vector &init_nodes,
23 const Vector &init_field)
24{
25 nodes0 = init_nodes;
26 field0 = init_field;
27}
28
30 Vector &new_field,
31 int new_nodes_ordering)
32{
34#ifdef MFEM_USE_MPI
35 if (pfes) { space = pfes; }
36#endif
37 int fes_ordering = space->GetOrdering(),
38 ncomp = space->GetVDim();
39
40 // TODO: Implement for AMR meshes.
41 const int pnt_cnt = field0.Size() / ncomp;
42
43 new_field = field0;
44 Vector new_field_temp;
45 for (int i = 0; i < ncomp; i++)
46 {
47 if (fes_ordering == Ordering::byNODES)
48 {
49 new_field_temp.MakeRef(new_field, i*pnt_cnt, pnt_cnt);
50 }
51 else
52 {
53 new_field_temp.SetSize(pnt_cnt);
54 for (int j = 0; j < pnt_cnt; j++)
55 {
56 new_field_temp(j) = new_field(i + j*ncomp);
57 }
58 }
59 ComputeAtNewPositionScalar(new_nodes, new_field_temp);
60 if (fes_ordering == Ordering::byVDIM)
61 {
62 for (int j = 0; j < pnt_cnt; j++)
63 {
64 new_field(i + j*ncomp) = new_field_temp(j);
65 }
66 }
67 }
68
69 field0 = new_field;
70 nodes0 = new_nodes;
71}
72
73void AdvectorCG::ComputeAtNewPositionScalar(const Vector &new_nodes,
74 Vector &new_field)
75{
76 Mesh *m = mesh;
77#ifdef MFEM_USE_MPI
78 if (pmesh) { m = pmesh; }
79#endif
80
81 MFEM_VERIFY(m != NULL, "No mesh has been given to the AdaptivityEvaluator.");
82
83 // This will be used to move the positions.
84 GridFunction *mesh_nodes = m->GetNodes();
85 *mesh_nodes = nodes0;
86 real_t minv = new_field.Min(), maxv = new_field.Max();
87
88 // Velocity of the positions.
89 GridFunction u(mesh_nodes->FESpace());
90 subtract(new_nodes, nodes0, u);
91
92 // Define a scalar FE space for the solution, and the advection operator.
93 TimeDependentOperator *oper = NULL;
94 FiniteElementSpace *fess = NULL;
95#ifdef MFEM_USE_MPI
96 ParFiniteElementSpace *pfess = NULL;
97#endif
98 if (fes)
99 {
100 fess = new FiniteElementSpace(fes->GetMesh(), fes->FEColl(), 1);
101 oper = new SerialAdvectorCGOper(nodes0, u, *fess, al);
102 }
103#ifdef MFEM_USE_MPI
104 else if (pfes)
105 {
106 pfess = new ParFiniteElementSpace(pfes->GetParMesh(), pfes->FEColl(), 1);
107 oper = new ParAdvectorCGOper(nodes0, u, *pfess, al, opt_mt);
108 }
109#endif
110 MFEM_VERIFY(oper != NULL,
111 "No FE space has been given to the AdaptivityEvaluator.");
112 ode_solver.Init(*oper);
113
114 // Compute some time step [mesh_size / speed].
115 real_t h_min = std::numeric_limits<real_t>::infinity();
116 for (int i = 0; i < m->GetNE(); i++)
117 {
118 h_min = std::min(h_min, m->GetElementSize(i));
119 }
120 real_t v_max = 0.0;
121 const int s = new_field.Size();
122
123 u.HostReadWrite();
124 for (int i = 0; i < s; i++)
125 {
126 real_t vel = 0.;
127 for (int j = 0; j < m->Dimension(); j++)
128 {
129 vel += u(i+j*s)*u(i+j*s);
130 }
131 v_max = std::max(v_max, vel);
132 }
133
134#ifdef MFEM_USE_MPI
135 if (pfes)
136 {
137 real_t v_loc = v_max, h_loc = h_min;
138 MPI_Allreduce(&v_loc, &v_max, 1, MPITypeMap<real_t>::mpi_type, MPI_MAX,
139 pfes->GetComm());
140 MPI_Allreduce(&h_loc, &h_min, 1, MPITypeMap<real_t>::mpi_type, MPI_MIN,
141 pfes->GetComm());
142 }
143#endif
144
145 if (v_max == 0.0) // No need to change the field.
146 {
147 delete oper;
148 delete fess;
149#ifdef MFEM_USE_MPI
150 delete pfess;
151#endif
152 return;
153 }
154
155 v_max = std::sqrt(v_max);
156 real_t dt = dt_scale * h_min / v_max;
157
158 real_t t = 0.0;
159 bool last_step = false;
160 while (!last_step)
161 {
162 if (t + dt >= 1.0)
163 {
164 dt = 1.0 - t;
165 last_step = true;
166 }
167 ode_solver.Step(new_field, t, dt);
168 }
169
170 real_t glob_minv = minv,
171 glob_maxv = maxv;
172#ifdef MFEM_USE_MPI
173 if (pfes)
174 {
175 MPI_Allreduce(&minv, &glob_minv, 1, MPITypeMap<real_t>::mpi_type, MPI_MIN,
176 pfes->GetComm());
177 MPI_Allreduce(&maxv, &glob_maxv, 1, MPITypeMap<real_t>::mpi_type, MPI_MAX,
178 pfes->GetComm());
179 }
180#endif
181
182 // Trim the overshoots and undershoots.
183 new_field.HostReadWrite();
184 for (int i = 0; i < s; i++)
185 {
186 if (new_field(i) < glob_minv) { new_field(i) = glob_minv; }
187 if (new_field(i) > glob_maxv) { new_field(i) = glob_maxv; }
188 }
189
190 delete oper;
191 delete fess;
192#ifdef MFEM_USE_MPI
193 delete pfess;
194#endif
195}
196
200 AssemblyLevel al)
201 : TimeDependentOperator(fes.GetVSize()),
202 x0(x_start), x_now(*fes.GetMesh()->GetNodes()),
203 u(vel), u_coeff(&u), M(&fes), K(&fes), al(al)
204{
206 K.AddDomainIntegrator(Kinteg);
208 K.Assemble(0);
209 K.Finalize(0);
210
211 MassIntegrator *Minteg = new MassIntegrator;
212 M.AddDomainIntegrator(Minteg);
214 M.Assemble(0);
215 M.Finalize(0);
216}
217
218void SerialAdvectorCGOper::Mult(const Vector &ind, Vector &di_dt) const
219{
220 // Move the mesh.
221 const real_t t = GetTime();
222 add(x0, t, u, x_now);
224
225 // Assemble on the new mesh.
226 K.BilinearForm::operator=(0.0);
227 K.Assemble();
228 Vector rhs(K.Size());
229 K.Mult(ind, rhs);
230 M.BilinearForm::operator=(0.0);
231 M.Assemble();
232
233 di_dt = 0.0;
234 CGSolver lin_solver;
235 Solver *prec = nullptr;
236 Array<int> ess_tdof_list;
238 {
239 prec = new OperatorJacobiSmoother(M, ess_tdof_list);
240 lin_solver.SetOperator(M);
241 }
242 else
243 {
244 prec = new DSmoother(M.SpMat());
245 lin_solver.SetOperator(M.SpMat());
246 }
247 lin_solver.SetPreconditioner(*prec);
248#ifdef MFEM_USE_SINGLE
249 const real_t rtol = 1e-4;
250#else
251 const real_t rtol = 1e-12;
252#endif
253 lin_solver.SetRelTol(rtol); lin_solver.SetAbsTol(0.0);
254 lin_solver.SetMaxIter(100);
255 lin_solver.SetPrintLevel(0);
256 lin_solver.Mult(rhs, di_dt);
257
258 delete prec;
259}
260
261#ifdef MFEM_USE_MPI
265 AssemblyLevel al,
266 MemoryType mt)
267 : TimeDependentOperator(pfes.GetVSize()),
268 x0(x_start), x_now(*pfes.GetMesh()->GetNodes()),
269 u(vel), u_coeff(&u), M(&pfes), K(&pfes), al(al)
270{
273 {
274 Kinteg->SetPAMemoryType(mt);
275 }
276 K.AddDomainIntegrator(Kinteg);
278 K.Assemble(0);
279 K.Finalize(0);
280
281 MassIntegrator *Minteg = new MassIntegrator;
283 {
284 Minteg->SetPAMemoryType(mt);
285 }
286 M.AddDomainIntegrator(Minteg);
288 M.Assemble(0);
289 M.Finalize(0);
290}
291
292void ParAdvectorCGOper::Mult(const Vector &ind, Vector &di_dt) const
293{
294 // Move the mesh.
295 const real_t t = GetTime();
296 add(x0, t, u, x_now);
298
299 // Assemble on the new mesh.
300 K.BilinearForm::operator=(0.0);
301 K.Assemble();
303 K.Mult(ind, rhs);
304 M.BilinearForm::operator=(0.0);
305 M.Assemble();
306
307 HypreParVector *RHS = rhs.ParallelAssemble();
309 X = 0.0;
310
311 OperatorHandle Mop;
312 Solver *prec = nullptr;
313 Array<int> ess_tdof_list;
315 {
316 M.FormSystemMatrix(ess_tdof_list, Mop);
317 prec = new OperatorJacobiSmoother(M, ess_tdof_list);
318 }
319 else
320 {
321 Mop.Reset(M.ParallelAssemble());
322 prec = new HypreSmoother;
323 static_cast<HypreSmoother*>(prec)->SetType(HypreSmoother::Jacobi, 1);
324 }
325
326 CGSolver lin_solver(M.ParFESpace()->GetParMesh()->GetComm());
327 lin_solver.SetPreconditioner(*prec);
328 lin_solver.SetOperator(*Mop);
329#ifdef MFEM_USE_SINGLE
330 const real_t rtol = 1e-4;
331#else
332 const real_t rtol = 1e-8;
333#endif
334 lin_solver.SetRelTol(rtol); lin_solver.SetAbsTol(0.0);
335 lin_solver.SetMaxIter(100);
336 lin_solver.SetPrintLevel(0);
337 lin_solver.Mult(*RHS, X);
338 K.ParFESpace()->GetProlongationMatrix()->Mult(X, di_dt);
339
340 delete RHS;
341 delete prec;
342}
343#endif
344
345#ifdef MFEM_USE_GSLIB
347 const Vector &init_field)
348{
349 nodes0 = init_nodes;
350 Mesh *m = mesh;
351#ifdef MFEM_USE_MPI
352 if (pmesh) { m = pmesh; }
353#endif
354 m->SetNodes(nodes0);
355
356 const real_t rel_bbox_el = 0.1;
357 const real_t newton_tol = 1.0e-12;
358 const int npts_at_once = 256;
359
360 if (finder)
361 {
362 finder->FreeData();
363 delete finder;
364 }
365
367#ifdef MFEM_USE_MPI
368 if (pfes)
369 {
370 f = pfes;
371 finder = new FindPointsGSLIB(pfes->GetComm());
372 }
373 else { finder = new FindPointsGSLIB(); }
374#else
375 finder = new FindPointsGSLIB();
376#endif
377 finder->Setup(*m, rel_bbox_el, newton_tol, npts_at_once);
378
379 field0_gf.SetSpace(f);
380 field0_gf = init_field;
381}
382
384 Vector &new_field,
385 int new_nodes_ordering)
386{
387 finder->Interpolate(new_nodes, field0_gf, new_field, new_nodes_ordering);
388}
389
390#endif
391
393 const Vector &b) const
394{
395 const FiniteElementSpace *fes = NULL;
396 real_t energy_in = 0.0;
397#ifdef MFEM_USE_MPI
398 const ParNonlinearForm *p_nlf = dynamic_cast<const ParNonlinearForm *>(oper);
399 MFEM_VERIFY(!(parallel && p_nlf == NULL), "Invalid Operator subclass.");
400 if (parallel)
401 {
402 fes = p_nlf->FESpace();
403 energy_in = p_nlf->GetEnergy(x);
404 }
405#endif
406 const bool serial = !parallel;
407 const NonlinearForm *nlf = dynamic_cast<const NonlinearForm *>(oper);
408 MFEM_VERIFY(!(serial && nlf == NULL), "Invalid Operator subclass.");
409 if (serial)
410 {
411 fes = nlf->FESpace();
412 energy_in = nlf->GetEnergy(x);
413 }
414
415 // Get the local prolongation of the solution vector.
416 Vector x_out_loc(fes->GetVSize(),
418 if (serial)
419 {
420 const SparseMatrix *cP = fes->GetConformingProlongation();
421 if (!cP) { x_out_loc = x; }
422 else { cP->Mult(x, x_out_loc); }
423 }
424#ifdef MFEM_USE_MPI
425 else
426 {
427 fes->GetProlongationMatrix()->Mult(x, x_out_loc);
428 }
429#endif
430
431 real_t scale = 1.0;
432 real_t avg_surf_fit_err, max_surf_fit_err = 0.0;
433 if (surf_fit_max_threshold > 0.0)
434 {
435 GetSurfaceFittingError(x_out_loc, avg_surf_fit_err, max_surf_fit_err);
436 if (max_surf_fit_err < surf_fit_max_threshold)
437 {
439 {
440 mfem::out << "TMOPNewtonSolver converged "
441 "based on the surface fitting error.\n";
442 }
443 scale = 0.0;
444 return scale;
445 }
446 }
448 {
450 {
451 mfem::out << "TMOPNewtonSolver converged "
452 "based on max number of times surface fitting weight can"
453 "be increased. \n";
454 }
455 scale = 0.0;
456 return scale;
457 }
458
459 // Check if the starting mesh (given by x) is inverted. Note that x hasn't
460 // been modified by the Newton update yet.
461 const real_t min_detT_in = ComputeMinDet(x_out_loc, *fes);
462 const bool untangling = (min_detT_in <= 0.0) ? true : false;
463 const real_t untangle_factor = 1.5;
464 if (untangling)
465 {
466 // Needed for the line search below. The untangling metrics see this
467 // reference to detect deteriorations.
468 MFEM_VERIFY(min_det_ptr != NULL, " Initial mesh was valid, but"
469 " intermediate mesh is invalid. Contact TMOP Developers.");
470 MFEM_VERIFY(min_detJ_threshold == 0.0,
471 "This setup is not supported. Contact TMOP Developers.");
472 *min_det_ptr = untangle_factor * min_detT_in;
473 }
474
475 const bool have_b = (b.Size() == Height());
476
477 Vector x_out(x.Size());
478 bool x_out_ok = false;
479 real_t energy_out = 0.0, min_detT_out;
480 const real_t norm_in = Norm(r);
481
482 const real_t detJ_factor = (solver_type == 1) ? 0.25 : 0.5;
484 // TODO:
485 // - Customized line search for worst-quality optimization.
486 // - What is the Newton exit criterion for worst-quality optimization?
487
488 // Perform the line search.
489 for (int i = 0; i < 12; i++)
490 {
491 // Update the mesh and get the L-vector in x_out_loc.
492 add(x, -scale, c, x_out);
493 if (serial)
494 {
495 const SparseMatrix *cP = fes->GetConformingProlongation();
496 if (!cP) { x_out_loc = x_out; }
497 else { cP->Mult(x_out, x_out_loc); }
498 }
499#ifdef MFEM_USE_MPI
500 else { fes->GetProlongationMatrix()->Mult(x_out, x_out_loc); }
501#endif
502
503 // Check the changes in detJ.
504 min_detT_out = ComputeMinDet(x_out_loc, *fes);
505 if (untangling == false && min_detT_out <= min_detJ_threshold)
506 {
507 // No untangling, and detJ got negative (or small) -- no good.
509 {
510 mfem::out << "Scale = " << scale << " Neg det(J) found.\n";
511 }
512 scale *= detJ_factor; continue;
513 }
514 if (untangling == true && min_detT_out < *min_det_ptr)
515 {
516 // Untangling, and detJ got even more negative -- no good.
518 {
519 mfem::out << "Scale = " << scale << " Neg det(J) decreased.\n";
520 }
521 scale *= detJ_factor; continue;
522 }
523
524 // Skip the energy and residual checks when we're untangling. The
525 // untangling metrics change their denominators, which can affect the
526 // energy and residual, so their increase/decrease is not relevant.
527 if (untangling) { x_out_ok = true; break; }
528
529 // Check the changes in total energy.
530 ProcessNewState(x_out);
531
532 real_t avg_fit_err, max_fit_err = 0.0;
533 if (surf_fit_max_threshold > 0.0)
534 {
535 GetSurfaceFittingError(x_out_loc, avg_fit_err, max_fit_err);
536 }
537 if (surf_fit_max_threshold > 0.0 && max_fit_err >= 1.2*max_surf_fit_err)
538 {
540 {
541 mfem::out << "Scale = " << scale << " Surf fit err increased.\n";
542 }
543 scale *= 0.5; continue;
544 }
545
546 if (serial)
547 {
548 energy_out = nlf->GetGridFunctionEnergy(x_out_loc);
549 }
550#ifdef MFEM_USE_MPI
551 else
552 {
553 energy_out = p_nlf->GetParGridFunctionEnergy(x_out_loc);
554 }
555#endif
556 if (energy_out > energy_in + 0.2*fabs(energy_in) ||
557 std::isnan(energy_out) != 0)
558 {
560 {
561 mfem::out << "Scale = " << scale << " Increasing energy: "
562 << energy_in << " --> " << energy_out << '\n';
563 }
564 scale *= 0.5; continue;
565 }
566
567 // Check the changes in the Newton residual.
568 oper->Mult(x_out, r);
569 if (have_b) { r -= b; }
570 real_t norm_out = Norm(r);
571
572 if (norm_out > 1.2*norm_in)
573 {
575 {
576 mfem::out << "Scale = " << scale << " Norm increased: "
577 << norm_in << " --> " << norm_out << '\n';
578 }
579 scale *= 0.5; continue;
580 }
581 else { x_out_ok = true; break; }
582 } // end line search
583
584 if (untangling)
585 {
586 // Update the global min detJ. Untangling metrics see this min_det_ptr.
587 if (min_detT_out > 0.0)
588 {
589 *min_det_ptr = 0.0;
592 { mfem::out << "The mesh has been untangled at the used points!\n"; }
593 }
594 else { *min_det_ptr = untangle_factor * min_detT_out; }
595 }
596
599 {
600 if (untangling)
601 {
602 mfem::out << "Min det(T) change: "
603 << min_detT_in << " -> " << min_detT_out
604 << " with " << scale << " scaling.\n";
605 }
606 else
607 {
608 mfem::out << "Energy decrease: "
609 << energy_in << " --> " << energy_out << " or "
610 << (energy_in - energy_out) / energy_in * 100.0
611 << "% with " << scale << " scaling.\n";
612 }
613 }
614
615 if (x_out_ok == false) { scale = 0.0; }
616
617 if (surf_fit_scale_factor > 0.0) { update_surf_fit_coeff = true; }
619
620 return scale;
621}
622
624{
625 const NonlinearForm *nlf = dynamic_cast<const NonlinearForm *>(oper);
626 const Array<NonlinearFormIntegrator*> &integs = *nlf->GetDNFI();
627 TMOP_Integrator *ti = NULL;
628 TMOPComboIntegrator *co = NULL;
629 for (int i = 0; i < integs.Size(); i++)
630 {
631 ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
632 if (ti)
633 {
634 ti->UpdateSurfaceFittingWeight(factor);
635 }
636 co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
637 if (co)
638 {
640 for (int j = 0; j < ati.Size(); j++)
641 {
642 ati[j]->UpdateSurfaceFittingWeight(factor);
643 }
644 }
645 }
646}
647
649{
650 const NonlinearForm *nlf = dynamic_cast<const NonlinearForm *>(oper);
651 const Array<NonlinearFormIntegrator*> &integs = *nlf->GetDNFI();
652 TMOP_Integrator *ti = NULL;
653 TMOPComboIntegrator *co = NULL;
654 weights.SetSize(0);
655 real_t weight;
656
657 for (int i = 0; i < integs.Size(); i++)
658 {
659 ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
660 if (ti)
661 {
662 weight = ti->GetSurfaceFittingWeight();
663 weights.Append(weight);
664 }
665 co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
666 if (co)
667 {
669 for (int j = 0; j < ati.Size(); j++)
670 {
671 weight = ati[j]->GetSurfaceFittingWeight();
672 weights.Append(weight);
673 }
674 }
675 }
676}
677
679 real_t &err_avg,
680 real_t &err_max) const
681{
682 const NonlinearForm *nlf = dynamic_cast<const NonlinearForm *>(oper);
683 const Array<NonlinearFormIntegrator*> &integs = *nlf->GetDNFI();
684 TMOP_Integrator *ti = NULL;
685 TMOPComboIntegrator *co = NULL;
686
687 err_avg = 0.0;
688 err_max = 0.0;
689 real_t err_avg_loc, err_max_loc;
690 for (int i = 0; i < integs.Size(); i++)
691 {
692 ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
693 if (ti)
694 {
695 if (ti->IsSurfaceFittingEnabled())
696 {
697 ti->GetSurfaceFittingErrors(x_loc, err_avg_loc, err_max_loc);
698 err_avg = std::max(err_avg_loc, err_avg);
699 err_max = std::max(err_max_loc, err_max);
700 }
701 }
702 co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
703 if (co)
704 {
706 for (int j = 0; j < ati.Size(); j++)
707 {
708 if (ati[j]->IsSurfaceFittingEnabled())
709 {
710 ati[j]->GetSurfaceFittingErrors(x_loc, err_avg_loc, err_max_loc);
711 err_avg = std::max(err_avg_loc, err_avg);
712 err_max = std::max(err_max_loc, err_max);
713 }
714 }
715 }
716 }
717}
718
720{
721 const NonlinearForm *nlf = dynamic_cast<const NonlinearForm *>(oper);
722 const Array<NonlinearFormIntegrator*> &integs = *nlf->GetDNFI();
723
724 // Reset the update flags of all TargetConstructors. This is done to avoid
725 // repeated updates of shared TargetConstructors.
726 TMOP_Integrator *ti = NULL;
727 TMOPComboIntegrator *co = NULL;
728 DiscreteAdaptTC *dtc = NULL;
729 for (int i = 0; i < integs.Size(); i++)
730 {
731 ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
732 if (ti)
733 {
734 dtc = ti->GetDiscreteAdaptTC();
735 if (dtc) { dtc->ResetUpdateFlags(); }
736 }
737 co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
738 if (co)
739 {
741 for (int j = 0; j < ati.Size(); j++)
742 {
743 dtc = ati[j]->GetDiscreteAdaptTC();
744 if (dtc) { dtc->ResetUpdateFlags(); }
745 }
746 }
747 }
748
749 Vector x_loc;
750 const FiniteElementSpace *x_fes = nullptr;
751 if (parallel)
752 {
753#ifdef MFEM_USE_MPI
754 const ParNonlinearForm *pnlf =
755 dynamic_cast<const ParNonlinearForm *>(oper);
756
757 x_fes = pnlf->ParFESpace();
758 x_loc.SetSize(x_fes->GetVSize());
759 x_fes->GetProlongationMatrix()->Mult(x, x_loc);
760#endif
761 }
762 else
763 {
764 x_fes = nlf->FESpace();
765 const Operator *P = nlf->GetProlongation();
766 if (P)
767 {
768 x_loc.SetSize(P->Height());
769 P->Mult(x,x_loc);
770 }
771 else { x_loc = x; }
772 }
773
774 for (int i = 0; i < integs.Size(); i++)
775 {
776 ti = dynamic_cast<TMOP_Integrator *>(integs[i]);
777 if (ti)
778 {
779 ti->UpdateAfterMeshPositionChange(x_loc, *x_fes);
781 {
782 ti->ComputeUntangleMetricQuantiles(x_loc, *x_fes);
783 }
784 }
785 co = dynamic_cast<TMOPComboIntegrator *>(integs[i]);
786 if (co)
787 {
789 for (int j = 0; j < ati.Size(); j++)
790 {
791 ati[j]->UpdateAfterMeshPositionChange(x_loc, *x_fes);
793 {
794 ati[j]->ComputeUntangleMetricQuantiles(x_loc, *x_fes);
795 }
796 }
797 }
798 }
799
800 // Constant coefficient associated with the surface fitting terms if
801 // adaptive surface fitting is enabled. The idea is to increase the
802 // coefficient if the surface fitting error does not sufficiently
803 // decrease between subsequent TMOPNewtonSolver iterations.
805 {
806 // Get surface fitting errors.
808 // Get array with surface fitting weights.
809 Array<real_t> weights;
811
813 {
814 mfem::out << "Avg/Max surface fitting error: " <<
815 surf_fit_err_avg << " " <<
816 surf_fit_err_max << "\n";
817 mfem::out << "Min/Max surface fitting weight: " <<
818 weights.Min() << " " << weights.Max() << "\n";
819 }
820
821 real_t change_surf_fit_err = surf_fit_err_avg_prvs-surf_fit_err_avg;
822 real_t rel_change_surf_fit_err = change_surf_fit_err/surf_fit_err_avg_prvs;
823 // Increase the surface fitting coefficient if the surface fitting error
824 // does not decrease sufficiently.
825 if (rel_change_surf_fit_err < surf_fit_rel_change_threshold)
826 {
828 adapt_inc_count += 1;
829 }
830 else
831 {
832 adapt_inc_count = 0;
833 }
835 update_surf_fit_coeff = false;
836 }
837}
838
840 const FiniteElementSpace &fes) const
841{
842 real_t min_detJ = infinity();
843 const int NE = fes.GetNE(), dim = fes.GetMesh()->Dimension();
844 Array<int> xdofs;
845 DenseMatrix Jpr(dim);
846 const bool mixed_mesh = fes.GetMesh()->GetNumGeometries(dim) > 1;
847 if (dim == 1 || mixed_mesh || UsesTensorBasis(fes) == false)
848 {
849 for (int i = 0; i < NE; i++)
850 {
851 const int dof = fes.GetFE(i)->GetDof();
852 DenseMatrix dshape(dof, dim), pos(dof, dim);
853 Vector posV(pos.Data(), dof * dim);
854
855 fes.GetElementVDofs(i, xdofs);
856 x_loc.GetSubVector(xdofs, posV);
857
858 const IntegrationRule &irule = GetIntegrationRule(*fes.GetFE(i));
859 const int nsp = irule.GetNPoints();
860 for (int j = 0; j < nsp; j++)
861 {
862 fes.GetFE(i)->CalcDShape(irule.IntPoint(j), dshape);
863 MultAtB(pos, dshape, Jpr);
864 min_detJ = std::min(min_detJ, Jpr.Det());
865 }
866 }
867 }
868 else
869 {
870 min_detJ = dim == 2 ? MinDetJpr_2D(&fes, x_loc) :
871 dim == 3 ? MinDetJpr_3D(&fes, x_loc) : 0.0;
872 }
873 real_t min_detT_all = min_detJ;
874#ifdef MFEM_USE_MPI
875 if (parallel)
876 {
877 auto p_nlf = dynamic_cast<const ParNonlinearForm *>(oper);
878 MPI_Allreduce(&min_detJ, &min_detT_all, 1, MPITypeMap<real_t>::mpi_type,
879 MPI_MIN,
880 p_nlf->ParFESpace()->GetComm());
881 }
882#endif
883 const DenseMatrix &Wideal =
885 min_detT_all /= Wideal.Det();
886
887 return min_detT_all;
888}
889
890#ifdef MFEM_USE_MPI
891// Metric values are visualized by creating an L2 finite element functions and
892// computing the metric values at the nodes.
894 const TargetConstructor &tc, ParMesh &pmesh,
895 char *title, int position)
896{
898 ParFiniteElementSpace fes(&pmesh, &fec, 1);
899 ParGridFunction metric(&fes);
900 InterpolateTMOP_QualityMetric(qm, tc, pmesh, metric);
901 socketstream sock;
902 if (pmesh.GetMyRank() == 0)
903 {
904 sock.open("localhost", 19916);
905 sock << "solution\n";
906 }
907 pmesh.PrintAsOne(sock);
908 metric.SaveAsOne(sock);
909 if (pmesh.GetMyRank() == 0)
910 {
911 sock << "window_title '"<< title << "'\n"
912 << "window_geometry "
913 << position << " " << 0 << " " << 600 << " " << 600 << "\n"
914 << "keys jRmclA\n";
915 }
916}
917#endif
918
919// Metric values are visualized by creating an L2 finite element functions and
920// computing the metric values at the nodes.
922 const TargetConstructor &tc, Mesh &mesh,
923 char *title, int position)
924{
926 FiniteElementSpace fes(&mesh, &fec, 1);
927 GridFunction metric(&fes);
928 InterpolateTMOP_QualityMetric(qm, tc, mesh, metric);
929 osockstream sock(19916, "localhost");
930 sock << "solution\n";
931 mesh.Print(sock);
932 metric.Save(sock);
933 sock.send();
934 sock << "window_title '"<< title << "'\n"
935 << "window_geometry "
936 << position << " " << 0 << " " << 600 << " " << 600 << "\n"
937 << "keys jRmclA\n";
938}
939
940}
FiniteElementSpace * fes
Definition tmop.hpp:1284
ParFiniteElementSpace * pfes
Definition tmop.hpp:1289
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field, int new_nodes_ordering=Ordering::byNODES)
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
T Min() const
Find the minimal element in the array, using the comparison operator < for class T.
Definition array.cpp:85
int Size() const
Return the logical size of the array.
Definition array.hpp:144
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:769
@ GaussLobatto
Closed type.
Definition fe_base.hpp:32
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY....
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication: .
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
int Size() const
Get the size of the BilinearForm as a square matrix.
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix: .
Conjugate gradient method.
Definition solvers.hpp:513
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition solvers.hpp:526
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:718
Data type for scaled Jacobi-type smoother of sparse matrix.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
real_t * Data() const
Returns the matrix data array.
Definition densemat.hpp:111
real_t Det() const
Definition densemat.cpp:535
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition device.hpp:274
void ResetUpdateFlags()
Used in combination with the Update methods to avoid extra computations.
Definition tmop.hpp:1614
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points....
Definition gslib.hpp:53
void Setup(Mesh &m, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
Definition gslib.cpp:109
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition gslib.cpp:856
virtual void FreeData()
Definition gslib.cpp:283
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition fespace.hpp:597
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
Definition fespace.cpp:287
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3168
Ordering::Type GetOrdering() const
Return the ordering method.
Definition fespace.hpp:725
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:740
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition fespace.cpp:1302
const FiniteElementCollection * FEColl() const
Definition fespace.hpp:727
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:559
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:713
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition fe_base.hpp:326
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:329
const DenseMatrix & GetGeomToPerfGeomJac(int GeomType) const
Definition geom.hpp:98
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
FiniteElementSpace * FESpace()
Definition gridfunc.hpp:696
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition gridfunc.cpp:195
Wrapper for hypre's parallel vector class.
Definition hypre.hpp:206
Parallel smoothers in hypre.
Definition hypre.hpp:1020
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
virtual void ComputeAtNewPosition(const Vector &new_nodes, Vector &new_field, int new_nodes_ordering=Ordering::byNODES)
virtual void SetInitialField(const Vector &init_nodes, const Vector &init_field)
PrintLevel print_options
Output behavior for the iterative solver.
Definition solvers.hpp:138
const Operator * oper
Definition solvers.hpp:121
void SetRelTol(real_t rtol)
Definition solvers.hpp:209
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.cpp:173
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition solvers.cpp:71
void SetMaxIter(int max_it)
Definition solvers.hpp:211
void SetAbsTol(real_t atol)
Definition solvers.hpp:210
real_t Norm(const Vector &x) const
Return the inner product norm of x, using the inner product defined by Dot()
Definition solvers.hpp:180
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
Mesh data type.
Definition mesh.hpp:56
void NodesUpdated()
This function should be called after the mesh node coordinates have been updated externally,...
Definition mesh.hpp:2075
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Definition mesh.hpp:2288
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
real_t GetElementSize(int i, int type=0)
Get the size of the i-th element relative to the perfect reference element.
Definition mesh.cpp:107
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:8973
void SetNodes(const Vector &node_coord)
Updates the vertex/node locations. Invokes NodesUpdated().
Definition mesh.cpp:8985
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition mesh.cpp:6961
void SetPAMemoryType(MemoryType mt)
FiniteElementSpace * FESpace()
virtual real_t GetEnergy(const Vector &x) const
Compute the energy corresponding to the state x.
Array< NonlinearFormIntegrator * > * GetDNFI()
Access all integrators added with AddDomainIntegrator().
real_t GetGridFunctionEnergy(const Vector &x) const
Compute the energy corresponding to the state x.
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
Pointer to an Operator of a specified type.
Definition handle.hpp:34
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition handle.hpp:145
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition solvers.hpp:313
Abstract operator.
Definition operator.hpp:25
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
Performs a single remap advection step in parallel.
ParAdvectorCGOper(const Vector &x_start, GridFunction &vel, ParFiniteElementSpace &pfes, AssemblyLevel al=AssemblyLevel::LEGACY, MemoryType mt=MemoryType::DEFAULT)
VectorGridFunctionCoefficient u_coeff
virtual void Mult(const Vector &ind, Vector &di_dt) const
Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x,...
const AssemblyLevel al
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void Assemble(int skip_zeros=1)
Assemble the local matrix.
ParFiniteElementSpace * ParFESpace() const
Return the parallel FE space associated with the ParBilinearForm.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
Abstract parallel finite element space.
Definition pfespace.hpp:29
MPI_Comm GetComm() const
Definition pfespace.hpp:273
const Operator * GetProlongationMatrix() const override
The returned Operator is owned by the FiniteElementSpace.
ParMesh * GetParMesh() const
Definition pfespace.hpp:277
Class for parallel grid function.
Definition pgridfunc.hpp:33
void ParallelAssemble(Vector &tv) const
Returns the vector assembled on the true dofs.
void SaveAsOne(const char *fname, int precision=16) const
Class for parallel meshes.
Definition pmesh.hpp:34
MPI_Comm GetComm() const
Definition pmesh.hpp:402
int GetMyRank() const
Definition pmesh.hpp:404
void PrintAsOne(std::ostream &out=mfem::out, const std::string &comments="") const
Definition pmesh.cpp:4968
Parallel non-linear operator on the true dofs.
real_t GetParGridFunctionEnergy(const Vector &x) const
Compute the energy corresponding to the state x.
real_t GetEnergy(const ParGridFunction &x) const
Compute the energy of a ParGridFunction.
ParFiniteElementSpace * ParFESpace() const
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:100
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:109
Performs a single remap advection step in serial.
VectorGridFunctionCoefficient u_coeff
SerialAdvectorCGOper(const Vector &x_start, GridFunction &vel, FiniteElementSpace &fes, AssemblyLevel al=AssemblyLevel::LEGACY)
const AssemblyLevel al
virtual void Mult(const Vector &ind, Vector &di_dt) const
Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x,...
Base class for solvers.
Definition operator.hpp:683
Data type sparse matrix.
Definition sparsemat.hpp:51
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
const Array< TMOP_Integrator * > & GetTMOPIntegrators() const
Definition tmop.hpp:2273
virtual void ProcessNewState(const Vector &x) const
virtual real_t ComputeScalingFactor(const Vector &x, const Vector &b) const
virtual void GetSurfaceFittingError(const Vector &x_loc, real_t &err_avg, real_t &err_max) const
real_t MinDetJpr_3D(const FiniteElementSpace *, const Vector &) const
const IntegrationRule & GetIntegrationRule(const FiniteElement &el) const
void UpdateSurfaceFittingWeight(real_t factor) const
Update surface fitting weight as surf_fit_weight *= factor.
real_t MinDetJpr_2D(const FiniteElementSpace *, const Vector &) const
real_t ComputeMinDet(const Vector &x_loc, const FiniteElementSpace &fes) const
void GetSurfaceFittingWeight(Array< real_t > &weights) const
Get the surface fitting weight for all the TMOP integrators.
A TMOP integrator class based on any given TMOP_QualityMetric and TargetConstructor.
Definition tmop.hpp:1738
void ComputeUntangleMetricQuantiles(const Vector &x, const FiniteElementSpace &fes)
Definition tmop.cpp:4683
real_t GetSurfaceFittingWeight()
Get the surface fitting weight.
Definition tmop.cpp:4237
bool IsSurfaceFittingEnabled()
Definition tmop.hpp:2160
void UpdateAfterMeshPositionChange(const Vector &x_new, const FiniteElementSpace &x_fes)
Definition tmop.cpp:4380
void GetSurfaceFittingErrors(const Vector &pos, real_t &err_avg, real_t &err_max)
Definition tmop.cpp:3079
DiscreteAdaptTC * GetDiscreteAdaptTC() const
Definition tmop.hpp:2221
void UpdateSurfaceFittingWeight(real_t factor)
Update the surface fitting weight as surf_fit_coeff *= factor;.
Definition tmop.cpp:4225
Abstract class for local mesh quality metrics in the target-matrix optimization paradigm (TMOP) by P....
Definition tmop.hpp:24
Base class representing target-matrix construction algorithms for mesh optimization via the target-ma...
Definition tmop.hpp:1334
Base abstract class for first order time dependent operators.
Definition operator.hpp:317
real_t t
Current time.
Definition operator.hpp:340
virtual real_t GetTime() const
Read the currently set time.
Definition operator.hpp:357
Vector data type.
Definition vector.hpp:80
real_t Max() const
Returns the maximal element of the vector.
Definition vector.cpp:922
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
virtual real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:494
real_t Min() const
Returns the minimal element of the vector.
Definition vector.cpp:1212
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition vector.cpp:578
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition vector.hpp:602
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
int dim
Definition ex24.cpp:53
Mesh * GetMesh(int type)
Definition ex29.cpp:218
real_t b
Definition lissajous.cpp:42
string space
real_t infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition vector.hpp:45
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
Geometry Geometries
Definition fe.cpp:49
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition vector.cpp:316
AssemblyLevel
Enumeration defining the assembly level for bilinear and nonlinear form classes derived from Operator...
void vis_tmop_metric_s(int order, TMOP_QualityMetric &qm, const TargetConstructor &tc, Mesh &mesh, char *title, int position)
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:4906
void vis_tmop_metric_p(int order, TMOP_QualityMetric &qm, const TargetConstructor &tc, ParMesh &pmesh, char *title, int position)
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Definition fespace.hpp:1345
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition vector.cpp:472
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
float real_t
Definition config.hpp:43
MemoryType
Memory types supported by MFEM.
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
void vel(const Vector &x, real_t t, Vector &u)
RefCoord t[3]
RefCoord s[3]
bool iterations
Detailed information about each iteration will be reported to mfem::out.
Definition solvers.hpp:88
bool first_and_last
Information about the first and last iteration will be printed to mfem::out.
Definition solvers.hpp:94
bool summary
A summary of the solver process will be reported after the last iteration to mfem::out.
Definition solvers.hpp:91
Helper struct to convert a C++ type to an MPI type.