20 #include <type_traits>
26 #include "serac/serac_config.hpp"
31 #include "serac/numerics/functional/integral.hpp"
32 #include "serac/numerics/functional/differentiate_wrt.hpp"
33 #include "serac/numerics/functional/element_restriction.hpp"
39 constexpr
int SOURCE = 0;
40 constexpr
int FLUX = 1;
41 constexpr
int VALUE = 0;
42 constexpr
int DERIVATIVE = 1;
59 template <
typename... T>
62 constexpr uint32_t n =
sizeof...(T);
63 bool matching[] = {std::is_same_v<T, differentiate_wrt_this>...};
64 for (uint32_t i = 0; i < n; i++) {
69 return NO_DIFFERENTIATION;
80 constexpr
operator int() {
return ind; }
86 if (mesh.GetNodes() ==
nullptr) {
89 the provided mesh does not have a nodal gridfunction.
90 If you created an mfem::Mesh manually, make sure that the
91 following member functions are invoked before use
93 > mfem::Mesh::EnsureNodes();
94 > mfem::ParMesh::ExchangeFaceNbrData();
96 or else the mfem::Mesh won't be fully initialized
104 int num_elements = mesh.GetNE();
105 for (
int e = 0; e < num_elements; e++) {
106 auto type = mesh.GetElementType(e);
107 if (
type == mfem::Element::POINT ||
type == mfem::Element::WEDGE ||
type == mfem::Element::PYRAMID) {
108 SLIC_ERROR_ROOT(
"Mesh contains unsupported element type");
121 template <
typename function_space>
122 inline std::pair<std::unique_ptr<mfem::ParFiniteElementSpace>, std::unique_ptr<mfem::FiniteElementCollection>>
125 const int dim = mesh->Dimension();
126 std::unique_ptr<mfem::FiniteElementCollection> fec;
128 switch (function_space::family) {
130 fec = std::make_unique<mfem::H1_FECollection>(function_space::order, dim);
133 fec = std::make_unique<mfem::ND_FECollection>(function_space::order, dim);
136 fec = std::make_unique<mfem::RT_FECollection>(function_space::order, dim);
140 fec = std::make_unique<mfem::L2_FECollection>(function_space::order, dim, mfem::BasisType::GaussLobatto);
143 return {
nullptr,
nullptr};
149 std::make_unique<mfem::ParFiniteElementSpace>(mesh, fec.get(), function_space::components, serac::ordering);
151 return std::pair(std::move(fes), std::move(fec));
160 inline void updateFaceNbrData(
const mfem::ParFiniteElementSpace* const_trial_space, mfem::ParGridFunction& trial_pgf,
161 mfem::Vector& trial_tdof_vals)
163 mfem::ParFiniteElementSpace* nonconst_trial_space =
const_cast<mfem::ParFiniteElementSpace*
>(const_trial_space);
166 trial_pgf.MakeRef(nonconst_trial_space, trial_tdof_vals.GetData());
170 trial_pgf.ExchangeFaceNbrData();
177 inline void appendFaceNbrData(
const mfem::ParFiniteElementSpace* trial_space,
const mfem::ParGridFunction& trial_pgf,
178 const int LSize, mfem::Vector& input_L)
180 int num_ghost_elem = trial_space->GetParMesh()->GetNFaceNeighborElements();
181 int components_per_node = trial_space->GetVDim();
182 const mfem::Vector& face_neighbor_data = trial_pgf.FaceNbrData();
185 for (
int n = 0; n < num_ghost_elem; ++n) {
186 mfem::Array<int> old_shared_elem_vdof_ids;
187 trial_space->GetFaceNbrElementVDofs(n, old_shared_elem_vdof_ids);
188 int dofs_per_element = old_shared_elem_vdof_ids.Size() / components_per_node;
191 for (
int m = 0; m < old_shared_elem_vdof_ids.Size(); ++m) {
193 int new_id = (m % dofs_per_element) * components_per_node + m / dofs_per_element + offset;
194 input_L[LSize + new_id] = face_neighbor_data(old_shared_elem_vdof_ids[m]);
197 offset += old_shared_elem_vdof_ids.Size();
205 mfem::Array<HYPRE_BigInt>& face_nbr_glob_vdof_map)
207 int components_per_node = trial_space->GetVDim();
208 if (components_per_node == 1) {
209 face_nbr_glob_vdof_map = trial_space->face_nbr_glob_dof_map;
213 if (components_per_node > 1 && trial_space->GetOrdering() == mfem::Ordering::byNODES) {
214 SLIC_ERROR_ROOT(
"Unsupported: L2 vector field ordered by nodes");
217 face_nbr_glob_vdof_map.SetSize(trial_space->face_nbr_glob_dof_map.Size());
218 int num_ghost_elem = trial_space->GetParMesh()->GetNFaceNeighborElements();
219 const HYPRE_BigInt* face_nbr_glob_dof_index = trial_space->face_nbr_glob_dof_map;
222 for (
int n = 0; n < num_ghost_elem; ++n) {
223 mfem::Array<int> old_shared_elem_vdof_ids;
224 trial_space->GetFaceNbrElementVDofs(n, old_shared_elem_vdof_ids);
225 int dofs_per_element = old_shared_elem_vdof_ids.Size() / components_per_node;
228 for (
int m = 0; m < old_shared_elem_vdof_ids.Size(); ++m) {
229 int new_id = (m % dofs_per_element) * components_per_node + m / dofs_per_element + offset;
230 face_nbr_glob_vdof_map[new_id] = face_nbr_glob_dof_index[old_shared_elem_vdof_ids[m]];
233 offset += old_shared_elem_vdof_ids.Size();
238 template <
typename T, ExecutionSpace exec = serac::default_execution_space>
283 class Functional<test(trials...), exec> {
284 static constexpr
tuple<trials...> trial_spaces{};
285 static constexpr uint32_t num_trial_spaces =
sizeof...(trials);
286 static constexpr
auto Q =
std::max({test::order, trials::order...}) + 1;
288 static constexpr mfem::Geometry::Type elem_geom[4] = {mfem::Geometry::INVALID, mfem::Geometry::SEGMENT,
289 mfem::Geometry::SQUARE, mfem::Geometry::CUBE};
290 static constexpr mfem::Geometry::Type simplex_geom[4] = {mfem::Geometry::INVALID, mfem::Geometry::SEGMENT,
291 mfem::Geometry::TRIANGLE, mfem::Geometry::TETRAHEDRON};
296 template <u
int32_t i>
297 struct operator_paren_return {
298 using type =
typename std::conditional<
299 i == NO_DIFFERENTIATION,
312 Functional(
const mfem::ParFiniteElementSpace* test_fes,
313 std::array<const mfem::ParFiniteElementSpace*, num_trial_spaces> trial_fes)
314 : update_qdata_(false), test_space_(test_fes), trial_space_(trial_fes), mem_type(mfem::Device::GetMemoryType())
318 test_function_space_ = {test::family, test::order, test::components};
320 std::array<Family, num_trial_spaces> trial_families = {trials::family...};
321 std::array<int, num_trial_spaces> trial_orders = {trials::order...};
322 std::array<int, num_trial_spaces> trial_components = {trials::components...};
324 for (uint32_t i = 0; i < num_trial_spaces; i++) {
325 trial_function_spaces_[i] = {trial_families[i], trial_orders[i], trial_components[i]};
327 P_trial_[i] = trial_space_[i]->GetProlongationMatrix();
332 if (trial_function_spaces_[i].family == Family::L2) {
334 input_ldof_values_[i].SetSize(P_trial_[i]->Height(), mem_type);
337 updateFaceNbrData(trial_space_[i], trial_pargrid_functions_[i], input_ldof_values_[i]);
340 input_L_[i].SetSize(input_ldof_values_[i].Size() + trial_pargrid_functions_[i].FaceNbrData().Size(), mem_type);
346 input_L_[i].SetSize(P_trial_[i]->Height(), mem_type);
350 input_E_.push_back({});
351 input_E_buffer_.push_back({});
358 P_test_ = test_space_->GetProlongationMatrix();
360 if (test_function_space_.family == Family::L2) {
362 mfem::ParGridFunction X_pgf;
363 mfem::Vector X_ldof_values(P_test_->Height(), mem_type);
366 output_L_.SetSize(X_ldof_values.Size() + X_pgf.FaceNbrData().Size(), mem_type);
368 output_L_.SetSize(P_test_->Height(), mem_type);
371 output_T_.SetSize(test_fes->GetTrueVSize(), mem_type);
376 for (uint32_t i = 0; i < num_trial_spaces; i++) {
377 grad_.emplace_back(*
this, i);
392 template <
int dim,
int... args,
typename lambda,
typename qpt_data_type = Nothing>
393 void AddDomainIntegral(Dimension<dim>, DependsOn<args...>,
const lambda& integrand, Domain& domain,
394 std::shared_ptr<QuadratureData<qpt_data_type>> qdata =
NoQData)
396 if (domain.mesh_.GetNE() == 0)
return;
398 SLIC_ERROR_ROOT_IF(dim != domain.mesh_.Dimension(),
"invalid mesh dimension for domain integral");
403 std::vector<uint32_t> arg_vec = {args...};
404 for (uint32_t i : arg_vec) {
405 domain.insert_restriction(trial_space_[i], trial_function_spaces_[i]);
407 domain.insert_restriction(test_space_, test_function_space_);
409 using signature = test(decltype(serac::type<args>(trial_spaces))...);
410 integrals_.push_back(
411 MakeDomainIntegral<signature, Q, dim>(domain, integrand, qdata, std::vector<uint32_t>{args...}));
423 template <
int dim,
int... args,
typename lambda>
424 void AddBoundaryIntegral(Dimension<dim>, DependsOn<args...>,
const lambda& integrand, Domain& domain)
426 auto num_bdr_elements = domain.mesh_.GetNBE();
427 if (num_bdr_elements == 0)
return;
429 SLIC_ERROR_ROOT_IF(dim != domain.dim_,
"invalid domain of integration for boundary integral");
433 std::vector<uint32_t> arg_vec = {args...};
434 for (uint32_t i : arg_vec) {
435 domain.insert_restriction(trial_space_[i], trial_function_spaces_[i]);
437 domain.insert_restriction(test_space_, test_function_space_);
439 using signature = test(decltype(serac::type<args>(trial_spaces))...);
440 integrals_.push_back(MakeBoundaryIntegral<signature, Q, dim>(domain, integrand, std::vector<uint32_t>{args...}));
446 template <
int dim,
int... args,
typename Integrand>
447 void AddInteriorFaceIntegral(Dimension<dim>, DependsOn<args...>,
const Integrand& integrand, Domain& domain)
451 std::vector<uint32_t> arg_vec = {args...};
452 for (uint32_t i : arg_vec) {
453 domain.insert_restriction(trial_space_[i], trial_function_spaces_[i]);
455 domain.insert_restriction(test_space_, test_function_space_);
457 using signature = test(decltype(serac::type<args>(trial_spaces))...);
458 integrals_.push_back(
459 MakeInteriorFaceIntegral<signature, Q, dim>(domain, integrand, std::vector<uint32_t>{args...}));
471 template <
int... args,
typename lambda,
typename qpt_data_type = Nothing>
472 void AddAreaIntegral(DependsOn<args...> which_args,
const lambda& integrand, Domain& domain,
473 std::shared_ptr<QuadratureData<qpt_data_type>> data =
NoQData)
475 AddDomainIntegral(Dimension<2>{}, which_args, integrand, domain, data);
487 template <
int... args,
typename lambda,
typename qpt_data_type =
Nothing>
491 AddDomainIntegral(
Dimension<3>{}, which_args, integrand, domain, data);
495 template <
int... args,
typename lambda>
498 AddBoundaryIntegral(
Dimension<2>{}, which_args, integrand, domain);
512 void ActionOfGradient(
const mfem::Vector& input_T, mfem::Vector& output_T, uint32_t which)
const
515 if (trial_function_spaces_[which].family == Family::L2) {
519 P_trial_[which]->Mult(input_T, input_ldof_values_[which]);
520 input_L_[which].SetVector(input_ldof_values_[which], 0);
521 updateFaceNbrData(trial_space_[which], trial_pargrid_functions_[which], input_ldof_values_[which]);
524 if (trial_space_[which]->GetVDim() == 1) {
525 input_L_[which].SetVector(trial_pargrid_functions_[which].FaceNbrData(), input_ldof_values_[which].Size());
527 if (trial_space_[which]->GetOrdering() == mfem::Ordering::byVDIM) {
528 int local_size = input_ldof_values_[which].Size();
529 appendFaceNbrData(trial_space_[which], trial_pargrid_functions_[which], local_size, input_L_[which]);
531 SLIC_ERROR_ROOT(
"Unsupported: L2 vector field ordered by nodes");
535 P_trial_[which]->Mult(input_T, input_L_[which]);
540 for (
auto& integral : integrals_) {
541 if (integral.DependsOn(which)) {
542 Domain& dom = integral.domain_;
545 input_E_buffer_[which].SetSize(
int(G_trial.
ESize()));
546 input_E_[which].Update(input_E_buffer_[which], G_trial.
bOffsets());
547 G_trial.
Gather(input_L_[which], input_E_[which]);
550 output_E_buffer_.SetSize(
int(G_test.
ESize()));
551 output_E_.Update(output_E_buffer_, G_test.
bOffsets());
553 integral.GradientMult(input_E_[which], output_E_, which);
561 if (test_function_space_.family == Family::L2) {
564 mfem::Vector output_ldof_values(P_test_->Height(), mem_type);
565 for (
int j = 0; j < output_ldof_values.Size(); ++j) {
566 output_ldof_values[j] = output_L_[j];
568 P_test_->MultTranspose(output_ldof_values, output_T);
570 P_test_->MultTranspose(output_L_, output_T);
586 template <uint32_t wrt,
typename... T>
589 const mfem::Vector* input_T[] = {&
static_cast<const mfem::Vector&
>(args)...};
592 for (uint32_t i = 0; i < num_trial_spaces; i++) {
593 if (trial_function_spaces_[i].family == Family::L2) {
597 P_trial_[i]->Mult(*input_T[i], input_ldof_values_[i]);
601 input_L_[i].SetVector(input_ldof_values_[i], 0);
604 updateFaceNbrData(trial_space_[i], trial_pargrid_functions_[i], input_ldof_values_[i]);
618 if (trial_space_[i]->GetVDim() == 1) {
620 input_L_[i].SetVector(trial_pargrid_functions_[i].FaceNbrData(), input_ldof_values_[i].Size());
622 if (trial_space_[i]->GetOrdering() == mfem::Ordering::byVDIM) {
625 int local_size = input_ldof_values_[i].Size();
626 appendFaceNbrData(trial_space_[i], trial_pargrid_functions_[i], local_size, input_L_[i]);
635 SLIC_ERROR_ROOT(
"Unsupported: L2 vector field ordered by nodes");
639 P_trial_[i]->Mult(*input_T[i], input_L_[i]);
645 for (
auto& integral : integrals_) {
646 Domain& dom = integral.domain_;
650 for (
auto i : integral.active_trial_spaces_) {
652 input_E_buffer_[i].SetSize(
int(G_trial.
ESize()));
653 input_E_[i].Update(input_E_buffer_[i], G_trial.
bOffsets());
654 G_trial.
Gather(input_L_[i], input_E_[i]);
657 output_E_buffer_.SetSize(
int(G_test.
ESize()));
658 output_E_.Update(output_E_buffer_, G_test.
bOffsets());
659 integral.Mult(t, input_E_, output_E_, wrt, update_qdata_);
666 if (test_function_space_.family == Family::L2) {
669 mfem::Vector output_ldof_values(P_test_->Height(), mem_type);
670 for (
int j = 0; j < output_ldof_values.Size(); ++j) {
671 output_ldof_values[j] = output_L_[j];
673 P_test_->MultTranspose(output_ldof_values, output_T_);
675 P_test_->MultTranspose(output_L_, output_T_);
678 if constexpr (wrt != NO_DIFFERENTIATION) {
685 return {output_T_, grad_[wrt]};
688 if constexpr (wrt == NO_DIFFERENTIATION) {
699 template <
typename... T>
700 auto operator()(
double t,
const T&... args)
702 constexpr
int num_differentiated_arguments = (std::is_same_v<T, differentiate_wrt_this> + ...);
703 static_assert(num_differentiated_arguments <= 1,
704 "Error: Functional::operator() can only differentiate w.r.t. 1 argument a time");
705 static_assert(
sizeof...(T) == num_trial_spaces,
706 "Error: Functional::operator() must take exactly as many arguments as trial spaces");
710 return (*
this)(DifferentiateWRT<i>{}, t, args...);
721 void updateQdata(
bool update_flag) { update_qdata_ = update_flag; }
732 class Gradient :
public mfem::Operator {
738 Gradient(Functional<test(trials...), exec>& f, uint32_t which = 0)
739 : mfem::Operator(f.test_space_->GetTrueVSize(), f.trial_space_[which]->GetTrueVSize()),
741 which_argument(which),
742 test_space_(f.test_space_),
743 trial_space_(f.trial_space_[which]),
744 df_(f.test_space_->GetTrueVSize())
754 virtual void Mult(
const mfem::Vector& dx, mfem::Vector& df)
const override
756 form_.ActionOfGradient(dx, df, which_argument);
760 mfem::Vector& operator()(
const mfem::Vector& dx)
762 form_.ActionOfGradient(dx, df_, which_argument);
766 void initialize_sparsity_pattern()
768 using row_col = std::tuple<int, int>;
770 std::set<row_col> nonzero_entries;
772 for (
auto& integral : form_.integrals_) {
773 if (integral.DependsOn(which_argument)) {
774 Domain& dom = integral.domain_;
775 const auto& G_test = dom.get_restriction(form_.test_function_space_);
776 const auto& G_trial = dom.get_restriction(form_.trial_function_spaces_[which_argument]);
777 for (
const auto& [geom, test_restriction] : G_test.
restrictions) {
778 const auto& trial_restriction = G_trial.
restrictions.at(geom);
781 std::vector<int> test_vdofs(test_restriction.nodes_per_elem * test_restriction.components);
782 std::vector<int> trial_vdofs(trial_restriction.nodes_per_elem * trial_restriction.components);
784 auto num_elements =
static_cast<uint32_t
>(test_restriction.num_elements);
785 for (uint32_t e = 0; e < num_elements; e++) {
786 for (uint32_t i = 0; i < test_restriction.nodes_per_elem; i++) {
787 auto test_dof = test_restriction.dof_info(e, i);
788 for (uint32_t j = 0; j < test_restriction.components; j++) {
789 test_vdofs[i * test_restriction.components + j] = int(test_restriction.GetVDof(test_dof, j).index());
793 for (uint32_t i = 0; i < trial_restriction.nodes_per_elem; i++) {
794 auto trial_dof = trial_restriction.dof_info(e, i);
795 for (uint32_t j = 0; j < trial_restriction.components; j++) {
796 trial_vdofs[i * trial_restriction.components + j] =
797 int(trial_restriction.GetVDof(trial_dof, j).index());
801 for (
int row : test_vdofs) {
802 for (
int col : trial_vdofs) {
803 nonzero_entries.insert({row, col});
811 uint64_t nnz = nonzero_entries.size();
812 int nrows = form_.output_L_.Size();
814 row_ptr.resize(uint32_t(nrows + 1));
819 for (
auto [row, col] : nonzero_entries) {
820 col_ind[uint32_t(nz)] = col;
821 for (
int i = last_row + 1; i <= row; i++) {
822 row_ptr[uint32_t(i)] = nz;
827 for (
int i = last_row + 1; i <= nrows; i++) {
828 row_ptr[uint32_t(i)] = nz;
832 uint64_t max_buffer_size()
834 uint64_t max_entries = 0;
835 for (
auto& integral : form_.integrals_) {
836 if (integral.DependsOn(which_argument)) {
837 Domain& dom = integral.domain_;
838 const auto& G_test = dom.get_restriction(form_.test_function_space_);
839 const auto& G_trial = dom.get_restriction(form_.trial_function_spaces_[which_argument]);
840 for (
const auto& [geom, test_restriction] : G_test.
restrictions) {
841 const auto& trial_restriction = G_trial.
restrictions.at(geom);
842 uint64_t nrows_per_element = test_restriction.nodes_per_elem * test_restriction.components;
843 uint64_t ncols_per_element = trial_restriction.nodes_per_elem * trial_restriction.components;
844 uint64_t entries_per_element = nrows_per_element * ncols_per_element;
845 uint64_t entries_needed = test_restriction.num_elements * entries_per_element;
846 max_entries =
std::max(entries_needed, max_entries);
853 std::unique_ptr<mfem::HypreParMatrix> assemble()
855 if (row_ptr.empty()) {
856 initialize_sparsity_pattern();
861 constexpr
bool sparse_matrix_frees_graph_ptrs =
false;
862 constexpr
bool sparse_matrix_frees_values_ptr =
false;
863 constexpr
bool col_ind_is_sorted =
true;
867 std::vector<int> col_ind_copy = col_ind;
869 int nnz = row_ptr.back();
870 std::vector<double> values(uint32_t(nnz), 0.0);
871 auto A_local = mfem::SparseMatrix(row_ptr.data(), col_ind_copy.data(), values.data(), form_.output_L_.Size(),
872 form_.input_L_[which_argument].Size(), sparse_matrix_frees_graph_ptrs,
873 sparse_matrix_frees_values_ptr, col_ind_is_sorted);
875 std::vector<double> K_elem_buffer(max_buffer_size());
877 for (
auto& integral : form_.integrals_) {
879 if (integral.functional_to_integral_index_.count(which_argument) > 0) {
880 Domain& dom = integral.domain_;
882 uint32_t
id = integral.functional_to_integral_index_.at(which_argument);
883 const auto& G_test = dom.get_restriction(form_.test_function_space_);
884 const auto& G_trial = dom.get_restriction(form_.trial_function_spaces_[which_argument]);
885 for (
const auto& [geom, calculate_element_matrices_func] : integral.element_gradient_[
id]) {
886 const auto& test_restriction = G_test.
restrictions.at(geom);
887 const auto& trial_restriction = G_trial.
restrictions.at(geom);
890 CPUArrayView<double, 3> K_e(K_elem_buffer.data(), test_restriction.num_elements,
891 trial_restriction.nodes_per_elem * trial_restriction.components,
892 test_restriction.nodes_per_elem * test_restriction.components);
893 detail::zero_out(K_e);
896 calculate_element_matrices_func(K_e);
898 const std::vector<int>& element_ids = integral.domain_.get(geom);
900 uint32_t rows_per_elem = uint32_t(test_restriction.nodes_per_elem * test_restriction.components);
901 uint32_t cols_per_elem = uint32_t(trial_restriction.nodes_per_elem * trial_restriction.components);
903 std::vector<DoF> test_vdofs(rows_per_elem);
904 std::vector<DoF> trial_vdofs(cols_per_elem);
906 for (uint32_t e = 0; e < element_ids.size(); e++) {
907 test_restriction.GetElementVDofs(
int(e), test_vdofs);
908 trial_restriction.GetElementVDofs(
int(e), trial_vdofs);
910 for (uint32_t i = 0; i < cols_per_elem; i++) {
911 int col = int(trial_vdofs[i].index());
913 for (uint32_t j = 0; j < rows_per_elem; j++) {
914 int row = int(test_vdofs[j].index());
915 A_local.SearchRow(row, col) += K_e(e, i, j);
923 auto* R = form_.test_space_->Dof_TrueDof_Matrix();
945 mfem::HypreParMatrix* A_hypre;
946 if (
dynamic_cast<const mfem::L2_FECollection*
>(test_space_->FEColl()) ||
947 dynamic_cast<const mfem::L2_FECollection*
>(trial_space_->FEColl())) {
948 int lrows = test_space_->GetVSize();
949 int lcols = trial_space_->GetVSize();
950 HYPRE_BigInt col_offset = trial_space_->GetMyDofOffset();
951 mfem::Array<HYPRE_BigInt> glob_J(A_local.NumNonZeroElems());
953 int* J = A_local.GetJ();
954 for (
int j = 0; j < glob_J.Size(); ++j) {
956 glob_J[j] = J[j] + col_offset;
958 glob_J[j] = form_.face_nbr_glob_vdof_maps_[which_argument][J[j] - lcols];
961 A_hypre =
new mfem::HypreParMatrix(test_space_->GetComm(), lrows, test_space_->GlobalVSize(),
962 trial_space_->GlobalVSize(), A_local.GetI(), glob_J, A_local.GetData(),
963 test_space_->GetDofOffsets(), trial_space_->GetDofOffsets());
967 new mfem::HypreParMatrix(test_space_->GetComm(), test_space_->GlobalVSize(), trial_space_->GlobalVSize(),
968 test_space_->GetDofOffsets(), trial_space_->GetDofOffsets(), &A_local);
971 auto* P = trial_space_->Dof_TrueDof_Matrix();
973 std::unique_ptr<mfem::HypreParMatrix> A(mfem::RAP(R, A_hypre, P));
980 friend auto assemble(Gradient& g) {
return g.assemble(); }
984 Functional<test(trials...), exec>& form_;
986 std::vector<int> row_ptr;
987 std::vector<int> col_ind;
999 uint32_t which_argument;
1002 const mfem::ParFiniteElementSpace* test_space_;
1005 const mfem::ParFiniteElementSpace* trial_space_;
1012 const mfem::ParFiniteElementSpace* test_space_;
1015 std::array<const mfem::ParFiniteElementSpace*, num_trial_spaces> trial_space_;
1017 std::array<FunctionSpace, num_trial_spaces> trial_function_spaces_;
1018 FunctionSpace test_function_space_;
1021 std::array<mfem::Array<HYPRE_BigInt>, num_trial_spaces> face_nbr_glob_vdof_maps_;
1024 mutable std::array<mfem::ParGridFunction, num_trial_spaces> trial_pargrid_functions_;
1027 mutable std::array<mfem::Vector, num_trial_spaces> input_ldof_values_;
1033 const mfem::Operator* P_trial_[num_trial_spaces];
1036 mutable mfem::Vector input_L_[num_trial_spaces];
1038 mutable std::vector<mfem::Vector> input_E_buffer_;
1039 mutable std::vector<mfem::BlockVector> input_E_;
1041 mutable std::vector<Integral> integrals_;
1043 mutable mfem::Vector output_E_buffer_;
1044 mutable mfem::BlockVector output_E_;
1047 mutable mfem::Vector output_L_;
1049 const mfem::Operator* P_test_;
1052 mutable mfem::Vector output_T_;
1055 mutable std::vector<Gradient> grad_;
1057 const mfem::MemoryType mem_type;
many of the functions in this file amount to extracting element indices from an mesh_t like
This file contains helper traits and enumerations for classifying finite elements.
a specialization of serac::Functional for quantities of interest
This file contains the all the necessary functions and macros required for logging as well as a helpe...
Accelerator functionality.
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.
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...
void check_for_missing_nodal_gridfunc(const mfem::Mesh &mesh)
function for verifying that the mesh has been fully initialized
constexpr uint32_t index_of_differentiation()
given a list of types, this function returns the index that corresponds to the type dual_vector.
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...
tuple(T...) -> tuple< T... >
Class template argument deduction rule for tuples.
constexpr SERAC_HOST_DEVICE auto type(const tuple< T... > &values)
a function intended to be used for extracting the ith type from a tuple.
ExecutionSpace
enum used for signalling whether or not to perform certain calculations on the CPU or GPU
void check_for_unsupported_elements(const mfem::Mesh &mesh)
function for verifying that there are no unsupported element types in the mesh
void rearrangeFaceNbrDofGlobalIndex(const mfem::ParFiniteElementSpace *trial_space, mfem::Array< HYPRE_BigInt > &face_nbr_glob_vdof_map)
helper functional to reorder the face_nbr_glob_dof_map for L2 space to byVDIM
std::pair< std::unique_ptr< mfem::ParFiniteElementSpace >, std::unique_ptr< mfem::FiniteElementCollection > > generateParFiniteElementSpace(mfem::ParMesh *mesh)
create an mfem::ParFiniteElementSpace from one of serac's tag types: H1, Hcurl, L2
#define SERAC_MARK_FUNCTION
Definitions of quadrature rules for quads and hexes.
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
void ScatterAdd(const mfem::BlockVector &E_block_vector, mfem::Vector &L_vector) const
"E->L" in mfem parlance, each element scatter-adds its local vector into the appropriate place in the...
uint64_t ESize() const
the size of the "E-vector" associated with this restriction operator
std::map< mfem::Geometry::Type, ElementRestriction > restrictions
the individual ElementRestriction operators for each element geometry
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.
a class for representing a geometric region that can be used for integration
Compile-time alias for index of differentiation.
these classes are a little confusing. These two special types represent the similar (but different) c...
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 ...
Implementation of the tensor class used by Functional.