MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
dist_solver.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 "dist_solver.hpp"
13
14#ifdef MFEM_USE_MPI
15
16namespace mfem
17{
18
19namespace common
20{
21
22void DiffuseField(ParGridFunction &field, int smooth_steps)
23{
24 // Setup the Laplacian operator.
25 ParBilinearForm *Lap = new ParBilinearForm(field.ParFESpace());
27 Lap->Assemble();
28 Lap->Finalize();
30
31 HypreSmoother *S = new HypreSmoother(*A,0,smooth_steps);
32 S->iterative_mode = true;
33
34 Vector tmp(A->Width());
35 field.SetTrueVector();
36 Vector fieldtrue = field.GetTrueVector();
37 tmp = 0.0;
38 S->Mult(tmp, fieldtrue);
39
40 field.SetFromTrueDofs(fieldtrue);
41
42 delete A;
43 delete S;
44 delete Lap;
45}
46
48{
49 // Compute average mesh size (assumes similar cells).
50 real_t dx, loc_area = 0.0;
51 for (int i = 0; i < pmesh.GetNE(); i++)
52 {
53 loc_area += pmesh.GetElementVolume(i);
54 }
55 real_t glob_area;
56 MPI_Allreduce(&loc_area, &glob_area, 1, MPITypeMap<real_t>::mpi_type,
57 MPI_SUM, pmesh.GetComm());
58 const int glob_zones = pmesh.GetGlobalNE();
59 switch (pmesh.GetElementBaseGeometry(0))
60 {
62 dx = glob_area / glob_zones; break;
64 dx = sqrt(glob_area / glob_zones); break;
66 dx = sqrt(2.0 * glob_area / glob_zones); break;
67 case Geometry::CUBE:
68 dx = pow(glob_area / glob_zones, 1.0/3.0); break;
70 dx = pow(6.0 * glob_area / glob_zones, 1.0/3.0); break;
71 default: MFEM_ABORT("Unknown zone type!"); dx = 0.0;
72 }
73
74 return dx;
75}
76
78 ParGridFunction &dist_v)
79{
80 ParFiniteElementSpace &pfes = *dist_s.ParFESpace();
81 MFEM_VERIFY(pfes.GetOrdering()==Ordering::byNODES,
82 "Only Ordering::byNODES is supported.");
83
84 const int dim = pfes.GetMesh()->Dimension();
85 const int size = dist_s.Size();
86
87 ParGridFunction der(&pfes);
88 Vector magn(size);
89 magn = 0.0;
90 for (int d = 0; d < dim; d++)
91 {
92 dist_s.GetDerivative(1, d, der);
93 for (int i = 0; i < size; i++)
94 {
95 magn(i) += der(i) * der(i);
96 // The vector must point towards the level zero set.
97 dist_v(i + d*size) = (dist_s(i) > 0.0) ? -der(i) : der(i);
98 }
99 }
100
101 for (int i = 0; i < size; i++)
102 {
103 const real_t vec_magn = std::sqrt(magn(i) + 1e-12);
104 for (int d = 0; d < dim; d++)
105 {
106 dist_v(i + d*size) *= fabs(dist_s(i)) / vec_magn;
107 }
108 }
109}
110
112 ParGridFunction &distance)
113{
114 ParFiniteElementSpace &pfes = *distance.ParFESpace();
115 MFEM_VERIFY(pfes.GetVDim() == pfes.GetMesh()->Dimension(),
116 "This function expects a vector ParGridFunction!");
117
118 ParFiniteElementSpace pfes_s(pfes.GetParMesh(), pfes.FEColl());
119 ParGridFunction dist_s(&pfes_s);
120 ComputeScalarDistance(zero_level_set, dist_s);
121 ScalarDistToVector(dist_s, distance);
122}
123
124
126 ParGridFunction &distance)
127{
128 ParFiniteElementSpace &pfes = *distance.ParFESpace();
129
130 auto check_h1 = dynamic_cast<const H1_FECollection *>(pfes.FEColl());
131 MFEM_VERIFY(check_h1 && pfes.GetVDim() == 1,
132 "This solver supports only scalar H1 spaces.");
133
134 // Compute average mesh size (assumes similar cells).
135 ParMesh &pmesh = *pfes.GetParMesh();
136
137 // Step 0 - transform the input level set into a source-type bump.
138 ParGridFunction source(&pfes);
139 source.ProjectCoefficient(zero_level_set);
140 // Optional smoothing of the initial level set.
142 // Transform so that the peak is at 0.
143 // Assumes range [-1, 1].
144 if (transform)
145 {
146 for (int i = 0; i < source.Size(); i++)
147 {
148 const real_t x = source(i);
149 source(i) = ((x < -1.0) || (x > 1.0)) ? 0.0 : (1.0 - x) * (1.0 + x);
150 }
151 }
152
153 int amg_print_level = 0;
154
155 // Solver.
156 CGSolver cg(MPI_COMM_WORLD);
157 cg.SetRelTol(1e-12);
158 cg.SetMaxIter(100);
160 OperatorPtr A;
161 Vector B, X;
162
163 // Step 1 - diffuse.
164 ParGridFunction diffused_source(&pfes);
165 for (int i = 0; i < diffuse_iter; i++)
166 {
167 // Set up RHS.
168 ParLinearForm b(&pfes);
170 b.AddDomainIntegrator(new DomainLFIntegrator(src_coeff));
171 b.Assemble();
172
173 // Diffusion and mass terms in the LHS.
174 ParBilinearForm a_d(&pfes);
178 a_d.Assemble();
179
180 // Solve with Dirichlet BC.
181 Array<int> ess_tdof_list;
182 if (pmesh.bdr_attributes.Size())
183 {
184 Array<int> ess_bdr(pmesh.bdr_attributes.Max());
185 ess_bdr = 1;
186 pfes.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
187 }
188 ParGridFunction u_dirichlet(&pfes);
189 u_dirichlet = 0.0;
190 a_d.FormLinearSystem(ess_tdof_list, u_dirichlet, b, A, X, B);
191 auto *prec = new HypreBoomerAMG;
192 prec->SetPrintLevel(amg_print_level);
193 cg.SetPreconditioner(*prec);
194 cg.SetOperator(*A);
195 cg.Mult(B, X);
196 a_d.RecoverFEMSolution(X, b, u_dirichlet);
197 delete prec;
198
199 // Diffusion and mass terms in the LHS.
200 ParBilinearForm a_n(&pfes);
203 a_n.Assemble();
204
205 // Solve with Neumann BC.
206 ParGridFunction u_neumann(&pfes);
207 ess_tdof_list.DeleteAll();
208 a_n.FormLinearSystem(ess_tdof_list, u_neumann, b, A, X, B);
209 auto *prec2 = new HypreBoomerAMG;
210 prec2->SetPrintLevel(amg_print_level);
211 cg.SetPreconditioner(*prec2);
212 cg.SetOperator(*A);
213 cg.Mult(B, X);
214 a_n.RecoverFEMSolution(X, b, u_neumann);
215 delete prec2;
216
217 for (int ii = 0; ii < diffused_source.Size(); ii++)
218 {
219 // This assumes that the magnitudes of the two solutions are somewhat
220 // similar; otherwise one of the solutions would dominate and the BC
221 // won't look correct. To avoid this, it's good to have the source
222 // away from the boundary (i.e. have more resolution).
223 diffused_source(ii) = 0.5 * (u_neumann(ii) + u_dirichlet(ii));
224 }
225 source = diffused_source;
226 }
227
228 // Step 2 - solve for the distance using the normalized gradient.
229 {
230 // RHS - normalized gradient.
231 ParLinearForm b2(&pfes);
232 NormalizedGradCoefficient grad_u(diffused_source, pmesh.Dimension());
234 b2.Assemble();
235
236 // LHS - diffusion.
237 ParBilinearForm a2(&pfes);
239 a2.Assemble();
240
241 // No BC.
242 Array<int> no_ess_tdofs;
243
244 a2.FormLinearSystem(no_ess_tdofs, distance, b2, A, X, B);
245
246 auto *prec = new HypreBoomerAMG;
247 prec->SetPrintLevel(amg_print_level);
248 OrthoSolver ortho(pfes.GetComm());
249 ortho.SetSolver(*prec);
250 cg.SetPreconditioner(ortho);
251 cg.SetOperator(*A);
252 cg.Mult(B, X);
253 a2.RecoverFEMSolution(X, b2, distance);
254 delete prec;
255 }
256
257 // Shift the distance values to have minimum at zero.
258 real_t d_min_loc = distance.Min();
259 real_t d_min_glob;
260 MPI_Allreduce(&d_min_loc, &d_min_glob, 1, MPITypeMap<real_t>::mpi_type,
261 MPI_MIN, pfes.GetComm());
262 distance -= d_min_glob;
263
264 if (vis_glvis)
265 {
266 char vishost[] = "localhost";
267 int visport = 19916;
268
269 ParFiniteElementSpace fespace_vec(&pmesh, pfes.FEColl(),
270 pmesh.Dimension());
271 NormalizedGradCoefficient grad_u(diffused_source, pmesh.Dimension());
272 ParGridFunction x(&fespace_vec);
273 x.ProjectCoefficient(grad_u);
274
275 socketstream sol_sock_x(vishost, visport);
276 sol_sock_x << "parallel " << pfes.GetNRanks() << " "
277 << pfes.GetMyRank() << "\n";
278 sol_sock_x.precision(8);
279 sol_sock_x << "solution\n" << pmesh << x;
280 sol_sock_x << "window_geometry " << 0 << " " << 0 << " "
281 << 500 << " " << 500 << "\n"
282 << "window_title '" << "Heat Directions" << "'\n"
283 << "keys evvRj*******A\n" << std::flush;
284 }
285}
286
287real_t NormalizationDistanceSolver::NormalizationCoeff::
289{
290 T.SetIntPoint(&ip);
291 Vector u_grad;
292 u.GetGradient(T, u_grad);
293 const real_t u_value = u.GetValue(T, ip);
294
295 return u_value / sqrt(u_value * u_value + u_grad * u_grad + 1e-12);
296}
297
299 ParGridFunction &dist)
300{
301 ParFiniteElementSpace &pfes = *dist.ParFESpace();
302
303 ParGridFunction u_gf(&pfes);
304 u_gf.ProjectCoefficient(u_coeff);
305
306 NormalizationCoeff rv_coeff(u_gf);
308}
309
311 ParGridFunction &fdist)
312{
313 ParFiniteElementSpace* fesd=fdist.ParFESpace();
314
315 auto check_h1 = dynamic_cast<const H1_FECollection *>(fesd->FEColl());
316 auto check_l2 = dynamic_cast<const L2_FECollection *>(fesd->FEColl());
317 MFEM_VERIFY((check_h1 || check_l2) && fesd->GetVDim() == 1,
318 "This solver supports only scalar H1 or L2 spaces.");
319
320 ParMesh* mesh=fesd->GetParMesh();
321 const int dim=mesh->Dimension();
322
323 MPI_Comm lcomm=fesd->GetComm();
324 int myrank;
325 MPI_Comm_rank(lcomm,&myrank);
326
327 const int order = fesd->GetOrder(0);
328 H1_FECollection fecp(order, dim);
329 ParFiniteElementSpace fesp(mesh, &fecp, 1, Ordering::byVDIM);
330
331 ParGridFunction wf(&fesp);
332 wf.ProjectCoefficient(func);
333 GradientGridFunctionCoefficient gf(&wf); //gradient of wf
334
335
336 ParGridFunction xf(&fesp);
337 HypreParVector *sv = xf.GetTrueDofs();
338 *sv=1.0;
339
340 ParNonlinearForm *nf = new ParNonlinearForm(&fesp);
341
342 PUMPLaplacian* pint = new PUMPLaplacian(&func,&gf,false);
343 nf->AddDomainIntegrator(pint);
344
345 pint->SetPower(2);
346
347 //define the solvers
348 HypreBoomerAMG *prec = new HypreBoomerAMG();
349 prec->SetPrintLevel(0);
350
351 GMRESSolver *gmres;
352 gmres = new GMRESSolver(lcomm);
353 gmres->SetAbsTol(newton_abs_tol/10);
354 gmres->SetRelTol(newton_rel_tol/10);
355 gmres->SetMaxIter(100);
356 gmres->SetPrintLevel(0);
357 gmres->SetPreconditioner(*prec);
358
359 NewtonSolver ns(lcomm);
360 ns.iterative_mode = true;
361 ns.SetSolver(*gmres);
362 ns.SetOperator(*nf);
364 ns.SetRelTol(newton_rel_tol);
365 ns.SetAbsTol(newton_abs_tol);
366 ns.SetMaxIter(newton_iter);
367
368 Vector b; // RHS is zero
369 ns.Mult(b, *sv);
370
371 for (int pp=3; pp<maxp; pp++)
372 {
373 if (myrank == 0 && (print_level.summary || print_level.iterations))
374 {
375 std::cout << "pp = " << pp << std::endl;
376 }
377 pint->SetPower(pp);
378 ns.Mult(b, *sv);
379 }
380
381 xf.SetFromTrueDofs(*sv);
383 PProductCoefficient tsol(func,gfx);
384 fdist.ProjectCoefficient(tsol);
385
386 delete gmres;
387 delete prec;
388 delete nf;
389 delete sv;
390}
391
392
395 const Vector &elfun)
396{
397 real_t energy = 0.0;
398 int ndof = el.GetDof();
399 int ndim = el.GetDim();
400 const IntegrationRule *ir = NULL;
401 int order = 2 * el.GetOrder() + trans.OrderGrad(&el);
402 ir = &IntRules.Get(el.GetGeomType(), order);
403
404 Vector shapef(ndof);
405 real_t fval;
406 real_t pval;
407 DenseMatrix B(ndof, ndim);
408 Vector qval(ndim);
409
410 B=0.0;
411
412 real_t w;
413 real_t ngrad2;
414
415 for (int i = 0; i < ir->GetNPoints(); i++)
416 {
417 const IntegrationPoint &ip = ir->IntPoint(i);
418 trans.SetIntPoint(&ip);
419 w = trans.Weight();
420 w = ip.weight * w;
421
422 fval=func->Eval(trans,ip);
423
424 el.CalcPhysDShape(trans, B);
425 el.CalcPhysShape(trans,shapef);
426
427 B.MultTranspose(elfun,qval);
428
429 ngrad2=0.0;
430 for (int jj=0; jj<ndim; jj++)
431 {
432 ngrad2 = ngrad2 + qval(jj)*qval(jj);
433 }
434
435 energy = energy + w * ngrad2 * diffcoef * 0.5;
436
437 // add the external load -1 if fval > 0.0; 1 if fval < 0.0;
438 pval=shapef*elfun;
439
440 energy = energy + w * pval * pval * 0.5;
441
442 if (fval>0.0)
443 {
444 energy = energy - w*pval;
445 }
446 else if (fval<0.0)
447 {
448 energy = energy + w*pval;
449 }
450 }
451
452 return energy;
453}
454
457 const Vector &elfun,
458 Vector &elvect)
459{
460 int ndof = el.GetDof();
461 int ndim = el.GetDim();
462 const IntegrationRule *ir = NULL;
463 int order = 2 * el.GetOrder() + trans.OrderGrad(&el);
464 ir = &IntRules.Get(el.GetGeomType(), order);
465
466 elvect.SetSize(ndof);
467 elvect=0.0;
468
469 Vector shapef(ndof);
470 real_t fval;
471 real_t pval;
472
473 DenseMatrix B(ndof, ndim); //[diff_x,diff_y,diff_z]
474
475 Vector qval(ndim); //[diff_x,diff_y,diff_z,u]
476 Vector lvec(ndof); //residual at ip
477
478 B=0.0;
479 qval=0.0;
480
481 real_t w;
482
483 for (int i = 0; i < ir->GetNPoints(); i++)
484 {
485 const IntegrationPoint &ip = ir->IntPoint(i);
486 trans.SetIntPoint(&ip);
487 w = trans.Weight();
488 w = ip.weight * w;
489
490 fval=func->Eval(trans,ip);
491
492 el.CalcPhysDShape(trans, B);
493 el.CalcPhysShape(trans,shapef);
494
495 B.MultTranspose(elfun,qval);
496 B.Mult(qval,lvec);
497
498 elvect.Add(w * diffcoef,lvec);
499
500 pval=shapef*elfun;
501
502 elvect.Add(w * pval, shapef);
503
504
505 //add the load
506 //add the external load -1 if fval > 0.0; 1 if fval < 0.0;
507 pval=shapef*elfun;
508 if (fval>0.0)
509 {
510 elvect.Add( -w, shapef);
511 }
512 else if (fval<0.0)
513 {
514 elvect.Add( w, shapef);
515 }
516 }
517}
518
521 const Vector &elfun,
522 DenseMatrix &elmat)
523{
524 int ndof = el.GetDof();
525 int ndim = el.GetDim();
526 const IntegrationRule *ir = NULL;
527 int order = 2 * el.GetOrder() + trans.OrderGrad(&el);
528 ir = &IntRules.Get(el.GetGeomType(), order);
529
530 elmat.SetSize(ndof,ndof);
531 elmat=0.0;
532
533 Vector shapef(ndof);
534
535 DenseMatrix B(ndof, ndim); //[diff_x,diff_y,diff_z]
536 B = 0.0;
537
538 real_t w;
539
540 for (int i = 0; i < ir->GetNPoints(); i++)
541 {
542 const IntegrationPoint &ip = ir->IntPoint(i);
543 trans.SetIntPoint(&ip);
544 w = trans.Weight();
545 w = ip.weight * w;
546
547 el.CalcPhysDShape(trans, B);
548 el.CalcPhysShape(trans,shapef);
549
550 AddMult_a_VVt(w, shapef, elmat);
551 AddMult_a_AAt(w * diffcoef, B, elmat);
552 }
553}
554
555
558 const Vector &elfun)
559{
560 real_t energy = 0.0;
561 int ndof = el.GetDof();
562 int ndim = el.GetDim();
563 const IntegrationRule *ir = NULL;
564 int order = 2 * el.GetOrder() + trans.OrderGrad(&el);
565 ir = &IntRules.Get(el.GetGeomType(), order);
566
567 Vector shapef(ndof);
568 real_t fval;
569 real_t pval;
570 real_t tval;
571 Vector vgrad(ndim);
572 DenseMatrix dshape(ndof, ndim);
573 DenseMatrix B(ndof, ndim);
574 Vector qval(ndim);
575 Vector tmpv(ndof);
576
577 B=0.0;
578
579 real_t w;
580 real_t ngrad2;
581
582 for (int i = 0; i < ir->GetNPoints(); i++)
583 {
584 const IntegrationPoint &ip = ir->IntPoint(i);
585 trans.SetIntPoint(&ip);
586 w = trans.Weight();
587 w = ip.weight * w;
588
589 fval=func->Eval(trans,ip);
590 fgrad->Eval(vgrad,trans,ip);
591 tval=fval;
592 if (fval<0.0)
593 {
594 fval=-fval;
595 vgrad*=-1.0;
596 }
597
598 el.CalcPhysDShape(trans, dshape);
599 el.CalcPhysShape(trans,shapef);
600
601 for (int jj=0; jj<ndim; jj++)
602 {
603 dshape.GetColumn(jj,tmpv);
604 tmpv*=fval;
605 tmpv.Add(vgrad[jj],shapef);
606 B.SetCol(jj,tmpv);
607 }
608 B.MultTranspose(elfun,qval);
609
610 ngrad2=0.0;
611 for (int jj=0; jj<ndim; jj++)
612 {
613 ngrad2 = ngrad2 + qval(jj)*qval(jj);
614 }
615
616 energy = energy + w * std::pow(ngrad2+ee*ee,pp/2.0)/pp;
617
618 // add the external load -1 if fval > 0.0; 1 if fval < 0.0;
619 pval=shapef*elfun;
620 if (tval>0.0)
621 {
622 energy = energy - w * pval * tval;
623 }
624 else if (tval<0.0)
625 {
626 energy = energy + w * pval * tval;
627 }
628 }
629
630 return energy;
631}
632
635 const Vector &elfun,
636 Vector &elvect)
637{
638 int ndof = el.GetDof();
639 int ndim = el.GetDim();
640 const IntegrationRule *ir = NULL;
641 int order = 2 * el.GetOrder() + trans.OrderGrad(&el)+1;
642 ir = &IntRules.Get(el.GetGeomType(), order);
643
644 elvect.SetSize(ndof);
645 elvect=0.0;
646
647 Vector shapef(ndof);
648 real_t fval;
649 real_t tval;
650 Vector vgrad(3);
651
652 DenseMatrix dshape(ndof, ndim);
653 DenseMatrix B(ndof, ndim); // [diff_x,diff_y,diff_z]
654
655 Vector qval(ndim); // [diff_x,diff_y,diff_z,u]
656 Vector lvec(ndof); // residual at ip
657 Vector tmpv(ndof);
658
659 B=0.0;
660 qval=0.0;
661
662 real_t w;
663 real_t ngrad2;
664 real_t aa;
665
666 for (int i = 0; i < ir->GetNPoints(); i++)
667 {
668 const IntegrationPoint &ip = ir->IntPoint(i);
669 trans.SetIntPoint(&ip);
670 w = trans.Weight();
671 w = ip.weight * w;
672
673 fval=func->Eval(trans,ip);
674 fgrad->Eval(vgrad,trans,ip);
675 tval=fval;
676 if (fval<0.0)
677 {
678 fval=-fval;
679 vgrad*=-1.0;
680 }
681
682 el.CalcPhysDShape(trans, dshape);
683 el.CalcPhysShape(trans,shapef);
684
685 for (int jj=0; jj<ndim; jj++)
686 {
687 dshape.GetColumn(jj,tmpv);
688 tmpv*=fval;
689 tmpv.Add(vgrad[jj],shapef);
690 B.SetCol(jj,tmpv);
691 }
692
693 B.MultTranspose(elfun,qval);
694
695 ngrad2=0.0;
696 for (int jj=0; jj<ndim; jj++)
697 {
698 ngrad2 = ngrad2 + qval(jj)*qval(jj);
699 }
700
701 aa = ngrad2 + ee*ee;
702 aa = std::pow(aa, (pp - 2.0) / 2.0);
703 B.Mult(qval,lvec);
704 elvect.Add(w * aa,lvec);
705
706 // add the load
707 // add the external load -1 if tval > 0.0; 1 if tval < 0.0;
708 if (tval>0.0)
709 {
710 elvect.Add( -w*fval, shapef);
711 }
712 else if (tval<0.0)
713 {
714 elvect.Add( w*fval, shapef);
715 }
716 }
717}
718
721 const Vector &elfun,
722 DenseMatrix &elmat)
723{
724 int ndof = el.GetDof();
725 int ndim = el.GetDim();
726 const IntegrationRule *ir = NULL;
727 int order = 2 * el.GetOrder() + trans.OrderGrad(&el)+1;
728 ir = &IntRules.Get(el.GetGeomType(), order);
729
730 elmat.SetSize(ndof,ndof);
731 elmat=0.0;
732
733 Vector shapef(ndof);
734 real_t fval;
735 Vector vgrad(ndim);
736
737 Vector qval(ndim); // [diff_x,diff_y,diff_z,u]
738 DenseMatrix dshape(ndof, ndim);
739 DenseMatrix B(ndof, ndim); // [diff_x,diff_y,diff_z]
740 Vector lvec(ndof);
741 Vector tmpv(ndof);
742
743 B=0.0;
744
745 real_t w;
746 real_t ngrad2;
747 real_t aa, aa0, aa1;
748
749 for (int i = 0; i < ir->GetNPoints(); i++)
750 {
751 const IntegrationPoint &ip = ir->IntPoint(i);
752 trans.SetIntPoint(&ip);
753 w = trans.Weight();
754 w = ip.weight * w;
755
756 fval=func->Eval(trans,ip);
757 fgrad->Eval(vgrad,trans,ip);
758 if (fval<0.0)
759 {
760 fval=-fval;
761 vgrad*=-1.0;
762 }
763
764 el.CalcPhysDShape(trans, dshape);
765 el.CalcPhysShape(trans,shapef);
766
767 for (int jj=0; jj<ndim; jj++)
768 {
769 dshape.GetColumn(jj,tmpv);
770 tmpv*=fval;
771 tmpv.Add(vgrad[jj],shapef);
772 B.SetCol(jj,tmpv);
773 }
774
775 B.MultTranspose(elfun,qval);
776 B.Mult(qval,lvec);
777
778 ngrad2=0.0;
779 for (int jj=0; jj<ndim; jj++)
780 {
781 ngrad2 = ngrad2 + qval(jj)*qval(jj);
782 }
783
784 aa = ngrad2 + ee * ee;
785 aa1 = std::pow(aa, (pp - 2.0) / 2.0);
786 aa0 = (pp-2.0) * std::pow(aa, (pp - 4.0) / 2.0);
787
788 AddMult_a_VVt(w * aa0, lvec, elmat);
789 AddMult_a_AAt(w * aa1, B, elmat);
790 }
791}
792
793
795{
796 if (sint == nullptr)
797 {
798 sint = new ScreenedPoisson(func, rr);
799 nf->AddDomainIntegrator(sint);
800 *sv = 0.0;
801 gmres->SetOperator(nf->GetGradient(*sv));
802 }
803 else { sint->SetInput(func); }
804
805 // form RHS
806 *sv = 0.0;
807 Vector rhs(sv->Size());
808 nf->Mult(*sv, rhs);
809 // filter the input field
810 gmres->Mult(rhs, *sv);
811
812 gf.SetFromTrueDofs(*sv);
813 gf.Neg();
814
816 ffield.ProjectCoefficient(gfc);
817}
818
819} // namespace common
820
821} // namespace mfem
822
823#endif
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
int Size() const
Return the logical size of the array.
Definition array.hpp:144
void DeleteAll()
Delete the whole array.
Definition array.hpp:864
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY....
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
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
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
A coefficient that is constant across space and time.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication.
Definition densemat.cpp:220
void MultTranspose(const real_t *x, real_t *y) const
Multiply a vector with the transpose matrix.
Definition densemat.cpp:262
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:105
void SetCol(int c, const real_t *col)
void GetColumn(int c, Vector &col) const
Class for domain integrator .
Definition lininteg.hpp:146
Class for domain integration .
Definition lininteg.hpp:109
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition eltrans.hpp:93
Ordering::Type GetOrdering() const
Return the ordering method.
Definition fespace.hpp:725
int GetOrder(int i) const
Returns the polynomial degree of the i'th finite element.
Definition fespace.hpp:696
const FiniteElementCollection * FEColl() const
Definition fespace.hpp:727
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:559
int GetVDim() const
Returns vector dimension.
Definition fespace.hpp:706
Abstract class for all finite elements.
Definition fe_base.hpp:239
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition fe_base.hpp:333
void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in physical space at the poi...
Definition fe_base.cpp:192
int GetDim() const
Returns the reference space dimension for the finite element.
Definition fe_base.hpp:316
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
void CalcPhysShape(ElementTransformation &Trans, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in physical space at the point ...
Definition fe_base.cpp:182
GMRES method.
Definition solvers.hpp:547
Vector coefficient defined as the Gradient of a scalar GridFunction.
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition gridfunc.hpp:144
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition gridfunc.hpp:130
void GetGradient(ElementTransformation &tr, Vector &grad) const
Gradient of a scalar function at a quadrature point.
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
The BoomerAMG solver in hypre.
Definition hypre.hpp:1691
void SetPrintLevel(int print_level)
Definition hypre.hpp:1771
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
Wrapper for hypre's parallel vector class.
Definition hypre.hpp:206
Parallel smoothers in hypre.
Definition hypre.hpp:1020
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Relax the linear system Ax=b.
Definition hypre.cpp:3777
Class for integration point with weight.
Definition intrules.hpp:35
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
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
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
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:282
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 GetElementVolume(int i)
Definition mesh.cpp:121
long long GetGlobalNE() const
Return the total (global) number of elements.
Definition mesh.hpp:1255
Geometry::Type GetElementBaseGeometry(int i) const
Definition mesh.hpp:1385
Newton's method for solving F(x)=b for a given operator F.
Definition solvers.hpp:667
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition solvers.cpp:1810
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition solvers.hpp:714
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition solvers.cpp:1821
void AddDomainIntegrator(NonlinearFormIntegrator *nlfi)
Adds new Domain Integrator.
Pointer to an Operator of a specified type.
Definition handle.hpp:34
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
Solver wrapper which orthogonalizes the input and output vector.
Definition solvers.hpp:1218
void SetSolver(Solver &s)
Set the solver used by the OrthoSolver.
Definition solvers.cpp:3436
Class for parallel bilinear form.
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(....
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.
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Abstract parallel finite element space.
Definition pfespace.hpp:29
MPI_Comm GetComm() const
Definition pfespace.hpp:273
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
ParMesh * GetParMesh() const
Definition pfespace.hpp:277
Class for parallel grid function.
Definition pgridfunc.hpp:33
void ProjectDiscCoefficient(VectorCoefficient &coeff) override
Project a discontinuous vector coefficient as a grid function on a continuous finite element space....
real_t GetValue(int i, const IntegrationPoint &ip, int vdim=1) const override
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
void SetFromTrueDofs(const Vector &tv) override
Set the GridFunction from the given true-dof vector.
ParFiniteElementSpace * ParFESpace() const
void GetDerivative(int comp, int der_comp, ParGridFunction &der)
Parallel version of GridFunction::GetDerivative(); see its documentation.
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Class for parallel linear form.
void Assemble()
Assembles the ParLinearForm i.e. sums over all domain/bdr integrators.
Class for parallel meshes.
Definition pmesh.hpp:34
MPI_Comm GetComm() const
Definition pmesh.hpp:402
Parallel non-linear operator on the true dofs.
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition operator.hpp:686
Vector data type.
Definition vector.hpp:80
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
real_t Min() const
Returns the minimal element of the vector.
Definition vector.cpp:1212
Vector & Add(const real_t a, const Vector &Va)
(*this) += a * Va
Definition vector.cpp:247
IterativeSolver::PrintLevel print_level
virtual void ComputeScalarDistance(Coefficient &zero_level_set, ParGridFunction &distance)=0
void ScalarDistToVector(ParGridFunction &dist_s, ParGridFunction &dist_v)
virtual void ComputeVectorDistance(Coefficient &zero_level_set, ParGridFunction &distance)
void ComputeScalarDistance(Coefficient &zero_level_set, ParGridFunction &distance)
void ComputeScalarDistance(Coefficient &u_coeff, ParGridFunction &dist)
void Filter(ParGridFunction &func, ParGridFunction &ffield)
void ComputeScalarDistance(Coefficient &func, ParGridFunction &fdist)
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat) override
Assemble the local gradient matrix.
virtual real_t GetElementEnergy(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun) override
Compute the local energy.
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect) override
Perform the local action of the NonlinearFormIntegrator.
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat) override
Assemble the local gradient matrix.
virtual real_t GetElementEnergy(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun) override
Compute the local energy.
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect) override
Perform the local action of the NonlinearFormIntegrator.
int dim
Definition ex24.cpp:53
void source(const Vector &x, Vector &f)
Definition ex25.cpp:620
void trans(const Vector &u, Vector &x)
Definition ex27.cpp:412
real_t b
Definition lissajous.cpp:42
real_t AvgElementSize(ParMesh &pmesh)
void DiffuseField(ParGridFunction &field, int smooth_steps)
const int visport
void AddMult_a_VVt(const real_t a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
void AddMult_a_AAt(real_t a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.
float real_t
Definition config.hpp:43
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:486
const char vishost[]
bool iterations
Detailed information about each iteration will be reported to mfem::out.
Definition solvers.hpp:88
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.