Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
finite_element_state.cpp
1 // Copyright (c) Lawrence Livermore National Security, LLC and
2 // other Smith Project Developers. See the top-level LICENSE file for
3 // details.
4 //
5 // SPDX-License-Identifier: (BSD-3-Clause)
6 
8 
11 
12 namespace smith {
13 
14 void FiniteElementState::project(mfem::VectorCoefficient& coef, mfem::Array<int>& dof_list)
15 {
16  mfem::ParGridFunction& grid_function = gridFunction();
17  grid_function.ProjectCoefficient(coef, dof_list);
18  setFromGridFunction(grid_function);
19 }
20 
21 void FiniteElementState::project(mfem::Coefficient& coef, mfem::Array<int>& dof_list, std::optional<int> component)
22 {
23  mfem::ParGridFunction& grid_function = gridFunction();
24 
25  if (component) {
26  grid_function.ProjectCoefficient(coef, dof_list, *component);
27  } else {
28  grid_function.ProjectCoefficient(coef, dof_list);
29  }
30 
31  setFromGridFunction(grid_function);
32 }
33 
35 {
36  mfem::ParGridFunction& grid_function = gridFunction();
37 
38  // The generic lambda parameter, auto&&, allows the component type (mfem::Coef or mfem::VecCoef)
39  // to be deduced, and the appropriate version of ProjectCoefficient is dispatched.
40  visit(
41  [this, &grid_function](auto&& concrete_coef) {
42  grid_function.ProjectCoefficient(*concrete_coef);
43  setFromGridFunction(grid_function);
44  },
45  coef);
46 }
47 
48 void FiniteElementState::project(mfem::Coefficient& coef)
49 {
50  mfem::ParGridFunction& grid_function = gridFunction();
51  grid_function.ProjectCoefficient(coef);
52  setFromGridFunction(grid_function);
53 }
54 
55 void FiniteElementState::project(mfem::VectorCoefficient& coef)
56 {
57  mfem::ParGridFunction& grid_function = gridFunction();
58  grid_function.ProjectCoefficient(coef);
59  setFromGridFunction(grid_function);
60 }
61 
62 void FiniteElementState::projectOnBoundary(mfem::Coefficient& coef, const mfem::Array<int>& markers)
63 {
64  mfem::ParGridFunction& grid_function = gridFunction();
65  // markers should be const param in mfem, but it's not
66  grid_function.ProjectBdrCoefficient(coef, const_cast<mfem::Array<int>&>(markers));
67  setFromGridFunction(grid_function);
68 }
69 
70 void FiniteElementState::projectOnBoundary(mfem::VectorCoefficient& coef, const mfem::Array<int>& markers)
71 {
72  mfem::ParGridFunction& grid_function = gridFunction();
73  // markers should be const param in mfem, but it's not
74  grid_function.ProjectBdrCoefficient(coef, const_cast<mfem::Array<int>&>(markers));
75  setFromGridFunction(grid_function);
76 }
77 
78 void FiniteElementState::project(mfem::Coefficient& coef, const Domain& domain)
79 {
80  mfem::Array<int> uniq_dof_ids = domain.dof_list(gridFunction().ParFESpace());
81  mfem::ParGridFunction& grid_function = gridFunction();
82  grid_function.ProjectCoefficient(coef, uniq_dof_ids);
83  setFromGridFunction(grid_function);
84 }
85 
86 void FiniteElementState::project(mfem::VectorCoefficient& coef, const Domain& domain)
87 {
88  mfem::Array<int> uniq_dof_ids = domain.dof_list(gridFunction().ParFESpace());
89  mfem::ParGridFunction& grid_function = gridFunction();
90  grid_function.ProjectCoefficient(coef, uniq_dof_ids);
91  setFromGridFunction(grid_function);
92 }
93 
94 mfem::ParGridFunction& FiniteElementState::gridFunction() const
95 {
96  if (!grid_func_) {
97  grid_func_ = std::make_unique<mfem::ParGridFunction>(space_.get());
98  }
99 
101  return *grid_func_;
102 }
103 
104 double norm(const FiniteElementState& state, const double p)
105 {
106  if (state.space().GetVDim() == 1) {
107  mfem::ConstantCoefficient zero(0.0);
108  return state.gridFunction().ComputeLpError(p, zero);
109  } else {
110  mfem::Vector zero(state.space().GetVDim());
111  zero = 0.0;
112  mfem::VectorConstantCoefficient zerovec(zero);
113  return state.gridFunction().ComputeLpError(p, zerovec);
114  }
115 }
116 
117 double computeL2Error(const FiniteElementState& state, mfem::VectorCoefficient& exact_solution)
118 {
119  return state.gridFunction().ComputeL2Error(exact_solution);
120 }
121 
122 double computeL2Error(const FiniteElementState& state, mfem::Coefficient& exact_solution)
123 {
124  return state.gridFunction().ComputeL2Error(exact_solution);
125 }
126 
127 } // namespace smith
Class for encapsulating the critical MFEM components of a primal finite element field.
mfem::ParGridFunction & gridFunction() const
Construct a grid function from the finite element state true vector.
void fillGridFunction(mfem::ParGridFunction &grid_function) const
Fill a user-provided grid function based on the underlying true vector.
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 project(mfem::VectorCoefficient &coef, mfem::Array< int > &dof_list)
Project a vector coefficient onto a set of dofs.
void setFromGridFunction(const mfem::ParGridFunction &grid_function)
Initialize the true vector in the FiniteElementState based on an input grid function.
std::unique_ptr< mfem::ParFiniteElementSpace > space_
Handle to the mfem::ParFiniteElementSpace, which is owned by MFEMSidreDataCollection.
mfem::ParFiniteElementSpace & space()
Returns a non-owning reference to the internal FESpace.
many of the functions in this file amount to extracting element indices from an mesh_t like
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: smith.cpp:36
constexpr SMITH_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
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...
a class for representing a geometric region that can be used for integration
Definition: domain.hpp:33
mfem::Array< int > dof_list(const fes_t *fes) const
get mfem degree of freedom list for a given FiniteElementSpace
Definition: domain.cpp:466
A sentinel struct for eliding no-op tensor operations.
Definition: tensor.hpp:122