MFEM  v4.6.0
Finite element discretization library
spde_solver.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 <algorithm>
13 #include <ctime>
14 
15 #include "../../examples/ex33.hpp"
16 #include "spde_solver.hpp"
17 
18 namespace mfem
19 {
20 namespace spde
21 {
22 
23 // Helper function that determines if output should be printed to the console.
24 // The output is printed if the rank is 0 and if the print level is greater than
25 // 0. The rank is retrieved via the fespace.
26 bool PrintOutput(const ParFiniteElementSpace *fespace_ptr, int print_level)
27 {
28  return (fespace_ptr->GetMyRank() == 0 && print_level > 0);
29 }
30 
31 void Boundary::PrintInfo(std::ostream &os) const
32 {
33  os << "\n<Boundary Info>\n"
34  << " Boundary Conditions:\n";
35  for (const auto &it : boundary_attributes)
36  {
37  os << " Boundary " << it.first << ": ";
38  switch (it.second)
39  {
41  os << "Neumann";
42  break;
44  os << "Dirichlet";
45  break;
47  os << "Robin, coefficient: " << robin_coefficient;
48  break;
50  os << "Periodic";
51  break;
52  default:
53  os << "Undefined";
54  break;
55  }
56  os << "\n";
57  }
58  bool first_print_statement = true;
59  // If the map is not empty
60  if (!dirichlet_coefficients.empty())
61  {
62  os << " Inhomogeneous Dirichlet defined on ";
63  for (const auto &it : dirichlet_coefficients)
64  {
65  if (!first_print_statement)
66  {
67  os << ", ";
68  }
69  else
70  {
71  first_print_statement = false;
72  }
73  os << it.first << "(=" << it.second << ")";
74  }
75  os << "\n";
76  }
77  os << "<Boundary Info>\n\n";
78 }
79 
80 void Boundary::VerifyDefinedBoundaries(const Mesh &mesh) const
81 {
82  // Verify that all defined boundaries are actually defined on the
83  // mesh, i.e. if the keys of boundary attributes appear in the boundary
84  // attributes of the mesh.
85  mfem::out << "\n<Boundary Verify>\n";
86  const Array<int> boundary(mesh.bdr_attributes);
87  for (const auto &it : boundary_attributes)
88  {
89  if (boundary.Find(it.first) == -1)
90  {
91  MFEM_ABORT(" Boundary "
92  << it.first
93  << " is not defined on the mesh but in Boundary class."
94  << "Exiting...")
95  }
96  }
97 
98  /// Verify if all boundary attributes appear as keys in the
99  /// boundary attributes, if not let the user know that we use Neumann by
100  /// default.
101  std::vector<int> boundary_attributes_keys;
102  for (int i = 0; i < boundary.Size(); i++)
103  {
104  if (boundary_attributes.find(boundary[i]) == boundary_attributes.end())
105  {
106  boundary_attributes_keys.push_back(boundary[i]);
107  }
108  }
109  if (!boundary_attributes_keys.empty())
110  {
111  mfem::out << " Boundaries (";
112  for (const auto &it : boundary_attributes_keys)
113  {
114  mfem::out << it << ", ";
115  }
116  mfem::out << ") are defined on the mesh but not in the";
117  mfem::out << " boundary attributes (Use Neumann).";
118  }
119 
120  /// Check if any periodic boundary is registered
121  for (const auto &it : boundary_attributes)
122  {
123  if (it.second == BoundaryType::kPeriodic)
124  {
125  MFEM_ABORT(" Periodic boundaries must be defined on the mesh"
126  << ", not in Boundaries. Exiting...")
127  }
128  }
129 
130  mfem::out << "\n<Boundary Verify>\n\n";
131 }
132 
134 {
135  const ParFiniteElementSpace &fes = *solution.ParFESpace();
136  const ParMesh &pmesh = *fes.GetParMesh();
137 
138  if (PrintOutput(&fes, 1))
139  {
140  mfem::out << "<Boundary::ComputeBoundaryError>"
141  << "\n";
142  mfem::out << " GetVDim: " << fes.GetVDim() << "\n";
143  }
144 
145  double alpha{0.0};
146  double beta{1.0};
147  double gamma{0.0};
148 
149  // Index i needs to be incremented by one to map to the boundary attributes
150  // in the mesh.
151  for (int i = 0; i < pmesh.bdr_attributes.Max(); i++)
152  {
153  double error{0};
154  double avg{0};
155  Array<int> bdr(pmesh.bdr_attributes.Max());
156  bdr = 0;
157  bdr[i] = 1;
158 
159  UpdateIntegrationCoefficients(i + 1, alpha, beta, gamma);
160  avg = IntegrateBC(solution, bdr, alpha, beta, gamma, error);
161  if (PrintOutput(&fes, 1))
162  {
163  mfem::out << "->Boundary " << i + 1 << "\n";
164  mfem::out << " Alpha : " << alpha << "\n";
165  mfem::out << " Beta : " << beta << "\n";
166  mfem::out << " Gamma : " << gamma << "\n";
167  mfem::out << " Average : " << avg << "\n";
168  mfem::out << " Error : " << error << "\n\n";
169  }
170  }
171 
172  if (PrintOutput(&fes, 1))
173  {
174  mfem::out << "<Boundary::ComputeBoundaryError>" << std::endl;
175  }
176 }
177 
179  double &gamma)
180 {
181  // Check if i is a key in boundary_attributes
182  if (boundary_attributes.find(i) != boundary_attributes.end())
183  {
184  switch (boundary_attributes[i])
185  {
187  alpha = 1.0;
188  beta = 0.0;
189  gamma = 0.0;
190  break;
192  alpha = 0.0;
193  beta = 1.0;
194  if (dirichlet_coefficients.find(i) != dirichlet_coefficients.end())
195  {
196  gamma = dirichlet_coefficients[i];
197  }
198  else
199  {
200  gamma = 0.0;
201  }
202  break;
204  alpha = 1.0;
206  gamma = 0.0;
207  break;
208  default:
209  alpha = 1.0;
210  beta = 0.0;
211  gamma = 0.0;
212  break;
213  }
214  }
215  else
216  {
217  // If i is not a key in boundary_attributes, it corresponds to Neumann.
218  alpha = 1.0;
219  beta = 0.0;
220  gamma = 0.0;
221  }
222 }
223 
225  BoundaryType type)
226 {
227  boundary_attributes[boundary] = type;
228 }
229 
231  double coefficient)
232 {
234  dirichlet_coefficients[boundary] = coefficient;
235 }
236 
237 void Boundary::SetRobinCoefficient(double coefficient)
238 {
239  robin_coefficient = coefficient;
240 };
241 
242 double IntegrateBC(const ParGridFunction &x, const Array<int> &bdr,
243  double alpha, double beta, double gamma, double &glb_err)
244 {
245  double loc_vals[3];
246  double &nrm = loc_vals[0];
247  double &avg = loc_vals[1];
248  double &error = loc_vals[2];
249 
250  nrm = 0.0;
251  avg = 0.0;
252  error = 0.0;
253 
254  const bool a_is_zero = alpha == 0.0;
255  const bool b_is_zero = beta == 0.0;
256 
257  const ParFiniteElementSpace &fes = *x.ParFESpace();
258  MFEM_ASSERT(fes.GetVDim() == 1, "");
259  ParMesh &mesh = *fes.GetParMesh();
260  Vector shape;
261  Vector loc_dofs;
262  Vector w_nor;
263  DenseMatrix dshape;
264  Array<int> dof_ids;
265  for (int i = 0; i < mesh.GetNBE(); i++)
266  {
267  if (bdr[mesh.GetBdrAttribute(i) - 1] == 0)
268  {
269  continue;
270  }
271 
272  FaceElementTransformations *FTr = mesh.GetBdrFaceTransformations(i);
273  if (FTr == nullptr)
274  {
275  continue;
276  }
277 
278  const FiniteElement &fe = *fes.GetFE(FTr->Elem1No);
279  MFEM_ASSERT(fe.GetMapType() == FiniteElement::VALUE, "");
280  const int int_order = 2 * fe.GetOrder() + 3;
281  const IntegrationRule &ir = IntRules.Get(FTr->FaceGeom, int_order);
282 
283  fes.GetElementDofs(FTr->Elem1No, dof_ids);
284  x.GetSubVector(dof_ids, loc_dofs);
285  if (!a_is_zero)
286  {
287  const int sdim = FTr->Face->GetSpaceDim();
288  w_nor.SetSize(sdim);
289  dshape.SetSize(fe.GetDof(), sdim);
290  }
291  if (!b_is_zero)
292  {
293  shape.SetSize(fe.GetDof());
294  }
295  for (int j = 0; j < ir.GetNPoints(); j++)
296  {
297  const IntegrationPoint &ip = ir.IntPoint(j);
298  IntegrationPoint eip;
299  FTr->Loc1.Transform(ip, eip);
300  FTr->Face->SetIntPoint(&ip);
301  double face_weight = FTr->Face->Weight();
302  double val = 0.0;
303  if (!a_is_zero)
304  {
305  FTr->Elem1->SetIntPoint(&eip);
306  fe.CalcPhysDShape(*FTr->Elem1, dshape);
307  CalcOrtho(FTr->Face->Jacobian(), w_nor);
308  val += alpha * dshape.InnerProduct(w_nor, loc_dofs) / face_weight;
309  }
310  if (!b_is_zero)
311  {
312  fe.CalcShape(eip, shape);
313  val += beta * (shape * loc_dofs);
314  }
315 
316  // Measure the length of the boundary
317  nrm += ip.weight * face_weight;
318 
319  // Integrate alpha * n.Grad(x) + beta * x
320  avg += val * ip.weight * face_weight;
321 
322  // Integrate |alpha * n.Grad(x) + beta * x - gamma|^2
323  val -= gamma;
324  error += (val * val) * ip.weight * face_weight;
325  }
326  }
327 
328  double glb_vals[3];
329  MPI_Allreduce(loc_vals, glb_vals, 3, MPI_DOUBLE, MPI_SUM, fes.GetComm());
330 
331  double glb_nrm = glb_vals[0];
332  double glb_avg = glb_vals[1];
333  glb_err = glb_vals[2];
334 
335  // Normalize by the length of the boundary
336  if (std::abs(glb_nrm) > 0.0)
337  {
338  glb_err /= glb_nrm;
339  glb_avg /= glb_nrm;
340  }
341 
342  // Compute l2 norm of the error in the boundary condition (negative
343  // quadrature weights may produce negative 'error')
344  glb_err = (glb_err >= 0.0) ? sqrt(glb_err) : -sqrt(-glb_err);
345 
346  // Return the average value of alpha * n.Grad(x) + beta * x
347  return glb_avg;
348 }
349 
350 SPDESolver::SPDESolver(double nu, const Boundary &bc,
351  ParFiniteElementSpace *fespace, double l1, double l2,
352  double l3, double e1, double e2, double e3)
353  : k_(fespace),
354  m_(fespace),
355  fespace_ptr_(fespace),
356  bc_(bc),
357  nu_(nu),
358  l1_(l1),
359  l2_(l2),
360  l3_(l3),
361  e1_(e1),
362  e2_(e2),
363  e3_(e3)
364 {
365  if (PrintOutput(fespace_ptr_, print_level_))
366  {
367  mfem::out << "<SPDESolver> Initialize Solver .." << std::endl;
368  }
369  StopWatch sw;
370  sw.Start();
371 
372  // Resize the marker arrays for the boundary conditions
373  // Number of boundary attributes in the mesh
374  int nbc{0};
375  const auto &bdr_attributes = fespace_ptr_->GetParMesh()->bdr_attributes;
376  if (bdr_attributes.Size() > 0)
377  {
378  // Assumes a contiguous range of boundary attributes (1, 2, 3, ...)
379  nbc = bdr_attributes.Max() - bdr_attributes.Min() + 1;
380  }
381  dbc_marker_.SetSize(nbc);
382  rbc_marker_.SetSize(nbc);
383  dbc_marker_ = 0;
384  rbc_marker_ = 0;
385 
386  // Fill the marker arrays for the boundary conditions. We decrement the number
387  // it.first by one because the boundary attributes in the mesh start at 1 and
388  // the marker arrays start at 0.
389  for (const auto &it : bc_.boundary_attributes)
390  {
391  switch (it.second)
392  {
394  dbc_marker_[it.first - 1] = 1;
395  break;
397  rbc_marker_[it.first - 1] = 1;
398  break;
399  default:
400  break;
401  }
402  }
403 
404  // Handle homogeneous Dirichlet boundary conditions
405  // Note: for non zero DBC we usually need to project the boundary onto the
406  // solution. This is not necessary in this case since the boundary is
407  // homogeneous. For inhomogeneous Dirichlet we consider a lifting scheme.
408  fespace_ptr_->GetEssentialTrueDofs(dbc_marker_, ess_tdof_list_);
409 
410  // Compute the rational approximation coefficients.
411  int dim = fespace_ptr_->GetParMesh()->Dimension();
412  int space_dim = fespace_ptr_->GetParMesh()->SpaceDimension();
413  alpha_ = (nu_ + dim / 2.0) / 2.0; // fractional exponent
414  integer_order_of_exponent_ = static_cast<int>(std::floor(alpha_));
415  double exponent_to_approximate = alpha_ - integer_order_of_exponent_;
416 
417  // Compute the rational approximation coefficients.
418  ComputeRationalCoefficients(exponent_to_approximate);
419 
420  // Set the bilinear forms.
421 
422  // Assemble stiffness matrix
423  auto diffusion_tensor =
424  ConstructMatrixCoefficient(l1_, l2_, l3_, e1_, e2_, e3_, nu_, space_dim);
425  MatrixConstantCoefficient diffusion_coefficient(diffusion_tensor);
426  k_.AddDomainIntegrator(new DiffusionIntegrator(diffusion_coefficient));
427  ConstantCoefficient robin_coefficient(bc_.robin_coefficient);
428  k_.AddBoundaryIntegrator(new MassIntegrator(robin_coefficient), rbc_marker_);
429  k_.Assemble(0);
430 
431  // Assemble mass matrix
432  ConstantCoefficient one(1.0);
433  m_.AddDomainIntegrator(new MassIntegrator(one));
434  m_.Assemble(0);
435 
436  // Form matrices for the linear system
437  Array<int> empty;
438  k_.FormSystemMatrix(empty, stiffness_);
439  m_.FormSystemMatrix(empty, mass_bc_);
440 
441  // Get the restriction and prolongation matrix for transformations
442  restriction_matrix_ = fespace->GetRestrictionMatrix();
443  prolongation_matrix_ = fespace->GetProlongationMatrix();
444 
445  // Resize the vectors B and X to the appropriate size
446  if (prolongation_matrix_)
447  {
448  B_.SetSize(prolongation_matrix_->Width());
449  }
450  else
451  {
452  mfem::err << "<SPDESolver> prolongation matrix is not defined" << std::endl;
453  }
454  if (restriction_matrix_)
455  {
456  X_.SetSize(restriction_matrix_->Height());
457  }
458  else
459  {
460  mfem::err << "<SPDESolver> restriction matrix is not defined" << std::endl;
461  }
462 
463  sw.Stop();
464  if (PrintOutput(fespace_ptr_, print_level_))
465  {
466  mfem::out << "<SPDESolver::Timing> matrix assembly " << sw.RealTime()
467  << " [s]" << std::endl;
468  }
469 }
470 
472 {
473  // ------------------------------------------------------------------------
474  // Solve the PDE (A)^N g = f, i.e. compute g = (A)^{-1}^N f iteratively.
475  // ------------------------------------------------------------------------
476 
477  StopWatch sw;
478  sw.Start();
479 
480  // Zero initialize x to avoid touching uninitialized memory
481  x = 0.0;
482 
483  ParGridFunction helper_gf(fespace_ptr_);
484  helper_gf = 0.0;
485 
486  if (integer_order_of_exponent_ > 0)
487  {
488  if (PrintOutput(fespace_ptr_, print_level_))
489  {
490  mfem::out << "<SPDESolver> Solving PDE (A)^" << integer_order_of_exponent_
491  << " u = f" << std::endl;
492  }
493  ActivateRepeatedSolve();
494  Solve(b, helper_gf, 1.0, 1.0, integer_order_of_exponent_);
495  if (integer_order_)
496  {
497  // If the exponent is an integer, we can directly add the solution to the
498  // final solution and return.
499  x += helper_gf;
500  if (!bc_.dirichlet_coefficients.empty())
501  {
502  LiftSolution(x);
503  }
504  return;
505  }
506  UpdateRHS(b);
507  DeactivateRepeatedSolve();
508  }
509 
510  // ------------------------------------------------------------------------
511  // Solve the (remaining) fractional PDE by solving M integer order PDEs and
512  // adding up the solutions.
513  // ------------------------------------------------------------------------
514  if (!integer_order_)
515  {
516  // Iterate over all expansion coefficient that contribute to the
517  // solution.
518  for (int i = 0; i < coeffs_.Size(); i++)
519  {
520  if (PrintOutput(fespace_ptr_, print_level_))
521  {
522  mfem::out << "\n<SPDESolver> Solving PDE -Δ u + " << -poles_[i]
523  << " u = " << coeffs_[i] << " g " << std::endl;
524  }
525  helper_gf = 0.0;
526  Solve(b, helper_gf, 1.0 - poles_[i], coeffs_[i]);
527  x += helper_gf;
528  }
529  }
530 
531  // Apply the inhomogeneous Dirichlet boundary conditions.
532  if (!bc_.dirichlet_coefficients.empty())
533  {
534  LiftSolution(x);
535  }
536 
537  sw.Stop();
538  if (PrintOutput(fespace_ptr_, print_level_))
539  {
540  mfem::out << "<SPDESolver::Timing> all PCG solves " << sw.RealTime()
541  << " [s]" << std::endl;
542  }
543 }
544 
546 {
547  delete b_wn;
548  integ =
549  new WhiteGaussianNoiseDomainLFIntegrator(fespace_ptr_->GetComm(), seed);
550  b_wn = new ParLinearForm(fespace_ptr_);
551  b_wn->AddDomainIntegrator(integ);
552 };
553 
555 {
556  if (!b_wn)
557  {
558  MFEM_ABORT("Need to call SPDESolver::SetupRandomFieldGenerator(...) first");
559  }
560  // Create stochastic load
561  b_wn->Assemble();
562  double normalization = ConstructNormalizationCoefficient(
563  nu_, l1_, l2_, l3_, fespace_ptr_->GetParMesh()->Dimension());
564  (*b_wn) *= normalization;
565 
566  // Call back to solve to generate the random field
567  Solve(*b_wn, x);
568 };
569 
571  double l2, double l3,
572  int dim)
573 {
574  // Computation considers squaring components, computing determinant, and
575  // squaring
576  double det = 0;
577  if (dim == 1)
578  {
579  det = l1;
580  }
581  else if (dim == 2)
582  {
583  det = l1 * l2;
584  }
585  else if (dim == 3)
586  {
587  det = l1 * l2 * l3;
588  }
589  const double gamma1 = tgamma(nu + static_cast<double>(dim) / 2.0);
590  const double gamma2 = tgamma(nu);
591  return sqrt(pow(2 * M_PI, dim / 2.0) * det * gamma1 /
592  (gamma2 * pow(nu, dim / 2.0)));
593 }
594 
596  double l3, double e1,
597  double e2, double e3,
598  double nu, int dim)
599 {
600  if (dim == 3)
601  {
602  // Compute cosine and sine of the angles e1, e2, e3
603  const double c1 = cos(e1);
604  const double s1 = sin(e1);
605  const double c2 = cos(e2);
606  const double s2 = sin(e2);
607  const double c3 = cos(e3);
608  const double s3 = sin(e3);
609 
610  // Fill the rotation matrix R with the Euler angles.
611  DenseMatrix R(3, 3);
612  R(0, 0) = c1 * c3 - c2 * s1 * s3;
613  R(0, 1) = -c1 * s3 - c2 * c3 * s1;
614  R(0, 2) = s1 * s2;
615  R(1, 0) = c3 * s1 + c1 * c2 * s3;
616  R(1, 1) = c1 * c2 * c3 - s1 * s3;
617  R(1, 2) = -c1 * s2;
618  R(2, 0) = s2 * s3;
619  R(2, 1) = c3 * s2;
620  R(2, 2) = c2;
621 
622  // Multiply the rotation matrix R with the translation vector.
623  Vector l(3);
624  l(0) = std::pow(l1, 2);
625  l(1) = std::pow(l2, 2);
626  l(2) = std::pow(l3, 2);
627  l *= (1 / (2.0 * nu));
628 
629  // Compute result = R^t diag(l) R
630  DenseMatrix res(3, 3);
631  R.Transpose();
632  MultADBt(R, l, R, res);
633  return res;
634  }
635  else if (dim == 2)
636  {
637  const double c1 = cos(e1);
638  const double s1 = sin(e1);
639  DenseMatrix Rt(2, 2);
640  Rt(0, 0) = c1;
641  Rt(0, 1) = s1;
642  Rt(1, 0) = -s1;
643  Rt(1, 1) = c1;
644  Vector l(2);
645  l(0) = std::pow(l1, 2);
646  l(1) = std::pow(l2, 2);
647  l *= (1 / (2.0 * nu));
648  DenseMatrix res(2, 2);
649  MultADAt(Rt,l,res);
650  return res;
651  }
652  else
653  {
654  DenseMatrix res(1, 1);
655  res(0, 0) = std::pow(l1, 2) / (2.0 * nu);
656  return res;
657  }
658 }
659 
660 void SPDESolver::Solve(const ParLinearForm &b, ParGridFunction &x, double alpha,
661  double beta, int exponent)
662 {
663  // Form system of equations. This is less general than
664  // BilinearForm::FormLinearSystem and kind of resembles the necessary subset
665  // of instructions that we need in this case.
666  if (prolongation_matrix_)
667  {
668  prolongation_matrix_->MultTranspose(b, B_);
669  }
670  else
671  {
672  B_ = b;
673  }
674  B_ *= beta;
675 
676  if (!apply_lift_)
677  {
678  // Initialize X_ to zero. Important! Might contain nan/inf -> crash.
679  X_ = 0.0;
680  }
681  else
682  {
683  restriction_matrix_->Mult(x, X_);
684  }
685 
686  HypreParMatrix *Op =
687  Add(1.0, stiffness_, alpha, mass_bc_); // construct Operator
688  HypreParMatrix *Ae = Op->EliminateRowsCols(ess_tdof_list_);
689  Op->EliminateBC(*Ae, ess_tdof_list_, X_, B_); // only for homogeneous BC
690 
691  for (int i = 0; i < exponent; i++)
692  {
693  // Solve the linear system Op X_ = B_
694  SolveLinearSystem(Op);
695  k_.RecoverFEMSolution(X_, b, x);
696  if (repeated_solve_)
697  {
698  // Prepare for next iteration. X is a primal and B is a dual vector. B_
699  // must be updated to represent X_ in the next step. Instead of copying
700  // it, we must transform it appropriately.
701  GridFunctionCoefficient gfc(&x);
702  ParLinearForm previous_solution(fespace_ptr_);
703  previous_solution.AddDomainIntegrator(new DomainLFIntegrator(gfc));
704  previous_solution.Assemble();
705  prolongation_matrix_->MultTranspose(previous_solution, B_);
706  Op->EliminateBC(*Ae, ess_tdof_list_, X_, B_);
707  }
708  }
709  delete Ae;
710  delete Op;
711 }
712 
713 void SPDESolver::LiftSolution(ParGridFunction &x)
714 {
715  // Set lifting flag
716  apply_lift_ = true;
717 
718  // Lifting of the solution takes care of inhomogeneous boundary conditions.
719  // See doi:10.1016/j.jcp.2019.109009; section 2.6
720  if (PrintOutput(fespace_ptr_, print_level_))
721  {
722  mfem::out << "\n<SPDESolver> Applying inhomogeneous DBC" << std::endl;
723  }
724 
725  // Define temporary grid function for lifting.
726  ParGridFunction helper_gf(fespace_ptr_);
727  helper_gf = 0.0;
728 
729  // Project the boundary conditions onto the solution space.
730  for (const auto &bc : bc_.dirichlet_coefficients)
731  {
732  Array<int> marker(fespace_ptr_->GetParMesh()->bdr_attributes.Max());
733  marker = 0;
734  marker[bc.first - 1] = 1;
735  ConstantCoefficient cc(bc.second);
736  helper_gf.ProjectBdrCoefficient(cc, marker);
737  }
738 
739  // Create linear form for the right hand side.
740  ParLinearForm b(fespace_ptr_);
741  ConstantCoefficient zero(0.0);
742  b.AddDomainIntegrator(new DomainLFIntegrator(zero));
743  b.Assemble();
744 
745  // Solve the PDE for the lifting.
746  Solve(b, helper_gf, 1.0, 1.0);
747 
748  // Add the lifting to the solution.
749  x += helper_gf;
750 
751  // Reset the lifting flag.
752  apply_lift_ = false;
753 }
754 
755 void SPDESolver::UpdateRHS(ParLinearForm &b) const
756 {
757  if (!repeated_solve_)
758  {
759  // This function is only relevant for repeated solves.
760  return;
761  }
762  if (restriction_matrix_)
763  {
764  // This effectively writes the solution of the previous iteration X_ to the
765  // linear form b. Note that at the end of solve we update B_ = Mass * X_.
766  restriction_matrix_->MultTranspose(B_, b);
767  }
768  else
769  {
770  b = B_;
771  }
772 }
773 
774 void SPDESolver::SolveLinearSystem(const HypreParMatrix *Op)
775 {
776  HypreBoomerAMG prec(*Op);
777  prec.SetPrintLevel(-1);
778  CGSolver cg(fespace_ptr_->GetComm());
779  cg.SetRelTol(1e-12);
780  cg.SetMaxIter(2000);
781  cg.SetPrintLevel(3);
782  cg.SetPreconditioner(prec);
783  cg.SetOperator(*Op);
784  cg.SetPrintLevel(std::max(0, print_level_ - 1));
785  cg.Mult(B_, X_);
786 }
787 
788 void SPDESolver::ComputeRationalCoefficients(double exponent)
789 {
790  if (abs(exponent) > 1e-12)
791  {
792  if (PrintOutput(fespace_ptr_, print_level_))
793  {
794  mfem::out << "<SPDESolver> Approximating the fractional exponent "
795  << exponent << std::endl;
796  }
797  ComputePartialFractionApproximation(exponent, coeffs_, poles_);
798 
799  // If the example is build without LAPACK, the exponent
800  // might be modified by the function call above.
801  alpha_ = exponent + integer_order_of_exponent_;
802  }
803  else
804  {
805  integer_order_ = true;
806  if (PrintOutput(fespace_ptr_, print_level_))
807  {
808  mfem::out << "<SPDESolver> Treating integer order PDE." << std::endl;
809  }
810  }
811 }
812 
813 SPDESolver::~SPDESolver() { delete b_wn; }
814 
815 } // namespace spde
816 } // namespace mfem
Abstract class for all finite elements.
Definition: fe_base.hpp:233
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: operator.hpp:93
A matrix coefficient that is constant in space and time.
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:733
void ComputeBoundaryError(const ParGridFunction &solution)
void SetRobinCoefficient(double coefficient)
Set the robin coefficient for the boundary.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1031
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:253
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:96
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:1152
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:980
void Solve(ParLinearForm &b, ParGridFunction &x)
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Definition: globals.hpp:71
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:517
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:480
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition: array.cpp:68
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
static DenseMatrix ConstructMatrixCoefficient(double l1, double l2, double l3, double e1, double e2, double e3, double nu, int dim)
double RealTime()
Return the number of real seconds elapsed since the stopwatch was started.
Definition: tic_toc.cpp:429
void UpdateIntegrationCoefficients(int i, double &alpha, double &beta, double &gamma)
Abstract parallel finite element space.
Definition: pfespace.hpp:28
Vector beta
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:579
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:379
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:2673
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:480
void Stop()
Stop the stopwatch.
Definition: tic_toc.cpp:419
void PrintInfo(std::ostream &os=mfem::out) const
Print the information specifying the boundary conditions.
Definition: spde_solver.cpp:31
std::map< int, double > dirichlet_coefficients
Coefficient for inhomogeneous Dirichlet boundary conditions.
Definition: spde_solver.hpp:60
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2320
double IntegrateBC(const ParGridFunction &x, const Array< int > &bdr, double alpha, double beta, double gamma, double &glb_err)
void GenerateRandomField(ParGridFunction &x)
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:256
Class for parallel linear form.
Definition: plinearform.hpp:26
~SPDESolver()
Destructor.
IntegrationPointTransformation Loc1
Definition: eltrans.hpp:524
void Assemble()
Assembles the ParLinearForm i.e. sums over all domain/bdr integrators.
Definition: plinearform.cpp:46
double b
Definition: lissajous.cpp:42
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:119
int Find(const T &el) const
Return the first index where &#39;el&#39; is found; return -1 if not found.
Definition: array.hpp:818
Timing object.
Definition: tic_toc.hpp:35
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
Definition: densemat.cpp:335
ParMesh * GetParMesh() const
Definition: pfespace.hpp:273
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
void MultADBt(const DenseMatrix &A, const Vector &D, const DenseMatrix &B, DenseMatrix &ADBt)
ADBt = A D B^t, where D is diagonal.
Definition: densemat.cpp:2826
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void SetupRandomFieldGenerator(int seed=0)
void AddInhomogeneousDirichletBoundaryCondition(int boundary, double coefficient)
Add a inhomogeneous Dirichlet boundary condition to the boundary.
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:131
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
bool PrintOutput(const ParFiniteElementSpace *fespace_ptr, int print_level)
Definition: spde_solver.cpp:26
void Start()
Start the stopwatch. The elapsed time is not cleared.
Definition: tic_toc.cpp:408
void AddHomogeneousBoundaryCondition(int boundary, BoundaryType type)
Add a homogeneous boundary condition to the boundary.
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:702
virtual DofTransformation * GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: pfespace.cpp:470
virtual void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:908
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:275
void Transpose()
(*this) = (*this)^t
Definition: densemat.cpp:1419
virtual const FiniteElement * GetFE(int i) const
Definition: pfespace.cpp:534
SPDESolver(double nu, const Boundary &bc, ParFiniteElementSpace *fespace, double l1=0.1, double l2=0.1, double l3=0.1, double e1=0.0, double e2=0.0, double e3=0.0)
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
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:113
void ComputePartialFractionApproximation(double &alpha, Array< double > &coeffs, Array< double > &poles, double lmax=1000., double tol=1e-10, int npoints=1000, int max_order=100)
Definition: ex33.hpp:295
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:687
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition: eltrans.cpp:546
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition: mesh.hpp:1023
void MultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt = A D A^t, where D is diagonal.
Definition: densemat.cpp:2745
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
virtual int GetSpaceDim() const =0
Get the dimension of the target (physical) space.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
Class for integration point with weight.
Definition: intrules.hpp:31
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds new Boundary Integrator. Assumes ownership of bfi.
int dim
Definition: ex24.cpp:53
static double ConstructNormalizationCoefficient(double nu, double l1, double l2, double l3, int dim)
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
std::map< int, BoundaryType > boundary_attributes
Map to assign homogeneous boundary conditions to defined boundary types.
Definition: spde_solver.hpp:58
ElementTransformation * Elem1
Definition: eltrans.hpp:522
ElementTransformation * Face
Definition: eltrans.hpp:523
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
Definition: fe_base.hpp:349
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
MPI_Comm GetComm() const
Definition: pfespace.hpp:269
const double alpha
Definition: ex15.cpp:369
Vector data type.
Definition: vector.hpp:58
Class for parallel grid function.
Definition: pgridfunc.hpp:32
void VerifyDefinedBoundaries(const Mesh &mesh) const
Definition: spde_solver.cpp:80
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
Class for parallel meshes.
Definition: pmesh.hpp:32
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe_base.hpp:327