Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
functional_qoi.inl
Go to the documentation of this file.
1 // Copyright (c) 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 
13 namespace serac {
14 
21  QoIProlongation() {}
22 
24  QoIProlongation(MPI_Comm c) : comm(c) {}
25 
27  void Mult(const mfem::Vector&, mfem::Vector&) const { SLIC_ERROR_ROOT("QoIProlongation::Mult() is not defined"); }
28 
30  void MultTranspose(const mfem::Vector& input, mfem::Vector& output) const
31  {
32  // const_cast to work around clang@14.0.6 compiler error:
33  // "argument type 'const double *' doesn't match specified 'MPI' type tag that requires 'double *'"
34  MPI_Allreduce(const_cast<double*>(&input[0]), &output[0], 1, MPI_DOUBLE, MPI_SUM, comm);
35  }
36 
37  MPI_Comm comm;
38 };
39 
52  void ScatterAdd(const mfem::Vector& input, mfem::Vector& output) const { output[0] += input.Sum(); }
53 };
54 
58 template <typename... trials, ExecutionSpace exec>
59 class Functional<double(trials...), exec> {
60  using test = QOI;
61  static constexpr tuple<trials...> trial_spaces{};
62  static constexpr uint32_t num_trial_spaces = sizeof...(trials);
63  static constexpr auto Q = std::max({test::order, trials::order...}) + 1;
64 
65  class Gradient;
66 
67  // clang-format off
68  template <uint32_t i>
69  struct operator_paren_return {
70  using type = typename std::conditional<
71  i == NO_DIFFERENTIATION, // if `i` is greater than or equal to zero,
72  double, // wise, we just return the value
73  serac::tuple<double&, Gradient&> // otherwise, we return the value and the derivative w.r.t arg `i`
74  >::type;
75  };
76  // clang-format on
77 
78  public:
83  Functional(std::array<const mfem::ParFiniteElementSpace*, num_trial_spaces> trial_fes)
84  : test_fec_(0, trial_fes[0]->GetMesh()->Dimension()),
85  test_space_(dynamic_cast<mfem::ParMesh*>(trial_fes[0]->GetMesh()), &test_fec_, 1, serac::ordering),
86  trial_space_(trial_fes),
87  mem_type(mfem::Device::GetMemoryType())
88  {
90 
91  std::array<Family, num_trial_spaces> trial_families = {trials::family...};
92  std::array<int, num_trial_spaces> trial_orders = {trials::order...};
93  std::array<int, num_trial_spaces> trial_components = {trials::components...};
94 
95  for (uint32_t i = 0; i < num_trial_spaces; i++) {
96  trial_function_spaces_[i] = {trial_families[i], trial_orders[i], trial_components[i]};
97 
98  P_trial_[i] = trial_space_[i]->GetProlongationMatrix();
99 
100  // For L2 spaces, configure a ParGridFunction with external data to exchange face neighbor data.
101  // This is for DG method where interior faces on the processor boundary need to access data on
102  // the neighrbor element
103  if (trial_function_spaces_[i].family == Family::L2) {
104  // Set the vector size to store tdof data
105  input_ldof_values_[i].SetSize(P_trial_[i]->Height(), mem_type);
106 
107  // Initialize face neighbor data vector
108  updateFaceNbrData(trial_space_[i], trial_pargrid_functions_[i], input_ldof_values_[i]);
109 
110  // Set the local input vector size to store local + ghost data
111  input_L_[i].SetSize(input_ldof_values_[i].Size() + trial_pargrid_functions_[i].FaceNbrData().Size(), mem_type);
112  } else {
113  input_L_[i].SetSize(P_trial_[i]->Height(), mem_type);
114  }
115 
116  // create the necessary number of empty mfem::Vectors, to be resized later
117  input_E_.push_back({});
118  input_E_buffer_.push_back({});
119  }
120 
121  G_test_ = QoIElementRestriction();
122  P_test_ = QoIProlongation(trial_fes[0]->GetParMesh()->GetComm());
123 
124  output_L_.SetSize(1, mem_type);
125 
126  output_T_.SetSize(1, mem_type);
127 
128  // gradient objects depend on some member variables in
129  // Functional, so we initialize the gradient objects last
130  // to ensure that those member variables are initialized first
131  for (uint32_t i = 0; i < num_trial_spaces; i++) {
132  grad_.emplace_back(*this, i);
133  }
134  }
135 
148  template <int dim, int... args, typename lambda, typename qpt_data_type = Nothing>
149  void AddDomainIntegral(Dimension<dim>, DependsOn<args...>, const lambda& integrand, Domain& domain,
150  std::shared_ptr<QuadratureData<qpt_data_type>> qdata = NoQData)
151  {
152  if (domain.mesh_.GetNE() == 0) return;
153 
154  SLIC_ERROR_ROOT_IF(dim != domain.mesh_.Dimension(), "invalid mesh dimension for domain integral");
155 
158 
159  std::vector<uint32_t> arg_vec = {args...};
160  for (uint32_t i : arg_vec) {
161  domain.insert_restriction(trial_space_[i], trial_function_spaces_[i]);
162  }
163 
164  using signature = test(decltype(serac::type<args>(trial_spaces))...);
165  integrals_.push_back(
166  MakeDomainIntegral<signature, Q, dim>(domain, integrand, qdata, std::vector<uint32_t>{args...}));
167  }
168 
180  template <int dim, int... args, typename lambda>
181  void AddBoundaryIntegral(Dimension<dim>, DependsOn<args...>, const lambda& integrand, Domain& domain)
182  {
183  auto num_bdr_elements = domain.mesh_.GetNBE();
184  if (num_bdr_elements == 0) return;
185 
186  SLIC_ERROR_ROOT_IF(dim != domain.dim_, "invalid domain of integration for boundary integral");
187 
189 
190  std::vector<uint32_t> arg_vec = {args...};
191  for (uint32_t i : arg_vec) {
192  domain.insert_restriction(trial_space_[i], trial_function_spaces_[i]);
193  }
194 
195  using signature = test(decltype(serac::type<args>(trial_spaces))...);
196  integrals_.push_back(MakeBoundaryIntegral<signature, Q, dim>(domain, integrand, arg_vec));
197  }
198 
210  template <int dim, int... args, typename Integrand>
211  void AddInteriorFaceIntegral(Dimension<dim>, DependsOn<args...>, const Integrand& integrand, Domain& domain)
212  {
214 
215  std::vector<uint32_t> arg_vec = {args...};
216  for (uint32_t i : arg_vec) {
217  domain.insert_restriction(trial_space_[i], trial_function_spaces_[i]);
218  }
219 
220  using signature = test(decltype(serac::type<args>(trial_spaces))...);
221  integrals_.push_back(MakeInteriorFaceIntegral<signature, Q, dim>(domain, integrand, arg_vec));
222  }
223 
234  template <int... args, typename lambda, typename qpt_data_type = Nothing>
235  void AddAreaIntegral(DependsOn<args...> which_args, const lambda& integrand, Domain& domain,
236  std::shared_ptr<QuadratureData<qpt_data_type>>& data = NoQData)
237  {
238  AddDomainIntegral(Dimension<2>{}, which_args, integrand, domain, data);
239  }
240 
251  template <int... args, typename lambda, typename qpt_data_type = Nothing>
252  void AddVolumeIntegral(DependsOn<args...> which_args, const lambda& integrand, Domain& domain,
253  std::shared_ptr<QuadratureData<qpt_data_type>>& data = NoQData)
254  {
255  AddDomainIntegral(Dimension<3>{}, which_args, integrand, domain, data);
256  }
257 
259  template <int... args, typename lambda>
260  void AddSurfaceIntegral(DependsOn<args...> which_args, const lambda& integrand, Domain& domain)
261  {
262  AddBoundaryIntegral(Dimension<2>{}, which_args, integrand, domain);
263  }
264 
275  double ActionOfGradient(const mfem::Vector& input_T, uint32_t which) const
276  {
277  // Please refer to 'ActionOfGradient' in 'functional.hpp' for detailed comments of this implementation
278  if (trial_function_spaces_[which].family == Family::L2) {
279  P_trial_[which]->Mult(input_T, input_ldof_values_[which]);
280  input_L_[which].SetVector(input_ldof_values_[which], 0);
281  updateFaceNbrData(trial_space_[which], trial_pargrid_functions_[which], input_ldof_values_[which]);
282 
283  if (trial_space_[which]->GetVDim() == 1) {
284  input_L_[which].SetVector(trial_pargrid_functions_[which].FaceNbrData(), input_ldof_values_[which].Size());
285  } else {
286  if (trial_space_[which]->GetOrdering() == mfem::Ordering::byVDIM) {
287  int local_size = input_ldof_values_[which].Size();
288  appendFaceNbrData(trial_space_[which], trial_pargrid_functions_[which], local_size, input_L_[which]);
289  } else {
290  SLIC_ERROR_ROOT("Unsupported: L2 vector field ordered by nodes");
291  }
292  }
293  } else {
294  P_trial_[which]->Mult(input_T, input_L_[which]);
295  }
296 
297  output_L_ = 0.0;
298 
299  for (auto& integral : integrals_) {
300  if (integral.DependsOn(which)) {
301  Domain& dom = integral.domain_;
302 
303  const serac::BlockElementRestriction& G_trial = dom.get_restriction(trial_function_spaces_[which]);
304  input_E_buffer_[which].SetSize(int(G_trial.ESize()));
305  input_E_[which].Update(input_E_buffer_[which], G_trial.bOffsets());
306  G_trial.Gather(input_L_[which], input_E_[which]);
307 
308  output_E_buffer_.SetSize(dom.total_elements());
309  output_E_.Update(output_E_buffer_, dom.bOffsets());
310 
311  integral.GradientMult(input_E_[which], output_E_, which);
312 
313  // make sure shared interior faces are integrated only once
314  if (dom.type_ == Domain::Type::InteriorFaces) {
315  for (int shared_id : dom.shared_interior_face_ids_) {
316  output_E_[shared_id] *= 0.5;
317  }
318  }
319 
320  // scatter-add to compute QoI value for the local processor
321  G_test_.ScatterAdd(output_E_, output_L_);
322  }
323  }
324 
325  // compute global QoI value by summing values from different processors
326  P_test_.MultTranspose(output_L_, output_T_);
327 
328  return output_T_[0];
329  }
330 
341  template <uint32_t wrt, typename... T>
342  typename operator_paren_return<wrt>::type operator()(DifferentiateWRT<wrt>, double t, const T&... args)
343  {
344  const mfem::Vector* input_T[] = {&static_cast<const mfem::Vector&>(args)...};
345 
346  // Please refer to 'operator()' in 'functional.hpp' for detailed comments of this implementation
347  // get the values for each local processor
348  for (uint32_t i = 0; i < num_trial_spaces; i++) {
349  if (trial_function_spaces_[i].family == Family::L2) {
350  P_trial_[i]->Mult(*input_T[i], input_ldof_values_[i]);
351  input_L_[i].SetVector(input_ldof_values_[i], 0);
352  updateFaceNbrData(trial_space_[i], trial_pargrid_functions_[i], input_ldof_values_[i]);
353 
354  if (trial_space_[i]->GetVDim() == 1) {
355  input_L_[i].SetVector(trial_pargrid_functions_[i].FaceNbrData(), input_ldof_values_[i].Size());
356  } else {
357  if (trial_space_[i]->GetOrdering() == mfem::Ordering::byVDIM) {
358  int local_size = input_ldof_values_[i].Size();
359  appendFaceNbrData(trial_space_[i], trial_pargrid_functions_[i], local_size, input_L_[i]);
360  } else {
361  SLIC_ERROR_ROOT("Unsupported: L2 vector field ordered by nodes");
362  }
363  }
364  } else {
365  P_trial_[i]->Mult(*input_T[i], input_L_[i]);
366  }
367  }
368 
369  output_L_ = 0.0;
370 
371  for (auto& integral : integrals_) {
372  Domain& dom = integral.domain_;
373 
374  for (auto i : integral.active_trial_spaces_) {
375  const serac::BlockElementRestriction& G_trial = dom.get_restriction(trial_function_spaces_[i]);
376  input_E_buffer_[i].SetSize(int(G_trial.ESize()));
377  input_E_[i].Update(input_E_buffer_[i], G_trial.bOffsets());
378  G_trial.Gather(input_L_[i], input_E_[i]);
379  }
380 
381  output_E_buffer_.SetSize(dom.total_elements());
382  output_E_.Update(output_E_buffer_, dom.bOffsets());
383 
384  const bool update_qdata = false;
385  integral.Mult(t, input_E_, output_E_, wrt, update_qdata);
386 
387  // make sure shared interior faces are integrated only once
388  if (dom.type_ == Domain::Type::InteriorFaces) {
389  for (int shared_id : dom.shared_interior_face_ids_) {
390  output_E_[shared_id] *= 0.5;
391  }
392  }
393 
394  // scatter-add to compute QoI value for the local processor
395  G_test_.ScatterAdd(output_E_, output_L_);
396  }
397 
398  // compute global QoI value by summing values from different processors
399  P_test_.MultTranspose(output_L_, output_T_);
400 
401  if constexpr (wrt != NO_DIFFERENTIATION) {
402  // if the user has indicated they'd like to evaluate and differentiate w.r.t.
403  // a specific argument, then we return both the value and gradient w.r.t. that argument
404  //
405  // mfem::Vector arg0 = ...;
406  // mfem::Vector arg1 = ...;
407  // e.g. auto [value, gradient_wrt_arg1] = my_functional(arg0, differentiate_wrt(arg1));
408  return {output_T_[0], grad_[wrt]};
409  }
410 
411  if constexpr (wrt == NO_DIFFERENTIATION) {
412  // if the user passes only `mfem::Vector`s then we assume they only want the output value
413  //
414  // mfem::Vector arg0 = ...;
415  // mfem::Vector arg1 = ...;
416  // e.g. mfem::Vector value = my_functional(arg0, arg1);
417  return output_T_[0];
418  }
419  }
420 
422  template <typename... T>
423  auto operator()(double t, const T&... args)
424  {
425  // below we add 0 so the number of differentiated arguments defaults to 0 if trial spaces are not provided
426  constexpr int num_differentiated_arguments = (std::is_same_v<T, differentiate_wrt_this> + ... + 0);
427  static_assert(num_differentiated_arguments <= 1,
428  "Error: Functional::operator() can only differentiate w.r.t. 1 argument a time");
429  static_assert(sizeof...(T) == num_trial_spaces,
430  "Error: Functional::operator() must take exactly as many arguments as trial spaces");
431 
432  [[maybe_unused]] constexpr uint32_t i = index_of_differentiation<T...>();
433 
434  return (*this)(DifferentiateWRT<i>{}, t, args...);
435  }
436 
437  private:
441  enum class Operation
442  {
443  Mult,
444  GradientMult
445  };
446 
450  class Gradient {
451  public:
456  Gradient(Functional<double(trials...)>& f, uint32_t which = 0)
457  : form_(f), which_argument(which), gradient_L_(f.trial_space_[which]->GetVSize())
458  {
459  }
460 
461  double operator()(const mfem::Vector& x) const { return form_.ActionOfGradient(x, which_argument); }
462 
463  uint64_t max_buffer_size()
464  {
465  uint64_t max_entries = 0;
466  for (auto& integral : form_.integrals_) {
467  if (integral.DependsOn(which_argument)) {
468  Domain& dom = integral.domain_;
469  const auto& G_trial = dom.get_restriction(form_.trial_function_spaces_[which_argument]);
470  for (const auto& [geom, test_restriction] : G_trial.restrictions) {
471  const auto& trial_restriction = G_trial.restrictions.at(geom);
472  uint64_t entries_per_element = trial_restriction.nodes_per_elem * trial_restriction.components;
473  max_entries = std::max(max_entries, trial_restriction.num_elements * entries_per_element);
474  }
475  }
476  }
477  return max_entries;
478  }
479 
480  std::unique_ptr<mfem::HypreParVector> assemble()
481  {
482  // The mfem method ParFiniteElementSpace.NewTrueDofVector should really be marked const
483  std::unique_ptr<mfem::HypreParVector> gradient_T(
484  const_cast<mfem::ParFiniteElementSpace*>(form_.trial_space_[which_argument])->NewTrueDofVector());
485 
486  gradient_L_ = 0.0;
487 
488  std::vector<double> K_elem_buffer(max_buffer_size());
489 
490  std::map<mfem::Geometry::Type, ExecArray<double, 3, exec>> element_gradients[Domain::num_types];
491 
492  int lcols = form_.trial_space_[which_argument]->GetVSize();
493 
495 
496  for (auto& integral : form_.integrals_) {
497  Domain& dom = integral.domain_;
498 
499  // if this integral's derivative isn't identically zero
500  if (integral.functional_to_integral_index_.count(which_argument) > 0) {
501  uint32_t id = integral.functional_to_integral_index_.at(which_argument);
502  const auto& G_trial = dom.get_restriction(form_.trial_function_spaces_[which_argument]);
503  for (const auto& [geom, calculate_element_gradients] : integral.element_gradient_[id]) {
504  const auto& trial_restriction = G_trial.restrictions.at(geom);
505 
506  // prepare a buffer to hold the element matrices
507  CPUArrayView<double, 3> K_e(K_elem_buffer.data(), trial_restriction.num_elements, 1,
508  trial_restriction.nodes_per_elem * trial_restriction.components);
509  detail::zero_out(K_e);
510 
511  calculate_element_gradients(K_e);
512 
513  const std::vector<int>& element_ids = integral.domain_.get(geom);
514 
515  uint32_t cols_per_elem = uint32_t(trial_restriction.nodes_per_elem * trial_restriction.components);
516  std::vector<DoF> trial_vdofs(cols_per_elem);
517 
518  for (uint32_t e = 0; e < element_ids.size(); e++) {
519  trial_restriction.GetElementVDofs(int(e), trial_vdofs);
520 
521  for (uint32_t i = 0; i < cols_per_elem; i++) {
522  int col = int(trial_vdofs[i].index());
523 
524  // only add local DG dof gradients (other FE spaces satisfy col < lcols inherently)
525  if (col < lcols) {
526  gradient_L_[col] += K_e(e, 0, i);
527  }
528  }
529  }
530  }
531  }
532  }
533 
535 
536  form_.P_trial_[which_argument]->MultTranspose(gradient_L_, *gradient_T);
537 
538  return gradient_T;
539  }
540 
541  friend auto assemble(Gradient& g) { return g.assemble(); }
542 
543  private:
547  Functional<double(trials...), exec>& form_;
548 
549  uint32_t which_argument;
550 
551  mfem::Vector gradient_L_;
552  };
553 
555  const mfem::L2_FECollection test_fec_;
556  const mfem::ParFiniteElementSpace test_space_;
557 
559  std::array<const mfem::ParFiniteElementSpace*, num_trial_spaces> trial_space_;
560 
561  std::array<FunctionSpace, num_trial_spaces> trial_function_spaces_;
562 
564  mutable std::array<mfem::ParGridFunction, num_trial_spaces> trial_pargrid_functions_;
565 
567  mutable std::array<mfem::Vector, num_trial_spaces> input_ldof_values_;
568 
573  const mfem::Operator* P_trial_[num_trial_spaces];
574 
576  mutable mfem::Vector input_L_[num_trial_spaces];
577 
578  mutable std::vector<mfem::Vector> input_E_buffer_;
579  mutable std::vector<mfem::BlockVector> input_E_;
580 
581  mutable std::vector<Integral> integrals_;
582 
583  mutable mfem::Vector output_E_buffer_;
584  mutable mfem::BlockVector output_E_;
585 
586  QoIElementRestriction G_test_;
587 
589  mutable mfem::Vector output_L_;
590 
591  QoIProlongation P_test_;
592 
594  mutable mfem::Vector output_T_;
595 
597  mutable std::vector<Gradient> grad_;
598 
599  const mfem::MemoryType mem_type;
600 };
601 
602 } // namespace serac
Functional(std::array< const mfem::ParFiniteElementSpace *, num_trial_spaces > trial_fes)
Constructs using a mfem::ParFiniteElementSpace object corresponding to the trial space.
auto operator()(double t, const T &... args)
This is an overloaded member function, provided for convenience. It differs from the above function o...
void AddInteriorFaceIntegral(Dimension< dim >, DependsOn< args... >, const Integrand &integrand, Domain &domain)
Adds a interior face integral term to the Functional object.
operator_paren_return< wrt >::type operator()(DifferentiateWRT< wrt >, double t, const T &... args)
this function lets the user evaluate the serac::Functional with the given trial space values
void AddAreaIntegral(DependsOn< args... > which_args, const lambda &integrand, Domain &domain, std::shared_ptr< QuadratureData< qpt_data_type >> &data=NoQData)
Adds an area integral, i.e., over 2D elements in R^2.
void AddDomainIntegral(Dimension< dim >, DependsOn< args... >, const lambda &integrand, Domain &domain, std::shared_ptr< QuadratureData< qpt_data_type >> qdata=NoQData)
Adds a domain integral term to the Functional object.
double ActionOfGradient(const mfem::Vector &input_T, uint32_t which) const
this function computes the directional derivative of the quantity of interest functional
void AddSurfaceIntegral(DependsOn< args... > which_args, const lambda &integrand, Domain &domain)
alias for Functional::AddBoundaryIntegral(Dimension<2>{}, integrand, domain);
void AddVolumeIntegral(DependsOn< args... > which_args, const lambda &integrand, Domain &domain, std::shared_ptr< QuadratureData< qpt_data_type >> &data=NoQData)
Adds a volume integral, i.e., over 3D elements in R^3.
void AddBoundaryIntegral(Dimension< dim >, DependsOn< args... >, const lambda &integrand, Domain &domain)
Adds a boundary integral term to the Functional object.
Accelerator functionality.
Definition: serac.cpp:36
std::shared_ptr< QuadratureData< Nothing > > NoQData
a single instance of a QuadratureData container of Nothings, since they are all interchangeable
SERAC_HOST_DEVICE auto max(dual< gradient_type > a, double b)
Implementation of max for dual numbers.
Definition: dual.hpp:229
void updateFaceNbrData(const mfem::ParFiniteElementSpace *const_trial_space, mfem::ParGridFunction &trial_pgf, mfem::Vector &trial_tdof_vals)
helper function to locally cast away const on FE space so we can update face neighbor data with Excha...
Definition: functional.hpp:151
void check_for_missing_nodal_gridfunc(const mfem::Mesh &mesh)
function for verifying that the mesh has been fully initialized
Definition: functional.hpp:84
constexpr uint32_t index_of_differentiation()
given a list of types, this function returns the index that corresponds to the type dual_vector.
Definition: functional.hpp:60
void appendFaceNbrData(const mfem::ParFiniteElementSpace *trial_space, const mfem::ParGridFunction &trial_pgf, const int LSize, mfem::Vector &input_L)
helper functional to reorder the ordering of FaceNbrData for L2 space to byVDIM and append this vecto...
Definition: functional.hpp:168
constexpr SERAC_HOST_DEVICE auto type(const tuple< T... > &values)
a function intended to be used for extracting the ith type from a tuple.
Definition: tuple.hpp:379
ExecutionSpace
enum used for signalling whether or not to perform certain calculations on the CPU or GPU
Definition: accelerator.hpp:73
void check_for_unsupported_elements(const mfem::Mesh &mesh)
function for verifying that there are no unsupported element types in the mesh
Definition: functional.hpp:93
#define SERAC_MARK_FUNCTION
Definition: profiling.hpp:90
a generalization of mfem::ElementRestriction that works with multiple kinds of element geometries....
mfem::Array< int > bOffsets() const
block offsets used when constructing mfem::HypreParVectors
uint64_t ESize() const
the size of the "E-vector" associated with this restriction operator
void Gather(const mfem::Vector &L_vector, mfem::BlockVector &E_block_vector) const
"L->E" in mfem parlance, each element gathers the values that belong to it, and stores them in the "E...
Compile-time alias for a dimension.
Definition: geometry.hpp:17
a class for representing a geometric region that can be used for integration
Definition: domain.hpp:33
void insert_restriction(const fes_t *fes, FunctionSpace space)
create a restriction operator over this domain, using its FunctionSpace as a key
Definition: domain.cpp:536
mfem::Array< int > bOffsets() const
returns an array of the prefix sum of element counts belonging to this domain. Primarily intended to ...
Definition: domain.hpp:206
static constexpr int num_types
the number of entries in the Type enum
Definition: domain.hpp:42
const BlockElementRestriction & get_restriction(FunctionSpace space)
getter for accessing a restriction operator by its function space
Definition: domain.cpp:544
const mesh_t & mesh_
the underyling mesh for this domain
Definition: domain.hpp:45
int total_elements() const
returns how many elements of any type belong to this domain
Definition: domain.hpp:197
Type type_
whether the elements in this domain are on the boundary or not
Definition: domain.hpp:51
int dim_
the geometric dimension of the domain
Definition: domain.hpp:48
std::vector< int > shared_interior_face_ids_
Ids of interior faces that lie on the boundary shared by two processors.
Definition: domain.hpp:106
these classes are a little confusing. These two special types represent the similar (but different) c...
"Quantity of Interest" elements (i.e. elements with a single shape function, 1)
this class behaves like a Restriction operator, except is specialized for the case of a quantity of i...
void ScatterAdd(const mfem::Vector &input, mfem::Vector &output) const
element-to-global ScatterAdd operation used in FEM assembly, for quantities of interest
this class behaves like a Prolongation operator, except is specialized for the case of a quantity of ...
void MultTranspose(const mfem::Vector &input, mfem::Vector &output) const
set the value of output to the distributed sum over input values from different processors
QoIProlongation(MPI_Comm c)
create a QoIProlongation for a Quantity of Interest
MPI_Comm comm
MPI communicator used to carry out the distributed reduction.
void Mult(const mfem::Vector &, mfem::Vector &) const
unimplemented: do not use
A class for storing and access user-defined types at quadrature points.
This is a class that mimics most of std::tuple's interface, except that it is usable in CUDA kernels ...
Definition: tuple.hpp:28