Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
finite_element_state.cpp
1 // Copyright (c) 2019-2024, Lawrence Livermore National Security, LLC and
2 // other Serac Project Developers. See the top-level LICENSE file for
3 // details.
4 //
5 // SPDX-License-Identifier: (BSD-3-Clause)
6 
9 
10 namespace serac {
11 
12 void FiniteElementState::project(mfem::VectorCoefficient& coef, mfem::Array<int>& dof_list)
13 {
14  mfem::ParGridFunction& grid_function = gridFunction();
15  grid_function.ProjectCoefficient(coef, dof_list);
16  setFromGridFunction(grid_function);
17 }
18 
19 void FiniteElementState::project(mfem::Coefficient& coef, mfem::Array<int>& dof_list, std::optional<int> component)
20 {
21  mfem::ParGridFunction& grid_function = gridFunction();
22 
23  if (component) {
24  grid_function.ProjectCoefficient(coef, dof_list, *component);
25  } else {
26  grid_function.ProjectCoefficient(coef, dof_list);
27  }
28 
29  setFromGridFunction(grid_function);
30 }
31 
33 {
34  mfem::ParGridFunction& grid_function = gridFunction();
35 
36  // The generic lambda parameter, auto&&, allows the component type (mfem::Coef or mfem::VecCoef)
37  // to be deduced, and the appropriate version of ProjectCoefficient is dispatched.
38  visit(
39  [this, &grid_function](auto&& concrete_coef) {
40  grid_function.ProjectCoefficient(*concrete_coef);
41  setFromGridFunction(grid_function);
42  },
43  coef);
44 }
45 
46 void FiniteElementState::project(mfem::Coefficient& coef)
47 {
48  mfem::ParGridFunction& grid_function = gridFunction();
49  grid_function.ProjectCoefficient(coef);
50  setFromGridFunction(grid_function);
51 }
52 
53 void FiniteElementState::project(mfem::VectorCoefficient& coef)
54 {
55  mfem::ParGridFunction& grid_function = gridFunction();
56  grid_function.ProjectCoefficient(coef);
57  setFromGridFunction(grid_function);
58 }
59 
60 void FiniteElementState::projectOnBoundary(mfem::Coefficient& coef, const mfem::Array<int>& markers)
61 {
62  mfem::ParGridFunction& grid_function = gridFunction();
63  // markers should be const param in mfem, but it's not
64  grid_function.ProjectBdrCoefficient(coef, const_cast<mfem::Array<int>&>(markers));
65  setFromGridFunction(grid_function);
66 }
67 
68 void FiniteElementState::projectOnBoundary(mfem::VectorCoefficient& coef, const mfem::Array<int>& markers)
69 {
70  mfem::ParGridFunction& grid_function = gridFunction();
71  // markers should be const param in mfem, but it's not
72  grid_function.ProjectBdrCoefficient(coef, const_cast<mfem::Array<int>&>(markers));
73  setFromGridFunction(grid_function);
74 }
75 
76 mfem::ParGridFunction& FiniteElementState::gridFunction() const
77 {
78  if (!grid_func_) {
79  grid_func_ = std::make_unique<mfem::ParGridFunction>(space_.get());
80  }
81 
83  return *grid_func_;
84 }
85 
86 double norm(const FiniteElementState& state, const double p)
87 {
88  if (state.space().GetVDim() == 1) {
89  mfem::ConstantCoefficient zero(0.0);
90  return state.gridFunction().ComputeLpError(p, zero);
91  } else {
92  mfem::Vector zero(state.space().GetVDim());
93  zero = 0.0;
94  mfem::VectorConstantCoefficient zerovec(zero);
95  return state.gridFunction().ComputeLpError(p, zerovec);
96  }
97 }
98 
99 double computeL2Error(const FiniteElementState& state, mfem::VectorCoefficient& exact_solution)
100 {
101  return state.gridFunction().ComputeL2Error(exact_solution);
102 }
103 
104 double computeL2Error(const FiniteElementState& state, mfem::Coefficient& exact_solution)
105 {
106  return state.gridFunction().ComputeL2Error(exact_solution);
107 }
108 
109 } // namespace serac
Class for encapsulating the critical MFEM components of a primal finite element field.
void projectOnBoundary(mfem::Coefficient &coef, const mfem::Array< int > &markers)
Project a coefficient on a specific set of marked boundaries.
std::unique_ptr< mfem::ParGridFunction > grid_func_
An optional container for a grid function (L-vector) view of the finite element state.
void setFromGridFunction(const mfem::ParGridFunction &grid_function)
Initialize the true vector in the FiniteElementState based on an input grid function.
void fillGridFunction(mfem::ParGridFunction &grid_function) const
Fill a user-provided grid function based on the underlying true vector.
mfem::ParGridFunction & gridFunction() const
Construct a grid function from the finite element state true vector.
void project(mfem::VectorCoefficient &coef, mfem::Array< int > &dof_list)
Project a vector coefficient onto a set of dofs.
mfem::ParFiniteElementSpace & space()
Returns a non-owning reference to the internal FESpace.
std::unique_ptr< mfem::ParFiniteElementSpace > space_
Handle to the mfem::ParFiniteElementSpace, which is owned by MFEMSidreDataCollection.
This file contains the declaration of structure that manages the MFEM objects that make up the state ...
This file contains the all the necessary functions and macros required for logging as well as a helpe...
Accelerator functionality.
Definition: serac.cpp:38
double computeL2Error(const FiniteElementState &state, mfem::VectorCoefficient &exact_solution)
Find the L2 norm of the error of a vector-valued finite element state with respect to an exact soluti...
constexpr SERAC_HOST_DEVICE auto norm(const isotropic_tensor< T, m, m > &I)
compute the Frobenius norm (sqrt(tr(dot(transpose(I), I)))) of an isotropic tensor
constexpr decltype(auto) visit(Visitor visitor, Variant &&v)
Applies a functor to the active variant element.
Definition: variant.hpp:365
A sentinel struct for eliding no-op tensor operations.
Definition: tensor.hpp:123