Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
functional_weak_form.hpp
Go to the documentation of this file.
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 
14 #pragma once
15 
17 #include "smith/physics/mesh.hpp"
21 
22 namespace smith {
23 
24 template <int spatial_dim, typename OutputSpace, typename inputs = Parameters<>,
25  typename input_indices = std::make_integer_sequence<int, inputs::n>>
27 
35 template <int spatial_dim, typename OutputSpace, typename... InputSpaces, int... input_indices>
36 class FunctionalWeakForm<spatial_dim, OutputSpace, Parameters<InputSpaces...>,
37  std::integer_sequence<int, input_indices...>> : public WeakForm {
38  public:
39  using SpacesT = std::vector<const mfem::ParFiniteElementSpace*>;
40 
42 
51  FunctionalWeakForm(std::string physics_name, std::shared_ptr<Mesh> mesh,
52  const mfem::ParFiniteElementSpace& output_mfem_space, const SpacesT& input_mfem_spaces)
53  : WeakForm(physics_name), mesh_(mesh)
54  {
55  std::array<const mfem::ParFiniteElementSpace*, sizeof...(InputSpaces)> trial_spaces;
56  std::array<const mfem::ParFiniteElementSpace*, sizeof...(InputSpaces) + 1> vector_residual_trial_spaces{
57  &output_mfem_space};
58 
59  SLIC_ERROR_ROOT_IF(
60  sizeof...(InputSpaces) != input_mfem_spaces.size(),
61  std::format("{} parameter spaces given in the template argument but {} input mfem spaces were supplied.",
62  sizeof...(InputSpaces), input_mfem_spaces.size()));
63 
64  // Validate output space
65  validateSpace<OutputSpace>(output_mfem_space, "output");
66 
67  // Validate input spaces
68  if constexpr (sizeof...(InputSpaces) > 0) {
69  // Validate each input space using a helper
70  validateInputSpaces<0>(input_mfem_spaces);
71 
72  for_constexpr<sizeof...(InputSpaces)>([&](auto i) { trial_spaces[i] = input_mfem_spaces[i]; });
73  for_constexpr<sizeof...(InputSpaces)>(
74  [&](auto i) { vector_residual_trial_spaces[i + 1] = input_mfem_spaces[i]; });
75  }
76 
77  const auto& shape_disp_space = mesh->shapeDisplacementSpace();
78 
79  weak_form_ = std::make_unique<ShapeAwareFunctional<ShapeDispSpace, OutputSpace(InputSpaces...)>>(
80  &shape_disp_space, &output_mfem_space, trial_spaces);
81 
82  v_dot_weak_form_residual_ =
83  std::make_unique<ShapeAwareFunctional<ShapeDispSpace, double(OutputSpace, InputSpaces...)>>(
84  &shape_disp_space, vector_residual_trial_spaces);
85  }
86 
110  template <int... active_parameters, typename BodyIntegralType>
111  void addBodyIntegral(DependsOn<active_parameters...>, std::string body_name, BodyIntegralType integrand)
112  {
113  weak_form_->AddDomainIntegral(Dimension<spatial_dim>{}, DependsOn<active_parameters...>{}, integrand,
114  mesh_->domain(body_name));
115 
116  v_dot_weak_form_residual_->AddDomainIntegral(
117  Dimension<spatial_dim>{}, DependsOn<0, 1 + active_parameters...>{},
118  [integrand](double time, auto X, auto V, auto... inputs) {
119  auto orig_tuple = integrand(time, X, inputs...);
120  return smith::inner(get<VALUE>(V), get<VALUE>(orig_tuple)) +
121  smith::inner(get<DERIVATIVE>(V), get<DERIVATIVE>(orig_tuple));
122  },
123  mesh_->domain(body_name));
124  }
125 
127  template <typename BodyForceType>
128  void addBodyIntegral(std::string body_name, BodyForceType body_integral)
129  {
130  addBodyIntegral(DependsOn<>{}, body_name, body_integral);
131  }
132 
153  template <int... active_parameters, typename BodyLoadType>
154  void addBodySource(DependsOn<active_parameters...> depends_on, std::string body_name, BodyLoadType load_function)
155  {
156  addBodyIntegral(depends_on, body_name, [load_function](double t, auto X, auto... inputs) {
157  return smith::tuple{-load_function(t, get<VALUE>(X), get<VALUE>(inputs)...), smith::zero{}};
158  });
159  }
160 
162  template <int... active_parameters, typename BodyLoadType>
163  void addBodySource(std::string body_name, BodyLoadType load_function)
164  {
165  return addBodySource(DependsOn<>{}, body_name, load_function);
166  }
167 
191  template <int... active_parameters, typename BoundaryIntegrandType>
192  void addBoundaryIntegral(DependsOn<active_parameters...>, std::string boundary_name, BoundaryIntegrandType integrand)
193  {
194  weak_form_->AddBoundaryIntegral(Dimension<spatial_dim - 1>{}, DependsOn<active_parameters...>{}, integrand,
195  mesh_->domain(boundary_name));
196 
197  v_dot_weak_form_residual_->AddBoundaryIntegral(
198  Dimension<spatial_dim - 1>{}, DependsOn<0, 1 + active_parameters...>{},
199  [integrand](double t, auto X, auto V, auto... params) {
200  auto orig_surface_flux = integrand(t, X, params...);
201  return smith::inner(get<VALUE>(V), orig_surface_flux);
202  },
203  mesh_->domain(boundary_name));
204  }
205 
207  template <typename BoundaryIntegrandType>
208  void addBoundaryIntegral(std::string boundary_name, const BoundaryIntegrandType& integrand)
209  {
210  addBoundaryIntegral(DependsOn<>{}, boundary_name, integrand);
211  }
212 
231  template <int... active_parameters, typename InteriorIntegrandType>
233  InteriorIntegrandType integrand)
234  {
235  weak_form_->AddInteriorFaceIntegral(Dimension<spatial_dim - 1>{}, DependsOn<active_parameters...>{}, integrand,
236  mesh_->domain(interior_name));
237 
238  v_dot_weak_form_residual_->AddInteriorFaceIntegral(
239  Dimension<spatial_dim - 1>{}, DependsOn<0, 1 + active_parameters...>{},
240  [integrand](double t, auto X, auto V, auto... params) {
241  auto [V1, V2] = V;
242  auto orig_surface_flux = integrand(t, X, params...);
243  auto [flux_pos, flux_neg] = orig_surface_flux;
244  return smith::inner(V1, flux_pos) + smith::inner(V2, flux_neg);
245  },
246  mesh_->domain(interior_name));
247  }
248 
250  template <typename InteriorIntegrandType>
251  void addInteriorBoundaryIntegral(std::string interior_name, const InteriorIntegrandType& integrand)
252  {
253  addInteriorBoundaryIntegral(DependsOn<>{}, interior_name, integrand);
254  }
255 
277  template <int... active_parameters, typename BoundaryFluxType>
278  void addBoundaryFlux(DependsOn<active_parameters...> depends_on, std::string boundary_name,
279  BoundaryFluxType flux_function)
280  {
281  addBoundaryIntegral(depends_on, boundary_name, [flux_function](double t, auto X, auto... inputs) {
282  auto n = cross(get<DERIVATIVE>(X));
283  return -flux_function(t, get<VALUE>(X), normalize(n), get<VALUE>(inputs)...);
284  });
285  }
286 
288  template <typename BoundaryFluxType>
289  void addBoundaryFlux(std::string boundary_name, const BoundaryFluxType& integrand)
290  {
291  addBoundaryFlux(DependsOn<>{}, boundary_name, integrand);
292  }
293 
295  mfem::Vector residual(TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector<ConstFieldPtr>& fields,
296  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& quad_fields = {}) const override
297  {
298  validateFields(fields, "residual");
299  dt_ = time_info.dt();
300  cycle_ = time_info.cycle();
301  auto ret = (*weak_form_)(time_info.time(), *shape_disp, *fields[input_indices]...);
302  return ret;
303  }
304 
306  std::unique_ptr<mfem::HypreParMatrix> jacobian(
307  TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector<ConstFieldPtr>& fields,
308  const std::vector<double>& jacobian_weights,
309  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& quad_fields = {}) const override
310  {
311  validateFields(fields, "jacobian");
312  dt_ = time_info.dt();
313  cycle_ = time_info.cycle();
314 
315  std::unique_ptr<mfem::HypreParMatrix> J;
316 
317  auto addToJ = [&J](double factor, std::unique_ptr<mfem::HypreParMatrix> jac_contrib) {
318  if (J) {
319  SLIC_ERROR_IF(J->N() != jac_contrib->N(),
320  "Multiple nonzero jacobian weights are being used on inconsistently sized input arguments.");
321  SLIC_ERROR_IF(J->M() != jac_contrib->M(),
322  "Multiple nonzero jacobian weights are being used on inconsistently sized input arguments.");
323  J->Add(factor, *jac_contrib);
324  } else {
325  J.reset(jac_contrib.release());
326  if (factor != 1.0) (*J) *= factor;
327  }
328  };
329 
330  auto jacs = jacobianFunctions(std::make_integer_sequence<int, sizeof...(input_indices)>{}, time_info.time(),
331  shape_disp, fields);
332 
333  for (size_t input_col = 0; input_col < jacobian_weights.size(); ++input_col) {
334  if (jacobian_weights[input_col] != 0.0) {
335  auto K = smith::get<DERIVATIVE>(jacs[input_col](time_info.time(), shape_disp, fields));
336  addToJ(jacobian_weights[input_col], assemble(K));
337  }
338  }
339 
340  return J;
341  }
342 
344  void jvp(TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector<ConstFieldPtr>& fields,
345  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& quad_fields,
346  [[maybe_unused]] ConstFieldPtr v_shape_disp, const std::vector<ConstFieldPtr>& v_fields,
347  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& v_quad_fields,
348  DualFieldPtr jvp_reaction) const override
349  {
350  validateFields(fields, "jvp");
351  SLIC_ERROR_IF(v_fields.size() != fields.size(),
352  "Invalid number of field sensitivities relative to the number of fields");
353 
354  dt_ = time_info.dt();
355  cycle_ = time_info.cycle();
356 
357  auto jacs = jacobianFunctions(std::make_integer_sequence<int, sizeof...(input_indices)>{}, time_info.time(),
358  shape_disp, fields);
359 
360  *jvp_reaction = 0.0;
361 
362  for (size_t input_col = 0; input_col < fields.size(); ++input_col) {
363  if (v_fields[input_col] != nullptr) {
364  auto K = smith::get<DERIVATIVE>(jacs[input_col](time_info.time(), shape_disp, fields));
365  K.AddMult(*v_fields[input_col], *jvp_reaction);
366  }
367  }
368  }
369 
371  void vjp(TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector<ConstFieldPtr>& fields,
372  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& quad_fields, ConstFieldPtr v_field,
373  DualFieldPtr vjp_shape_disp_sensitivity, const std::vector<DualFieldPtr>& vjp_sensitivities,
374  [[maybe_unused]] const std::vector<QuadratureFieldPtr>& vjp_quad_field_sensitivities) const override
375  {
376  validateFields(fields, "vjp");
377  SLIC_ERROR_IF(vjp_sensitivities.size() != fields.size(),
378  "Invalid number of field sensitivities relative to the number of fields");
379 
380  dt_ = time_info.dt();
381  cycle_ = time_info.cycle();
382 
383  auto vecJacs = vectorJacobianFunctions(std::make_integer_sequence<int, sizeof...(input_indices)>{},
384  time_info.time(), shape_disp, v_field, fields);
385  {
386  auto shape_vjp = smith::get<DERIVATIVE>((*v_dot_weak_form_residual_)(
387  DifferentiateWRT<0>{}, time_info.time(), *shape_disp, *v_field, *fields[input_indices]...));
388  auto shape_vjp_vector = assemble(shape_vjp);
389  *vjp_shape_disp_sensitivity += *shape_vjp_vector;
390  }
391 
392  for (size_t input_col = 0; input_col < fields.size(); ++input_col) {
393  if (vjp_sensitivities[input_col] != nullptr) {
394  auto vec_jac = smith::get<DERIVATIVE>(vecJacs[input_col](time_info.time(), shape_disp, v_field, fields));
395  auto vec_jac_mfem_vector = assemble(vec_jac);
396  *vjp_sensitivities[input_col] += *vec_jac_mfem_vector;
397  }
398  }
399  }
400 
403  ShapeAwareFunctional<ShapeDispSpace, OutputSpace(InputSpaces...)>& getShapeAwareResidual() { return *weak_form_; }
404 
408  ShapeAwareFunctional<ShapeDispSpace, double(OutputSpace, InputSpaces...)>& getShapeAwareVectorTimesResidual()
409  {
410  return *v_dot_weak_form_residual_;
411  }
412 
413  protected:
415  template <size_t I>
416  void validateInputSpaces(const SpacesT& input_mfem_spaces) const
417  {
418  if constexpr (I < sizeof...(InputSpaces)) {
419  using Space = typename std::tuple_element<I, std::tuple<InputSpaces...>>::type;
420  validateSpace<Space>(*input_mfem_spaces[I], axom::fmt::format("input[{}]", I));
421  validateInputSpaces<I + 1>(input_mfem_spaces);
422  }
423  }
424 
426  template <size_t I>
427  void validateFieldsRecursive(const std::vector<ConstFieldPtr>& fields, const std::string& method_name) const
428  {
429  if constexpr (I < sizeof...(InputSpaces)) {
430  using Space = typename std::tuple_element<I, std::tuple<InputSpaces...>>::type;
431  validateSpace<Space>(fields[I]->space(),
432  axom::fmt::format("{}(): field[{}] ('{}')", method_name, I, fields[I]->name()));
433  validateFieldsRecursive<I + 1>(fields, method_name);
434  }
435  }
436 
438  void validateFields(const std::vector<ConstFieldPtr>& fields, const std::string& method_name) const
439  {
440  SLIC_ERROR_ROOT_IF(fields.size() != sizeof...(InputSpaces),
441  axom::fmt::format("{}(): fields.size()={} but weak form expects {} InputSpaces", method_name,
442  fields.size(), sizeof...(InputSpaces)));
443 
444  if constexpr (sizeof...(InputSpaces) > 0) {
445  validateFieldsRecursive<0>(fields, method_name);
446  }
447  }
448 
450  template <typename Space>
451  static void validateSpace(const mfem::ParFiniteElementSpace& mfem_space, const std::string& space_name)
452  {
453  const auto* fec = mfem_space.FEColl();
454 
455  // Note: We check using the name string rather than dynamic_cast because MFEM
456  // has multiple FECollection types that represent H1 and L2 spaces
457  if constexpr (Space::family == Family::H1) {
458  std::string fec_name = fec->Name();
459  bool is_h1 =
460  (fec_name.find("H1") != std::string::npos || fec_name.find("ND_") != std::string::npos || // Nedelec elements
461  fec_name.find("Linear") != std::string::npos);
462  SLIC_ERROR_ROOT_IF(!is_h1, axom::fmt::format("Space '{}': Template specifies H1 family but mfem space uses '{}'",
463  space_name, fec_name));
464  } else if constexpr (Space::family == Family::L2) {
465  std::string fec_name = fec->Name();
466  bool is_l2 = (fec_name.find("L2") != std::string::npos || fec_name.find("DG") != std::string::npos ||
467  fec_name.find("Const") != std::string::npos);
468  SLIC_ERROR_ROOT_IF(
469  !is_l2, axom::fmt::format("Space '{}': Template specifies L2/DG family but mfem space uses '{}'", space_name,
470  fec_name));
471  }
472 
473  // Validate order
474  SLIC_ERROR_ROOT_IF(fec->GetOrder() != Space::order,
475  axom::fmt::format("Space '{}': Template specifies order {} but mfem space has order {}",
476  space_name, Space::order, fec->GetOrder()));
477 
478  // Validate vector dimension
479  SLIC_ERROR_ROOT_IF(
480  mfem_space.GetVDim() != Space::components,
481  axom::fmt::format("Space '{}': Template specifies {} components but mfem space has {} components (VDim)",
482  space_name, Space::components, mfem_space.GetVDim()));
483  }
484 
486  template <int... i>
487  auto jacobianFunctions(std::integer_sequence<int, i...>, double time, ConstFieldPtr shape_disp,
488  const std::vector<ConstFieldPtr>& fs) const
489  {
490  using JacFuncType = std::function<decltype((*weak_form_)(DifferentiateWRT<1>{}, time, *shape_disp, *fs[i]...))(
491  double, ConstFieldPtr, const std::vector<ConstFieldPtr>&)>;
492  return std::array<JacFuncType, sizeof...(i)>{
493  [this](double _time, ConstFieldPtr _shape_disp, const std::vector<ConstFieldPtr>& _fs) {
494  return (*weak_form_)(DifferentiateWRT<i + 1>{}, _time, *_shape_disp, *_fs[i]...);
495  }...};
496  }
497 
499  template <int... i>
500  auto vectorJacobianFunctions(std::integer_sequence<int, i...>, double time, ConstFieldPtr shape_disp, ConstFieldPtr v,
501  const std::vector<ConstFieldPtr>& fs) const
502  {
503  using GradFuncType =
504  std::function<decltype((*v_dot_weak_form_residual_)(DifferentiateWRT<1>{}, time, *shape_disp, *v, *fs[i]...))(
505  double, ConstFieldPtr, ConstFieldPtr, const std::vector<ConstFieldPtr>&)>;
506  return std::array<GradFuncType, sizeof...(i)>{
507  [this](double _time, ConstFieldPtr _shape_disp, ConstFieldPtr _v, const std::vector<ConstFieldPtr>& _fs) {
508  return (*v_dot_weak_form_residual_)(DifferentiateWRT<i + 2>{}, _time, *_shape_disp, *_v, *_fs[i]...);
509  }...};
510  }
511 
513  mutable double dt_ = std::numeric_limits<double>::max();
514 
516  mutable size_t cycle_ = 0;
517 
519  std::shared_ptr<Mesh> mesh_;
520 
522  std::unique_ptr<ShapeAwareFunctional<ShapeDispSpace, OutputSpace(InputSpaces...)>> weak_form_;
523 
525  std::unique_ptr<ShapeAwareFunctional<ShapeDispSpace, double(OutputSpace, InputSpaces...)>> v_dot_weak_form_residual_;
526 };
527 
530 inline std::vector<const mfem::ParFiniteElementSpace*> getSpaces(const std::vector<smith::FiniteElementState>& states)
531 {
532  std::vector<const mfem::ParFiniteElementSpace*> spaces;
533  for (auto& f : states) {
534  spaces.push_back(&f.space());
535  }
536  return spaces;
537 }
538 
539 } // namespace smith
Class for encapsulating the dual vector space of a finite element space (i.e. the space of linear for...
void jvp(TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector< ConstFieldPtr > &fields, [[maybe_unused]] const std::vector< ConstQuadratureFieldPtr > &quad_fields, [[maybe_unused]] ConstFieldPtr v_shape_disp, const std::vector< ConstFieldPtr > &v_fields, [[maybe_unused]] const std::vector< ConstQuadratureFieldPtr > &v_quad_fields, DualFieldPtr jvp_reaction) const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
auto vectorJacobianFunctions(std::integer_sequence< int, i... >, double time, ConstFieldPtr shape_disp, ConstFieldPtr v, const std::vector< ConstFieldPtr > &fs) const
Utility to get array of jvp functions, one for each input field in fs.
void validateFields(const std::vector< ConstFieldPtr > &fields, const std::string &method_name) const
Validate that fields vector matches the templated input spaces.
std::unique_ptr< ShapeAwareFunctional< ShapeDispSpace, double(OutputSpace, InputSpaces...)> > v_dot_weak_form_residual_
functional residual times and arbitrary vector v (same space as residual) evaluator,...
void addBodyIntegral(std::string body_name, BodyForceType body_integral)
This is an overloaded member function, provided for convenience. It differs from the above function o...
std::unique_ptr< ShapeAwareFunctional< ShapeDispSpace, OutputSpace(InputSpaces...)> > weak_form_
functional residual evaluator, shape aware
void addInteriorBoundaryIntegral(std::string interior_name, const InteriorIntegrandType &integrand)
This is an overloaded member function, provided for convenience. It differs from the above function o...
ShapeAwareFunctional< ShapeDispSpace, OutputSpace(InputSpaces...)> & getShapeAwareResidual()
Accessor to get a reference to the underlying ShapeAwareFunctional in case more direct access is need...
void validateFieldsRecursive(const std::vector< ConstFieldPtr > &fields, const std::string &method_name) const
Helper to validate fields passed to residual/jacobian/vjp/jvp.
FunctionalWeakForm(std::string physics_name, std::shared_ptr< Mesh > mesh, const mfem::ParFiniteElementSpace &output_mfem_space, const SpacesT &input_mfem_spaces)
Construct a new FunctionalWeakForm object.
void addBodySource(std::string body_name, BodyLoadType load_function)
This is an overloaded member function, provided for convenience. It differs from the above function o...
void addBoundaryIntegral(DependsOn< active_parameters... >, std::string boundary_name, BoundaryIntegrandType integrand)
Add a boundary integral term to the weak form.
auto jacobianFunctions(std::integer_sequence< int, i... >, double time, ConstFieldPtr shape_disp, const std::vector< ConstFieldPtr > &fs) const
Utility to get array of jacobian functions, one for each input field in fs.
void addBoundaryIntegral(std::string boundary_name, const BoundaryIntegrandType &integrand)
This is an overloaded member function, provided for convenience. It differs from the above function o...
std::unique_ptr< mfem::HypreParMatrix > jacobian(TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector< ConstFieldPtr > &fields, const std::vector< double > &jacobian_weights, [[maybe_unused]] const std::vector< ConstQuadratureFieldPtr > &quad_fields={}) const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
void addBoundaryFlux(std::string boundary_name, const BoundaryFluxType &integrand)
This is an overloaded member function, provided for convenience. It differs from the above function o...
void addBodySource(DependsOn< active_parameters... > depends_on, std::string body_name, BodyLoadType load_function)
Add a body source (body load) to the weak form.
void addBoundaryFlux(DependsOn< active_parameters... > depends_on, std::string boundary_name, BoundaryFluxType flux_function)
Add a boundary flux term to the weak form.
mfem::Vector residual(TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector< ConstFieldPtr > &fields, [[maybe_unused]] const std::vector< ConstQuadratureFieldPtr > &quad_fields={}) const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
void addInteriorBoundaryIntegral(DependsOn< active_parameters... >, std::string interior_name, InteriorIntegrandType integrand)
Add a interior boundary integral term to the weak form.
void vjp(TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector< ConstFieldPtr > &fields, [[maybe_unused]] const std::vector< ConstQuadratureFieldPtr > &quad_fields, ConstFieldPtr v_field, DualFieldPtr vjp_shape_disp_sensitivity, const std::vector< DualFieldPtr > &vjp_sensitivities, [[maybe_unused]] const std::vector< QuadratureFieldPtr > &vjp_quad_field_sensitivities) const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
static void validateSpace(const mfem::ParFiniteElementSpace &mfem_space, const std::string &space_name)
Validate that an mfem space matches the template space parameters.
void addBodyIntegral(DependsOn< active_parameters... >, std::string body_name, BodyIntegralType integrand)
Add a body integral contribution to the residual.
void validateInputSpaces(const SpacesT &input_mfem_spaces) const
Helper to validate input spaces recursively (for constructor)
ShapeAwareFunctional< ShapeDispSpace, double(OutputSpace, InputSpaces...)> & getShapeAwareVectorTimesResidual()
Accessor to get a reference to the underlying ShapeAwareFunctional vector-residual in case more direc...
Abstract WeakForm class.
Definition: weak_form.hpp:36
This contains a class that represents the dual of a finite element vector space, i....
This file contains the declaration of structure that manages the MFEM objects that make up the state ...
Smith mesh class which assists in constructing the appropriate parallel mfem meshes and registering a...
constexpr SMITH_HOST_DEVICE void for_constexpr(const lambda &f)
multidimensional loop tool that evaluates the lambda body inside the innermost loop.
Accelerator functionality.
Definition: smith.cpp:36
std::vector< const mfem::ParFiniteElementSpace * > getSpaces(const std::vector< smith::FiniteElementState > &states)
Helper function to construct vector of spaces from an existing vector of FiniteElementState.
std::vector< const mfem::ParFiniteElementSpace * > spaces(const std::vector< FieldState > &states, const std::vector< FieldState > &params={})
Get the spaces from the primal fields of a vector of field states.
SMITH_HOST_DEVICE auto cross(const tensor< T, 3, 2 > &A)
compute the cross product of the columns of A: A(:,1) x A(:,2)
Definition: tensor.hpp:966
SMITH_HOST_DEVICE auto max(dual< gradient_type > a, double b)
Implementation of max for dual numbers.
Definition: dual.hpp:229
mfem::future::tuple< T... > tuple
Expose MFEM tuple in the Smith namespace.
Definition: tuple.hpp:241
constexpr SMITH_HOST_DEVICE auto inner(const dual< S > &A, const dual< T > &B)
Definition: dual.hpp:281
mfem::future::tuple_element< I, T > tuple_element
Alias for the MFEM tuple element trait.
Definition: tuple.hpp:258
mfem::ParFiniteElementSpace & space(FieldState field)
Get the space from the primal field of a field states.
SMITH_HOST_DEVICE auto normalize(const tensor< T, n... > &A)
Normalizes the tensor Each element is divided by the Frobenius norm of the tensor,...
Definition: tensor.hpp:1122
FiniteElementState const * ConstFieldPtr
using
Definition: field_types.hpp:36
Wrapper of smith::Functional for evaluating integrals and derivatives of quantities with shape displa...
Compile-time alias for a dimension.
Definition: geometry.hpp:17
H1 elements of order p.
a struct that is used in the physics modules to clarify which template arguments are user-controlled ...
Definition: common.hpp:45
struct storing time and timestep information
Definition: common.hpp:18
double dt() const
accessor for dt
Definition: common.hpp:29
size_t cycle() const
accessor for cycle
Definition: common.hpp:32
double time() const
accessor for the current time
Definition: common.hpp:26
A sentinel struct for eliding no-op tensor operations.
Definition: tensor.hpp:122
Specifies interface for evaluating weak form residuals and their gradients.