xref: /libCEED/examples/deal.II/bps.h (revision 0850f9959daaba99c47e08ca6af4589f67130272)
18c81f8b0SPeter Munch // ---------------------------------------------------------------------
28c81f8b0SPeter Munch //
38c81f8b0SPeter Munch // Copyright (C) 2023 by the deal.II authors
48c81f8b0SPeter Munch //
58c81f8b0SPeter Munch // This file is part of the deal.II library.
68c81f8b0SPeter Munch //
78c81f8b0SPeter Munch // The deal.II library is free software; you can use it, redistribute
88c81f8b0SPeter Munch // it, and/or modify it under the terms of the GNU Lesser General
98c81f8b0SPeter Munch // Public License as published by the Free Software Foundation; either
108c81f8b0SPeter Munch // version 2.1 of the License, or (at your option) any later version.
118c81f8b0SPeter Munch // The full text of the license can be found in the file LICENSE.md at
128c81f8b0SPeter Munch // the top level directory of deal.II.
138c81f8b0SPeter Munch //
148c81f8b0SPeter Munch //  Authors: Peter Munch, Martin Kronbichler
158c81f8b0SPeter Munch //
168c81f8b0SPeter Munch // ---------------------------------------------------------------------
178c81f8b0SPeter Munch 
188c81f8b0SPeter Munch // deal.II includes
198c81f8b0SPeter Munch #include <deal.II/dofs/dof_tools.h>
208c81f8b0SPeter Munch 
218c81f8b0SPeter Munch #include <deal.II/fe/mapping.h>
228c81f8b0SPeter Munch 
238c81f8b0SPeter Munch #include <deal.II/lac/la_parallel_vector.h>
248c81f8b0SPeter Munch 
258c81f8b0SPeter Munch #include <deal.II/matrix_free/fe_evaluation.h>
268c81f8b0SPeter Munch #include <deal.II/matrix_free/matrix_free.h>
278c81f8b0SPeter Munch #include <deal.II/matrix_free/tools.h>
288c81f8b0SPeter Munch 
298c81f8b0SPeter Munch // libCEED includes
302097acd5SJeremy L Thompson #include <ceed.h>
312097acd5SJeremy L Thompson #include <ceed/backend.h>
328c81f8b0SPeter Munch 
338c81f8b0SPeter Munch // QFunction source
348c81f8b0SPeter Munch #include "bps-qfunctions.h"
358c81f8b0SPeter Munch 
368c81f8b0SPeter Munch using namespace dealii;
378c81f8b0SPeter Munch 
388c81f8b0SPeter Munch /**
398c81f8b0SPeter Munch  * BP types. For more details, see https://ceed.exascaleproject.org/bps/.
408c81f8b0SPeter Munch  */
418c81f8b0SPeter Munch enum class BPType : unsigned int
428c81f8b0SPeter Munch {
438c81f8b0SPeter Munch   BP1,
448c81f8b0SPeter Munch   BP2,
458c81f8b0SPeter Munch   BP3,
468c81f8b0SPeter Munch   BP4,
478c81f8b0SPeter Munch   BP5,
488c81f8b0SPeter Munch   BP6
498c81f8b0SPeter Munch };
508c81f8b0SPeter Munch 
518c81f8b0SPeter Munch 
528c81f8b0SPeter Munch 
538c81f8b0SPeter Munch /**
548c81f8b0SPeter Munch  * Struct storing relevant information regarding each BP.
558c81f8b0SPeter Munch  */
568c81f8b0SPeter Munch struct BPInfo
578c81f8b0SPeter Munch {
588c81f8b0SPeter Munch   BPInfo(const BPType type, const int dim, const int fe_degree)
598c81f8b0SPeter Munch     : type(type)
608c81f8b0SPeter Munch     , dim(dim)
618c81f8b0SPeter Munch     , fe_degree(fe_degree)
628c81f8b0SPeter Munch   {
638c81f8b0SPeter Munch     if (type == BPType::BP1)
648c81f8b0SPeter Munch       type_string = "BP1";
658c81f8b0SPeter Munch     else if (type == BPType::BP2)
668c81f8b0SPeter Munch       type_string = "BP2";
678c81f8b0SPeter Munch     else if (type == BPType::BP3)
688c81f8b0SPeter Munch       type_string = "BP3";
698c81f8b0SPeter Munch     else if (type == BPType::BP4)
708c81f8b0SPeter Munch       type_string = "BP4";
718c81f8b0SPeter Munch     else if (type == BPType::BP5)
728c81f8b0SPeter Munch       type_string = "BP5";
738c81f8b0SPeter Munch     else if (type == BPType::BP6)
748c81f8b0SPeter Munch       type_string = "BP6";
758c81f8b0SPeter Munch 
768c81f8b0SPeter Munch     this->n_q_points_1d = (type <= BPType::BP4) ? (fe_degree + 2) : (fe_degree + 1);
778c81f8b0SPeter Munch 
788c81f8b0SPeter Munch     this->n_components =
798c81f8b0SPeter Munch       (type == BPType::BP1 || type == BPType::BP3 || type == BPType::BP5) ? 1 : dim;
808c81f8b0SPeter Munch   }
818c81f8b0SPeter Munch 
828c81f8b0SPeter Munch 
838c81f8b0SPeter Munch   BPType       type;
848c81f8b0SPeter Munch   std::string  type_string;
858c81f8b0SPeter Munch   unsigned int dim;
868c81f8b0SPeter Munch   unsigned int fe_degree;
878c81f8b0SPeter Munch   unsigned int n_q_points_1d;
888c81f8b0SPeter Munch   unsigned int n_components;
898c81f8b0SPeter Munch };
908c81f8b0SPeter Munch 
918c81f8b0SPeter Munch 
928c81f8b0SPeter Munch 
938c81f8b0SPeter Munch /**
948c81f8b0SPeter Munch  * Base class of operators.
958c81f8b0SPeter Munch  */
968c81f8b0SPeter Munch template <typename Number>
978c81f8b0SPeter Munch class OperatorBase
988c81f8b0SPeter Munch {
998c81f8b0SPeter Munch public:
1008c81f8b0SPeter Munch   /**
1018c81f8b0SPeter Munch    * deal.II vector type
1028c81f8b0SPeter Munch    */
1038c81f8b0SPeter Munch   using VectorType = LinearAlgebra::distributed::Vector<Number>;
1048c81f8b0SPeter Munch 
1058c81f8b0SPeter Munch   /**
1068c81f8b0SPeter Munch    * Initialize vector.
1078c81f8b0SPeter Munch    */
1088c81f8b0SPeter Munch   virtual void
1098c81f8b0SPeter Munch   reinit() = 0;
1108c81f8b0SPeter Munch 
1118c81f8b0SPeter Munch   /**
1128c81f8b0SPeter Munch    * Perform matrix-vector product
1138c81f8b0SPeter Munch    */
1148c81f8b0SPeter Munch   virtual void
1158c81f8b0SPeter Munch   vmult(VectorType &dst, const VectorType &src) const = 0;
1168c81f8b0SPeter Munch 
1178c81f8b0SPeter Munch   /**
1188c81f8b0SPeter Munch    * Initialize vector.
1198c81f8b0SPeter Munch    */
1208c81f8b0SPeter Munch   virtual void
1218c81f8b0SPeter Munch   initialize_dof_vector(VectorType &vec) const = 0;
1228c81f8b0SPeter Munch 
1238c81f8b0SPeter Munch   /**
1248c81f8b0SPeter Munch    * Compute inverse of diagonal.
1258c81f8b0SPeter Munch    */
1268c81f8b0SPeter Munch   virtual void
1278c81f8b0SPeter Munch   compute_inverse_diagonal(VectorType &diagonal) const = 0;
1288c81f8b0SPeter Munch };
1298c81f8b0SPeter Munch 
1308c81f8b0SPeter Munch 
1318c81f8b0SPeter Munch /**
1328c81f8b0SPeter Munch  * Operator implementation using libCEED.
1338c81f8b0SPeter Munch  */
1348c81f8b0SPeter Munch template <int dim, typename Number>
1358c81f8b0SPeter Munch class OperatorCeed : public OperatorBase<Number>
1368c81f8b0SPeter Munch {
1378c81f8b0SPeter Munch public:
1388c81f8b0SPeter Munch   using VectorType = typename OperatorBase<Number>::VectorType;
1398c81f8b0SPeter Munch 
1408c81f8b0SPeter Munch   /**
1418c81f8b0SPeter Munch    * Constructor.
1428c81f8b0SPeter Munch    */
1438c81f8b0SPeter Munch   OperatorCeed(const Mapping<dim>              &mapping,
1448c81f8b0SPeter Munch                const DoFHandler<dim>           &dof_handler,
1458c81f8b0SPeter Munch                const AffineConstraints<Number> &constraints,
1468c81f8b0SPeter Munch                const Quadrature<dim>           &quadrature,
1478c81f8b0SPeter Munch                const BPType                    &bp,
1488c81f8b0SPeter Munch                const std::string               &resource)
1498c81f8b0SPeter Munch     : mapping(mapping)
1508c81f8b0SPeter Munch     , dof_handler(dof_handler)
1518c81f8b0SPeter Munch     , constraints(constraints)
1528c81f8b0SPeter Munch     , quadrature(quadrature)
1538c81f8b0SPeter Munch     , bp(bp)
1548c81f8b0SPeter Munch     , resource(resource)
1558c81f8b0SPeter Munch   {
1568c81f8b0SPeter Munch     reinit();
1578c81f8b0SPeter Munch   }
1588c81f8b0SPeter Munch 
1598c81f8b0SPeter Munch   /**
1608c81f8b0SPeter Munch    * Destructor.
1618c81f8b0SPeter Munch    */
1628c81f8b0SPeter Munch   ~OperatorCeed()
1638c81f8b0SPeter Munch   {
1642097acd5SJeremy L Thompson     CeedVectorDestroy(&src_ceed);
1652097acd5SJeremy L Thompson     CeedVectorDestroy(&dst_ceed);
16670dc8078SJeremy L Thompson     CeedOperatorDestroy(&op_apply);
1678c81f8b0SPeter Munch     CeedDestroy(&ceed);
1688c81f8b0SPeter Munch   }
1698c81f8b0SPeter Munch 
1708c81f8b0SPeter Munch   /**
1718c81f8b0SPeter Munch    * Initialized internal data structures, particularly, libCEED.
1728c81f8b0SPeter Munch    */
1738c81f8b0SPeter Munch   void
1748c81f8b0SPeter Munch   reinit() override
1758c81f8b0SPeter Munch   {
1768efac696SJeremy L Thompson     CeedVector           q_data;
1778efac696SJeremy L Thompson     CeedBasis            sol_basis;
1788efac696SJeremy L Thompson     CeedElemRestriction  sol_restriction;
1798efac696SJeremy L Thompson     CeedElemRestriction  q_data_restriction;
1808efac696SJeremy L Thompson     BuildContext         build_ctx_data;
1818efac696SJeremy L Thompson     CeedQFunctionContext build_ctx;
1828efac696SJeremy L Thompson     CeedQFunction        qf_apply;
1838efac696SJeremy L Thompson 
1848c81f8b0SPeter Munch     const auto &tria = dof_handler.get_triangulation();
1858c81f8b0SPeter Munch     const auto &fe   = dof_handler.get_fe();
1868c81f8b0SPeter Munch 
1878c81f8b0SPeter Munch     const auto n_components = fe.n_components();
1888c81f8b0SPeter Munch 
1898c81f8b0SPeter Munch     if (bp == BPType::BP1 || bp == BPType::BP3 || bp == BPType::BP5)
1908c81f8b0SPeter Munch       {
1918c81f8b0SPeter Munch         AssertThrow(n_components == 1, ExcInternalError());
1928c81f8b0SPeter Munch       }
1938c81f8b0SPeter Munch     else
1948c81f8b0SPeter Munch       {
1958c81f8b0SPeter Munch         AssertThrow(n_components == dim, ExcInternalError());
1968c81f8b0SPeter Munch       }
1978c81f8b0SPeter Munch 
1988c81f8b0SPeter Munch     // 1) create CEED instance -> "MatrixFree"
1998c81f8b0SPeter Munch     const char *ceed_spec = resource.c_str();
2008c81f8b0SPeter Munch     CeedInit(ceed_spec, &ceed);
2018c81f8b0SPeter Munch 
2028c81f8b0SPeter Munch     // 2) create shape functions -> "ShapeInfo"
2038c81f8b0SPeter Munch     const unsigned int fe_degree  = fe.tensor_degree();
2048c81f8b0SPeter Munch     const unsigned int n_q_points = quadrature.get_tensor_basis()[0].size();
2058c81f8b0SPeter Munch     CeedBasisCreateTensorH1Lagrange(
2068c81f8b0SPeter Munch       ceed, dim, n_components, fe_degree + 1, n_q_points, CEED_GAUSS, &sol_basis);
2078c81f8b0SPeter Munch 
2088c81f8b0SPeter Munch     // 3) create restriction matrix -> DoFInfo
2098c81f8b0SPeter Munch     unsigned int n_local_active_cells = 0;
2108c81f8b0SPeter Munch 
2118c81f8b0SPeter Munch     for (const auto &cell : dof_handler.active_cell_iterators())
2128c81f8b0SPeter Munch       if (cell->is_locally_owned())
2138c81f8b0SPeter Munch         n_local_active_cells++;
2148c81f8b0SPeter Munch 
2158c81f8b0SPeter Munch     partitioner =
2168c81f8b0SPeter Munch       std::make_shared<Utilities::MPI::Partitioner>(dof_handler.locally_owned_dofs(),
2178c81f8b0SPeter Munch                                                     DoFTools::extract_locally_active_dofs(
2188c81f8b0SPeter Munch                                                       dof_handler),
2198c81f8b0SPeter Munch                                                     dof_handler.get_communicator());
2208c81f8b0SPeter Munch 
2218c81f8b0SPeter Munch     std::vector<CeedInt> indices;
2228c81f8b0SPeter Munch     indices.reserve(n_local_active_cells * fe.n_dofs_per_cell() / n_components);
2238c81f8b0SPeter Munch 
2248c81f8b0SPeter Munch     const auto dof_mapping = FETools::lexicographic_to_hierarchic_numbering<dim>(fe_degree);
2258c81f8b0SPeter Munch 
2268c81f8b0SPeter Munch     std::vector<types::global_dof_index> local_indices(fe.n_dofs_per_cell());
2278c81f8b0SPeter Munch 
2288c81f8b0SPeter Munch     for (const auto &cell : dof_handler.active_cell_iterators())
2298c81f8b0SPeter Munch       if (cell->is_locally_owned())
2308c81f8b0SPeter Munch         {
2318c81f8b0SPeter Munch           cell->get_dof_indices(local_indices);
2328c81f8b0SPeter Munch 
2338c81f8b0SPeter Munch           for (const auto i : dof_mapping)
2348c81f8b0SPeter Munch             indices.emplace_back(
2358c81f8b0SPeter Munch               partitioner->global_to_local(local_indices[fe.component_to_system_index(0, i)]) /
2368c81f8b0SPeter Munch               n_components);
2378c81f8b0SPeter Munch         }
2388c81f8b0SPeter Munch 
2398c81f8b0SPeter Munch     CeedElemRestrictionCreate(ceed,
2408c81f8b0SPeter Munch                               n_local_active_cells,
2418c81f8b0SPeter Munch                               fe.n_dofs_per_cell() / n_components,
2428c81f8b0SPeter Munch                               n_components,
2438c81f8b0SPeter Munch                               std::max<unsigned int>(this->extended_local_size() / n_components, 1),
2448c81f8b0SPeter Munch                               this->extended_local_size(),
2458c81f8b0SPeter Munch                               CEED_MEM_HOST,
2468c81f8b0SPeter Munch                               CEED_COPY_VALUES,
2478c81f8b0SPeter Munch                               indices.data(),
2488c81f8b0SPeter Munch                               &sol_restriction);
2498c81f8b0SPeter Munch 
2508c81f8b0SPeter Munch     // 4) create mapping -> MappingInfo
2518c81f8b0SPeter Munch     const unsigned int n_components_metric = (bp <= BPType::BP2) ? 1 : (dim * (dim + 1) / 2);
2528c81f8b0SPeter Munch 
2538c81f8b0SPeter Munch     this->weights = compute_metric_data(ceed, mapping, tria, quadrature, bp);
2548c81f8b0SPeter Munch 
2558c81f8b0SPeter Munch     strides = {{1,
2568c81f8b0SPeter Munch                 static_cast<int>(quadrature.size()),
2578c81f8b0SPeter Munch                 static_cast<int>(quadrature.size() * n_components_metric)}};
2588c81f8b0SPeter Munch     CeedVectorCreate(ceed, weights.size(), &q_data);
2598c81f8b0SPeter Munch     CeedVectorSetArray(q_data, CEED_MEM_HOST, CEED_USE_POINTER, weights.data());
2608c81f8b0SPeter Munch     CeedElemRestrictionCreateStrided(ceed,
2618c81f8b0SPeter Munch                                      n_local_active_cells,
2628c81f8b0SPeter Munch                                      quadrature.size(),
2638c81f8b0SPeter Munch                                      n_components_metric,
2648c81f8b0SPeter Munch                                      weights.size(),
2658c81f8b0SPeter Munch                                      strides.data(),
2668c81f8b0SPeter Munch                                      &q_data_restriction);
2678c81f8b0SPeter Munch 
2688c81f8b0SPeter Munch     build_ctx_data.dim       = dim;
2698c81f8b0SPeter Munch     build_ctx_data.space_dim = dim;
2708c81f8b0SPeter Munch 
2718c81f8b0SPeter Munch     CeedQFunctionContextCreate(ceed, &build_ctx);
2728c81f8b0SPeter Munch     CeedQFunctionContextSetData(
2738efac696SJeremy L Thompson       build_ctx, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(build_ctx_data), &build_ctx_data);
2748c81f8b0SPeter Munch 
2758c81f8b0SPeter Munch     // 5) create q operation
2768c81f8b0SPeter Munch     if (bp == BPType::BP1)
2778c81f8b0SPeter Munch       CeedQFunctionCreateInterior(ceed, 1, f_apply_mass, f_apply_mass_loc, &qf_apply);
2788c81f8b0SPeter Munch     else if (bp == BPType::BP2)
2798c81f8b0SPeter Munch       CeedQFunctionCreateInterior(ceed, 1, f_apply_mass_vec, f_apply_mass_vec_loc, &qf_apply);
2808c81f8b0SPeter Munch     else if (bp == BPType::BP3 || bp == BPType::BP5)
2818c81f8b0SPeter Munch       CeedQFunctionCreateInterior(ceed, 1, f_apply_poisson, f_apply_poisson_loc, &qf_apply);
2828c81f8b0SPeter Munch     else if (bp == BPType::BP4 || bp == BPType::BP6)
2838c81f8b0SPeter Munch       CeedQFunctionCreateInterior(ceed, 1, f_apply_poisson_vec, f_apply_poisson_vec_loc, &qf_apply);
2848c81f8b0SPeter Munch     else
2858c81f8b0SPeter Munch       AssertThrow(false, ExcInternalError());
2868c81f8b0SPeter Munch 
2878c81f8b0SPeter Munch     if (bp <= BPType::BP2)
2888c81f8b0SPeter Munch       CeedQFunctionAddInput(qf_apply, "u", n_components, CEED_EVAL_INTERP);
2898c81f8b0SPeter Munch     else
2908c81f8b0SPeter Munch       CeedQFunctionAddInput(qf_apply, "u", dim * n_components, CEED_EVAL_GRAD);
2918c81f8b0SPeter Munch 
2928c81f8b0SPeter Munch     CeedQFunctionAddInput(qf_apply, "qdata", n_components_metric, CEED_EVAL_NONE);
2938c81f8b0SPeter Munch 
2948c81f8b0SPeter Munch     if (bp <= BPType::BP2)
2958c81f8b0SPeter Munch       CeedQFunctionAddOutput(qf_apply, "v", n_components, CEED_EVAL_INTERP);
2968c81f8b0SPeter Munch     else
2978c81f8b0SPeter Munch       CeedQFunctionAddOutput(qf_apply, "v", dim * n_components, CEED_EVAL_GRAD);
2988c81f8b0SPeter Munch 
2998c81f8b0SPeter Munch     CeedQFunctionSetContext(qf_apply, build_ctx);
3008c81f8b0SPeter Munch 
3018c81f8b0SPeter Munch     // 6) put everything together
3028c81f8b0SPeter Munch     CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_apply);
3038c81f8b0SPeter Munch 
3048c81f8b0SPeter Munch     CeedOperatorSetField(op_apply, "u", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
305fab7d8a4SJeremy L Thompson     CeedOperatorSetField(op_apply, "qdata", q_data_restriction, CEED_BASIS_NONE, q_data);
3068c81f8b0SPeter Munch     CeedOperatorSetField(op_apply, "v", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
3078efac696SJeremy L Thompson 
3082097acd5SJeremy L Thompson     // 7) libCEED vectors
3092097acd5SJeremy L Thompson     CeedElemRestrictionCreateVector(sol_restriction, &src_ceed, NULL);
3102097acd5SJeremy L Thompson     CeedElemRestrictionCreateVector(sol_restriction, &dst_ceed, NULL);
3112097acd5SJeremy L Thompson 
3122097acd5SJeremy L Thompson     // 8) cleanup
3138efac696SJeremy L Thompson     CeedVectorDestroy(&q_data);
3148efac696SJeremy L Thompson     CeedElemRestrictionDestroy(&q_data_restriction);
3158efac696SJeremy L Thompson     CeedElemRestrictionDestroy(&sol_restriction);
3168efac696SJeremy L Thompson     CeedBasisDestroy(&sol_basis);
3178efac696SJeremy L Thompson     CeedQFunctionContextDestroy(&build_ctx);
3188efac696SJeremy L Thompson     CeedQFunctionDestroy(&qf_apply);
3198c81f8b0SPeter Munch   }
3208c81f8b0SPeter Munch 
3218c81f8b0SPeter Munch   /**
3228c81f8b0SPeter Munch    * Perform matrix-vector product.
3238c81f8b0SPeter Munch    */
3248c81f8b0SPeter Munch   void
3258c81f8b0SPeter Munch   vmult(VectorType &dst, const VectorType &src) const override
3268c81f8b0SPeter Munch   {
3278c81f8b0SPeter Munch     // communicate: update ghost values
3288c81f8b0SPeter Munch     src.update_ghost_values();
3298c81f8b0SPeter Munch 
3308c81f8b0SPeter Munch     if (dof_handler.get_fe().n_components() == 1)
3318c81f8b0SPeter Munch       {
3322097acd5SJeremy L Thompson         // pass memory buffers to libCEED
3332097acd5SJeremy L Thompson         VectorTypeCeed x(src_ceed);
3342097acd5SJeremy L Thompson         VectorTypeCeed y(dst_ceed);
335*0850f995SJeremy L Thompson         x.import_array(src, CEED_MEM_HOST);
336*0850f995SJeremy L Thompson         y.import_array(dst, CEED_MEM_HOST);
3378c81f8b0SPeter Munch 
3388c81f8b0SPeter Munch         // apply operator
3392097acd5SJeremy L Thompson         CeedOperatorApply(op_apply, x(), y(), CEED_REQUEST_IMMEDIATE);
3402097acd5SJeremy L Thompson 
3412097acd5SJeremy L Thompson         // pull arrays back to deal.II
342*0850f995SJeremy L Thompson         x.sync_array();
343*0850f995SJeremy L Thompson         y.sync_array();
3448c81f8b0SPeter Munch       }
3458c81f8b0SPeter Munch     else // TODO: needed for multiple components
3468c81f8b0SPeter Munch       {
3478c81f8b0SPeter Munch         // allocate space for block vectors
3488c81f8b0SPeter Munch         src_tmp.reinit(this->extended_local_size(), true);
3498c81f8b0SPeter Munch         dst_tmp.reinit(this->extended_local_size(), true);
3508c81f8b0SPeter Munch 
3512097acd5SJeremy L Thompson         // copy to block vector
3522097acd5SJeremy L Thompson         copy_to_block_vector(src_tmp, src);
3538c81f8b0SPeter Munch 
3542097acd5SJeremy L Thompson         // pass memory buffers to libCEED
3552097acd5SJeremy L Thompson         VectorTypeCeed x(src_ceed);
3562097acd5SJeremy L Thompson         VectorTypeCeed y(dst_ceed);
357*0850f995SJeremy L Thompson         x.import_array(src_tmp, CEED_MEM_HOST);
358*0850f995SJeremy L Thompson         y.import_array(dst_tmp, CEED_MEM_HOST);
3598c81f8b0SPeter Munch 
3608c81f8b0SPeter Munch         // apply operator
3612097acd5SJeremy L Thompson         CeedOperatorApply(op_apply, x(), y(), CEED_REQUEST_IMMEDIATE);
3628c81f8b0SPeter Munch 
3632097acd5SJeremy L Thompson         // pull arrays back to deal.II
364*0850f995SJeremy L Thompson         x.sync_array();
365*0850f995SJeremy L Thompson         y.sync_array();
3662097acd5SJeremy L Thompson 
3672097acd5SJeremy L Thompson         // copy from block vector
3682097acd5SJeremy L Thompson         copy_from_block_vector(dst, dst_tmp);
3698c81f8b0SPeter Munch       }
3708c81f8b0SPeter Munch 
3718c81f8b0SPeter Munch     // communicate: compress
3728c81f8b0SPeter Munch     src.zero_out_ghost_values();
3738c81f8b0SPeter Munch     dst.compress(VectorOperation::add);
3748c81f8b0SPeter Munch 
3758c81f8b0SPeter Munch     // apply constraints: we assume homogeneous DBC
3768c81f8b0SPeter Munch     constraints.set_zero(dst);
3778c81f8b0SPeter Munch   }
3788c81f8b0SPeter Munch 
3798c81f8b0SPeter Munch   /**
3808c81f8b0SPeter Munch    * Initialized vector.
3818c81f8b0SPeter Munch    */
3828c81f8b0SPeter Munch   void
3838c81f8b0SPeter Munch   initialize_dof_vector(VectorType &vec) const override
3848c81f8b0SPeter Munch   {
3858c81f8b0SPeter Munch     vec.reinit(partitioner);
3868c81f8b0SPeter Munch   }
3878c81f8b0SPeter Munch 
3888c81f8b0SPeter Munch   /**
3898c81f8b0SPeter Munch    * Compute inverse of diagonal.
3908c81f8b0SPeter Munch    */
3918c81f8b0SPeter Munch   void
3928c81f8b0SPeter Munch   compute_inverse_diagonal(VectorType &diagonal) const override
3938c81f8b0SPeter Munch   {
3948c81f8b0SPeter Munch     this->initialize_dof_vector(diagonal);
3958c81f8b0SPeter Munch 
3962097acd5SJeremy L Thompson     // pass memory buffer to libCEED
3972097acd5SJeremy L Thompson     VectorTypeCeed y(dst_ceed);
398*0850f995SJeremy L Thompson     y.import_array(diagonal, CEED_MEM_HOST);
3998c81f8b0SPeter Munch 
4002097acd5SJeremy L Thompson     CeedOperatorLinearAssembleDiagonal(op_apply, y(), CEED_REQUEST_IMMEDIATE);
4012097acd5SJeremy L Thompson 
4022097acd5SJeremy L Thompson     // pull array back to deal.II
403*0850f995SJeremy L Thompson     y.sync_array();
4048c81f8b0SPeter Munch 
4058c81f8b0SPeter Munch     const unsigned int n_components = dof_handler.get_fe().n_components();
4068c81f8b0SPeter Munch 
4078c81f8b0SPeter Munch     if (n_components > 1) // TODO: needed for multiple components
4088c81f8b0SPeter Munch       {
4098c81f8b0SPeter Munch         VectorType tmp(diagonal);
4108c81f8b0SPeter Munch 
4118c81f8b0SPeter Munch         copy_from_block_vector(tmp, diagonal);
4128c81f8b0SPeter Munch 
4138c81f8b0SPeter Munch         std::swap(tmp, diagonal);
4148c81f8b0SPeter Munch       }
4158c81f8b0SPeter Munch 
4168c81f8b0SPeter Munch     diagonal.compress(VectorOperation::add);
4178c81f8b0SPeter Munch 
4188c81f8b0SPeter Munch     for (auto &i : diagonal)
4198c81f8b0SPeter Munch       i = (std::abs(i) > 1.0e-10) ? (1.0 / i) : 1.0;
4208c81f8b0SPeter Munch   }
4218c81f8b0SPeter Munch 
4228c81f8b0SPeter Munch private:
4238c81f8b0SPeter Munch   /**
4248c81f8b0SPeter Munch    * Wrapper around a deal.II vector to create a libCEED vector view.
4258c81f8b0SPeter Munch    */
4268c81f8b0SPeter Munch   class VectorTypeCeed
4278c81f8b0SPeter Munch   {
4288c81f8b0SPeter Munch   public:
4298c81f8b0SPeter Munch     /**
4308c81f8b0SPeter Munch      * Constructor.
4318c81f8b0SPeter Munch      */
4322097acd5SJeremy L Thompson     VectorTypeCeed(const CeedVector &vec_orig)
4338c81f8b0SPeter Munch     {
4342097acd5SJeremy L Thompson       vec_ceed = NULL;
4352097acd5SJeremy L Thompson       CeedVectorReferenceCopy(vec_orig, &vec_ceed);
4368c81f8b0SPeter Munch     }
4378c81f8b0SPeter Munch 
4388c81f8b0SPeter Munch     /**
4398c81f8b0SPeter Munch      * Return libCEED vector view.
4408c81f8b0SPeter Munch      */
4418c81f8b0SPeter Munch     CeedVector &
4428c81f8b0SPeter Munch     operator()()
4438c81f8b0SPeter Munch     {
4448c81f8b0SPeter Munch       return vec_ceed;
4458c81f8b0SPeter Munch     }
4468c81f8b0SPeter Munch 
4478c81f8b0SPeter Munch     /**
4482097acd5SJeremy L Thompson      * Set deal.II memory in libCEED vector.
4492097acd5SJeremy L Thompson      */
4502097acd5SJeremy L Thompson     void
451*0850f995SJeremy L Thompson     import_array(const VectorType &vec, const CeedMemType space)
4522097acd5SJeremy L Thompson     {
453*0850f995SJeremy L Thompson       mem_space = space;
454*0850f995SJeremy L Thompson       CeedVectorSetArray(vec_ceed, mem_space, CEED_USE_POINTER, vec.get_values());
4552097acd5SJeremy L Thompson     }
4562097acd5SJeremy L Thompson 
4572097acd5SJeremy L Thompson     /**
4588c81f8b0SPeter Munch      * Sync memory from device to host.
4598c81f8b0SPeter Munch      */
4608c81f8b0SPeter Munch     void
461*0850f995SJeremy L Thompson     sync_array()
4628c81f8b0SPeter Munch     {
463*0850f995SJeremy L Thompson       CeedVectorSyncArray(vec_ceed, mem_space);
4648c81f8b0SPeter Munch     }
4658c81f8b0SPeter Munch 
4668c81f8b0SPeter Munch     /**
4678c81f8b0SPeter Munch      * Destructor: destroy vector view.
4688c81f8b0SPeter Munch      */
4698c81f8b0SPeter Munch     ~VectorTypeCeed()
4708c81f8b0SPeter Munch     {
4712097acd5SJeremy L Thompson       bool has_array;
472*0850f995SJeremy L Thompson       CeedVectorHasBorrowedArrayOfType(vec_ceed, mem_space, &has_array);
4732097acd5SJeremy L Thompson       if (has_array)
4742097acd5SJeremy L Thompson         {
4758c81f8b0SPeter Munch           CeedScalar *ptr;
476*0850f995SJeremy L Thompson           CeedVectorTakeArray(vec_ceed, mem_space, &ptr);
4772097acd5SJeremy L Thompson         }
4788c81f8b0SPeter Munch       CeedVectorDestroy(&vec_ceed);
4798c81f8b0SPeter Munch     }
4808c81f8b0SPeter Munch 
4818c81f8b0SPeter Munch   private:
4828c81f8b0SPeter Munch     /**
4838c81f8b0SPeter Munch      * libCEED vector view.
4848c81f8b0SPeter Munch      */
485*0850f995SJeremy L Thompson     CeedMemType mem_space;
4868c81f8b0SPeter Munch     CeedVector  vec_ceed;
4878c81f8b0SPeter Munch   };
4888c81f8b0SPeter Munch 
4898c81f8b0SPeter Munch   /**
4908c81f8b0SPeter Munch    * Copy from block vector.
4918c81f8b0SPeter Munch    *
4928c81f8b0SPeter Munch    * @note Only needed for multiple components.
4938c81f8b0SPeter Munch    */
4948c81f8b0SPeter Munch   void
4958c81f8b0SPeter Munch   copy_from_block_vector(VectorType &dst, const VectorType &src) const
4968c81f8b0SPeter Munch   {
4978c81f8b0SPeter Munch     const unsigned int scalar_size = this->extended_local_size() / dim;
4988c81f8b0SPeter Munch 
4998c81f8b0SPeter Munch     for (unsigned int i = 0; i < scalar_size; ++i)
5008c81f8b0SPeter Munch       for (unsigned int j = 0; j < dim; ++j)
5018c81f8b0SPeter Munch         dst.get_values()[j + i * dim] = src.get_values()[j * scalar_size + i];
5028c81f8b0SPeter Munch   }
5038c81f8b0SPeter Munch 
5048c81f8b0SPeter Munch   /**
5058c81f8b0SPeter Munch    * Copy to block vector.
5068c81f8b0SPeter Munch    *
5078c81f8b0SPeter Munch    * @note Only needed for multiple components.
5088c81f8b0SPeter Munch    */
5098c81f8b0SPeter Munch   void
5108c81f8b0SPeter Munch   copy_to_block_vector(VectorType &dst, const VectorType &src) const
5118c81f8b0SPeter Munch   {
5128c81f8b0SPeter Munch     const unsigned int scalar_size = this->extended_local_size() / dim;
5138c81f8b0SPeter Munch 
5148c81f8b0SPeter Munch     for (unsigned int i = 0; i < scalar_size; ++i)
5158c81f8b0SPeter Munch       for (unsigned int j = 0; j < dim; ++j)
5168c81f8b0SPeter Munch         dst.get_values()[j * scalar_size + i] = src.get_values()[j + i * dim];
5178c81f8b0SPeter Munch   }
5188c81f8b0SPeter Munch 
5198c81f8b0SPeter Munch   /**
5208c81f8b0SPeter Munch    * Number of locally active DoFs.
5218c81f8b0SPeter Munch    */
5228c81f8b0SPeter Munch   unsigned int
5238c81f8b0SPeter Munch   extended_local_size() const
5248c81f8b0SPeter Munch   {
5258c81f8b0SPeter Munch     return partitioner->locally_owned_size() + partitioner->n_ghost_indices();
5268c81f8b0SPeter Munch   }
5278c81f8b0SPeter Munch 
5288c81f8b0SPeter Munch   /**
5298c81f8b0SPeter Munch    * Compute metric data: Jacobian, ...
5308c81f8b0SPeter Munch    */
5318c81f8b0SPeter Munch   static std::vector<double>
5328c81f8b0SPeter Munch   compute_metric_data(const Ceed               &ceed,
5338c81f8b0SPeter Munch                       const Mapping<dim>       &mapping,
5348c81f8b0SPeter Munch                       const Triangulation<dim> &tria,
5358c81f8b0SPeter Munch                       const Quadrature<dim>    &quadrature,
5368c81f8b0SPeter Munch                       const BPType              bp)
5378c81f8b0SPeter Munch   {
5388c81f8b0SPeter Munch     std::vector<double> weights;
5398c81f8b0SPeter Munch 
5408c81f8b0SPeter Munch     if (false)
5418c81f8b0SPeter Munch       {
5428c81f8b0SPeter Munch         FE_Nothing<dim> dummy_fe;
5438c81f8b0SPeter Munch         FEValues<dim>   fe_values(mapping, dummy_fe, quadrature, update_JxW_values);
5448c81f8b0SPeter Munch 
5458c81f8b0SPeter Munch         for (const auto &cell : tria.active_cell_iterators())
5468c81f8b0SPeter Munch           if (cell->is_locally_owned())
5478c81f8b0SPeter Munch             {
5488c81f8b0SPeter Munch               fe_values.reinit(cell);
5498c81f8b0SPeter Munch 
5508c81f8b0SPeter Munch               for (const auto q : fe_values.quadrature_point_indices())
5518c81f8b0SPeter Munch                 weights.emplace_back(fe_values.JxW(q));
5528c81f8b0SPeter Munch             }
5538c81f8b0SPeter Munch 
5548c81f8b0SPeter Munch         return weights;
5558c81f8b0SPeter Munch       }
5568c81f8b0SPeter Munch 
5578c81f8b0SPeter Munch     CeedBasis            geo_basis;
5588c81f8b0SPeter Munch     CeedVector           q_data;
5598c81f8b0SPeter Munch     CeedElemRestriction  q_data_restriction;
5608c81f8b0SPeter Munch     CeedVector           node_coords;
5618c81f8b0SPeter Munch     CeedElemRestriction  geo_restriction;
5628c81f8b0SPeter Munch     CeedQFunctionContext build_ctx;
5638c81f8b0SPeter Munch     CeedQFunction        qf_build;
5648c81f8b0SPeter Munch     CeedOperator         op_build;
5658c81f8b0SPeter Munch 
5668c81f8b0SPeter Munch     const unsigned int n_q_points = quadrature.get_tensor_basis()[0].size();
5678c81f8b0SPeter Munch 
5688c81f8b0SPeter Munch     const unsigned int n_components = (bp <= BPType::BP2) ? 1 : (dim * (dim + 1) / 2);
5698c81f8b0SPeter Munch 
5708c81f8b0SPeter Munch     const auto mapping_q = dynamic_cast<const MappingQ<dim> *>(&mapping);
5718c81f8b0SPeter Munch 
5728c81f8b0SPeter Munch     AssertThrow(mapping_q, ExcMessage("Wrong mapping!"));
5738c81f8b0SPeter Munch 
5748c81f8b0SPeter Munch     const unsigned int fe_degree = mapping_q->get_degree();
5758c81f8b0SPeter Munch 
5768c81f8b0SPeter Munch     CeedBasisCreateTensorH1Lagrange(
5778c81f8b0SPeter Munch       ceed, dim, dim, fe_degree + 1, n_q_points, CEED_GAUSS, &geo_basis);
5788c81f8b0SPeter Munch 
5798c81f8b0SPeter Munch     unsigned int n_local_active_cells = 0;
5808c81f8b0SPeter Munch 
5818c81f8b0SPeter Munch     for (const auto &cell : tria.active_cell_iterators())
5828c81f8b0SPeter Munch       if (cell->is_locally_owned())
5838c81f8b0SPeter Munch         n_local_active_cells++;
5848c81f8b0SPeter Munch 
5858c81f8b0SPeter Munch     std::vector<double>  geo_support_points;
5868c81f8b0SPeter Munch     std::vector<CeedInt> geo_indices;
5878c81f8b0SPeter Munch 
5888c81f8b0SPeter Munch     FE_Q<dim> geo_fe(fe_degree);
5898c81f8b0SPeter Munch 
5908c81f8b0SPeter Munch     DoFHandler<dim> geo_dof_handler(tria);
5918c81f8b0SPeter Munch     geo_dof_handler.distribute_dofs(geo_fe);
5928c81f8b0SPeter Munch 
5938c81f8b0SPeter Munch     const auto geo_partitioner =
5948c81f8b0SPeter Munch       std::make_shared<Utilities::MPI::Partitioner>(geo_dof_handler.locally_owned_dofs(),
5958c81f8b0SPeter Munch                                                     DoFTools::extract_locally_active_dofs(
5968c81f8b0SPeter Munch                                                       geo_dof_handler),
5978c81f8b0SPeter Munch                                                     geo_dof_handler.get_communicator());
5988c81f8b0SPeter Munch 
5998c81f8b0SPeter Munch     geo_indices.reserve(n_local_active_cells * geo_fe.n_dofs_per_cell());
6008c81f8b0SPeter Munch 
6018c81f8b0SPeter Munch     const auto dof_mapping = FETools::lexicographic_to_hierarchic_numbering<dim>(fe_degree);
6028c81f8b0SPeter Munch 
6038c81f8b0SPeter Munch     FEValues<dim> fe_values(mapping,
6048c81f8b0SPeter Munch                             geo_fe,
6058c81f8b0SPeter Munch                             geo_fe.get_unit_support_points(),
6068c81f8b0SPeter Munch                             update_quadrature_points);
6078c81f8b0SPeter Munch 
6088c81f8b0SPeter Munch     std::vector<types::global_dof_index> local_indices(geo_fe.n_dofs_per_cell());
6098c81f8b0SPeter Munch 
6108c81f8b0SPeter Munch     const unsigned int n_points =
6118c81f8b0SPeter Munch       geo_partitioner->locally_owned_size() + geo_partitioner->n_ghost_indices();
6128c81f8b0SPeter Munch 
6138c81f8b0SPeter Munch     geo_support_points.resize(dim * n_points);
6148c81f8b0SPeter Munch 
6158c81f8b0SPeter Munch     for (const auto &cell : geo_dof_handler.active_cell_iterators())
6168c81f8b0SPeter Munch       if (cell->is_locally_owned())
6178c81f8b0SPeter Munch         {
6188c81f8b0SPeter Munch           fe_values.reinit(cell);
6198c81f8b0SPeter Munch           cell->get_dof_indices(local_indices);
6208c81f8b0SPeter Munch 
6218c81f8b0SPeter Munch           for (const auto i : dof_mapping)
6228c81f8b0SPeter Munch             {
6238c81f8b0SPeter Munch               const auto index = geo_partitioner->global_to_local(local_indices[i]);
6248c81f8b0SPeter Munch               geo_indices.emplace_back(index);
6258c81f8b0SPeter Munch 
6268c81f8b0SPeter Munch               const auto point = fe_values.quadrature_point(i);
6278c81f8b0SPeter Munch 
6288c81f8b0SPeter Munch               for (unsigned int d = 0; d < dim; ++d)
6298c81f8b0SPeter Munch                 geo_support_points[index + d * n_points] = point[d];
6308c81f8b0SPeter Munch             }
6318c81f8b0SPeter Munch         }
6328c81f8b0SPeter Munch 
6338c81f8b0SPeter Munch     weights.resize(n_local_active_cells * quadrature.size() * n_components);
6348c81f8b0SPeter Munch 
6358c81f8b0SPeter Munch     CeedInt strides[3] = {1,
6368c81f8b0SPeter Munch                           static_cast<int>(quadrature.size()),
6378c81f8b0SPeter Munch                           static_cast<int>(quadrature.size() * n_components)};
6388c81f8b0SPeter Munch 
6398c81f8b0SPeter Munch     CeedVectorCreate(ceed, weights.size(), &q_data);
6408c81f8b0SPeter Munch     CeedVectorSetArray(q_data, CEED_MEM_HOST, CEED_USE_POINTER, weights.data());
6418c81f8b0SPeter Munch     CeedElemRestrictionCreateStrided(ceed,
6428c81f8b0SPeter Munch                                      n_local_active_cells,
6438c81f8b0SPeter Munch                                      quadrature.size(),
6448c81f8b0SPeter Munch                                      n_components,
6458c81f8b0SPeter Munch                                      weights.size(),
6468c81f8b0SPeter Munch                                      strides,
6478c81f8b0SPeter Munch                                      &q_data_restriction);
6488c81f8b0SPeter Munch 
6498c81f8b0SPeter Munch     CeedVectorCreate(ceed, geo_support_points.size(), &node_coords);
6508c81f8b0SPeter Munch     CeedVectorSetArray(node_coords, CEED_MEM_HOST, CEED_USE_POINTER, geo_support_points.data());
6518c81f8b0SPeter Munch 
6528c81f8b0SPeter Munch     CeedElemRestrictionCreate(ceed,
6538c81f8b0SPeter Munch                               n_local_active_cells,
6548c81f8b0SPeter Munch                               geo_fe.n_dofs_per_cell(),
6558c81f8b0SPeter Munch                               dim,
6568c81f8b0SPeter Munch                               std::max<unsigned int>(geo_support_points.size() / dim, 1),
6578c81f8b0SPeter Munch                               geo_support_points.size(),
6588c81f8b0SPeter Munch                               CEED_MEM_HOST,
6598c81f8b0SPeter Munch                               CEED_COPY_VALUES,
6608c81f8b0SPeter Munch                               geo_indices.data(),
6618c81f8b0SPeter Munch                               &geo_restriction);
6628c81f8b0SPeter Munch 
6638c81f8b0SPeter Munch     BuildContext build_ctx_data;
6648c81f8b0SPeter Munch     build_ctx_data.dim       = dim;
6658c81f8b0SPeter Munch     build_ctx_data.space_dim = dim;
6668c81f8b0SPeter Munch 
6678c81f8b0SPeter Munch     CeedQFunctionContextCreate(ceed, &build_ctx);
6688c81f8b0SPeter Munch     CeedQFunctionContextSetData(
6698c81f8b0SPeter Munch       build_ctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(build_ctx_data), &build_ctx_data);
6708c81f8b0SPeter Munch 
6718c81f8b0SPeter Munch     // 5) create q operation
6728c81f8b0SPeter Munch     if (bp <= BPType::BP2)
6738c81f8b0SPeter Munch       CeedQFunctionCreateInterior(ceed, 1, f_build_mass, f_build_mass_loc, &qf_build);
6748c81f8b0SPeter Munch     else
6758c81f8b0SPeter Munch       CeedQFunctionCreateInterior(ceed, 1, f_build_poisson, f_build_poisson_loc, &qf_build);
6768c81f8b0SPeter Munch 
6778c81f8b0SPeter Munch     CeedQFunctionAddInput(qf_build, "geo", dim * dim, CEED_EVAL_GRAD);
6788c81f8b0SPeter Munch     CeedQFunctionAddInput(qf_build, "weights", 1, CEED_EVAL_WEIGHT);
6798c81f8b0SPeter Munch     CeedQFunctionAddOutput(qf_build, "qdata", n_components, CEED_EVAL_NONE);
6808c81f8b0SPeter Munch     CeedQFunctionSetContext(qf_build, build_ctx);
6818c81f8b0SPeter Munch 
6828c81f8b0SPeter Munch     // 6) put everything together
6838c81f8b0SPeter Munch     CeedOperatorCreate(ceed, qf_build, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_build);
6848c81f8b0SPeter Munch     CeedOperatorSetField(op_build, "geo", geo_restriction, geo_basis, CEED_VECTOR_ACTIVE);
6858c81f8b0SPeter Munch     CeedOperatorSetField(
6868c81f8b0SPeter Munch       op_build, "weights", CEED_ELEMRESTRICTION_NONE, geo_basis, CEED_VECTOR_NONE);
6878c81f8b0SPeter Munch     CeedOperatorSetField(
688fab7d8a4SJeremy L Thompson       op_build, "qdata", q_data_restriction, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
6898c81f8b0SPeter Munch 
6908c81f8b0SPeter Munch     CeedOperatorApply(op_build, node_coords, q_data, CEED_REQUEST_IMMEDIATE);
6918c81f8b0SPeter Munch 
6928c81f8b0SPeter Munch     CeedVectorDestroy(&node_coords);
6938c81f8b0SPeter Munch     CeedVectorSyncArray(q_data, CEED_MEM_HOST);
6948c81f8b0SPeter Munch     CeedVectorDestroy(&q_data);
69570dc8078SJeremy L Thompson     CeedElemRestrictionDestroy(&geo_restriction);
69670dc8078SJeremy L Thompson     CeedElemRestrictionDestroy(&q_data_restriction);
6978c81f8b0SPeter Munch     CeedBasisDestroy(&geo_basis);
69870dc8078SJeremy L Thompson     CeedQFunctionContextDestroy(&build_ctx);
69970dc8078SJeremy L Thompson     CeedQFunctionDestroy(&qf_build);
70070dc8078SJeremy L Thompson     CeedOperatorDestroy(&op_build);
7018c81f8b0SPeter Munch 
7028c81f8b0SPeter Munch     return weights;
7038c81f8b0SPeter Munch   }
7048c81f8b0SPeter Munch 
7058c81f8b0SPeter Munch   /**
7068c81f8b0SPeter Munch    * Mapping object passed to the constructor.
7078c81f8b0SPeter Munch    */
7088c81f8b0SPeter Munch   const Mapping<dim> &mapping;
7098c81f8b0SPeter Munch 
7108c81f8b0SPeter Munch   /**
7118c81f8b0SPeter Munch    * DoFHandler object passed to the constructor.
7128c81f8b0SPeter Munch    */
7138c81f8b0SPeter Munch   const DoFHandler<dim> &dof_handler;
7148c81f8b0SPeter Munch 
7158c81f8b0SPeter Munch   /**
7168c81f8b0SPeter Munch    * Constraints object passed to the constructor.
7178c81f8b0SPeter Munch    */
7188c81f8b0SPeter Munch   const AffineConstraints<Number> &constraints;
7198c81f8b0SPeter Munch 
7208c81f8b0SPeter Munch   /**
7218c81f8b0SPeter Munch    * Quadrature rule object passed to the constructor.
7228c81f8b0SPeter Munch    */
7238c81f8b0SPeter Munch   const Quadrature<dim> &quadrature;
7248c81f8b0SPeter Munch 
7258c81f8b0SPeter Munch   /**
7268c81f8b0SPeter Munch    * Selected BP.
7278c81f8b0SPeter Munch    */
7288c81f8b0SPeter Munch   const BPType bp;
7298c81f8b0SPeter Munch 
7308c81f8b0SPeter Munch   /**
7318c81f8b0SPeter Munch    * Resource name.
7328c81f8b0SPeter Munch    */
7338c81f8b0SPeter Munch   const std::string resource;
7348c81f8b0SPeter Munch 
7358c81f8b0SPeter Munch   /**
7368c81f8b0SPeter Munch    * Partitioner for distributed vectors.
7378c81f8b0SPeter Munch    */
7388c81f8b0SPeter Munch   std::shared_ptr<Utilities::MPI::Partitioner> partitioner;
7398c81f8b0SPeter Munch 
7408c81f8b0SPeter Munch   /**
7418c81f8b0SPeter Munch    * libCEED data structures.
7428c81f8b0SPeter Munch    */
7438c81f8b0SPeter Munch   Ceed                   ceed;
7448c81f8b0SPeter Munch   std::vector<double>    weights;
7458c81f8b0SPeter Munch   std::array<CeedInt, 3> strides;
7462097acd5SJeremy L Thompson   CeedVector             src_ceed;
7472097acd5SJeremy L Thompson   CeedVector             dst_ceed;
7488c81f8b0SPeter Munch   CeedOperator           op_apply;
7498c81f8b0SPeter Munch 
7508c81f8b0SPeter Munch   /**
7518c81f8b0SPeter Munch    * Temporal (tempral) vectors.
7528c81f8b0SPeter Munch    *
7538c81f8b0SPeter Munch    * @note Only needed for multiple components.
7548c81f8b0SPeter Munch    */
7558c81f8b0SPeter Munch   mutable VectorType src_tmp;
7568c81f8b0SPeter Munch   mutable VectorType dst_tmp;
7578c81f8b0SPeter Munch };
7588c81f8b0SPeter Munch 
7598c81f8b0SPeter Munch 
7608c81f8b0SPeter Munch 
7618c81f8b0SPeter Munch template <int dim, typename Number>
7628c81f8b0SPeter Munch class OperatorDealii : public OperatorBase<Number>
7638c81f8b0SPeter Munch {
7648c81f8b0SPeter Munch public:
7658c81f8b0SPeter Munch   using VectorType = typename OperatorBase<Number>::VectorType;
7668c81f8b0SPeter Munch 
7678c81f8b0SPeter Munch   /**
7688c81f8b0SPeter Munch    * Constructor.
7698c81f8b0SPeter Munch    */
7708c81f8b0SPeter Munch   OperatorDealii(const Mapping<dim>              &mapping,
7718c81f8b0SPeter Munch                  const DoFHandler<dim>           &dof_handler,
7728c81f8b0SPeter Munch                  const AffineConstraints<Number> &constraints,
7738c81f8b0SPeter Munch                  const Quadrature<dim>           &quadrature,
7748c81f8b0SPeter Munch                  const BPType                    &bp)
7758c81f8b0SPeter Munch     : mapping(mapping)
7768c81f8b0SPeter Munch     , dof_handler(dof_handler)
7778c81f8b0SPeter Munch     , constraints(constraints)
7788c81f8b0SPeter Munch     , quadrature(quadrature)
7798c81f8b0SPeter Munch     , bp(bp)
7808c81f8b0SPeter Munch   {
7818c81f8b0SPeter Munch     reinit();
7828c81f8b0SPeter Munch   }
7838c81f8b0SPeter Munch 
7848c81f8b0SPeter Munch   /**
7858c81f8b0SPeter Munch    * Destructor.
7868c81f8b0SPeter Munch    */
7878c81f8b0SPeter Munch   ~OperatorDealii() = default;
7888c81f8b0SPeter Munch 
7898c81f8b0SPeter Munch   /**
7908c81f8b0SPeter Munch    * Initialized internal data structures, particularly, MatrixFree.
7918c81f8b0SPeter Munch    */
7928c81f8b0SPeter Munch   void
7938c81f8b0SPeter Munch   reinit() override
7948c81f8b0SPeter Munch   {
7958c81f8b0SPeter Munch     // configure MatrixFree
7968c81f8b0SPeter Munch     typename MatrixFree<dim, Number>::AdditionalData additional_data;
7978c81f8b0SPeter Munch     additional_data.tasks_parallel_scheme =
7988c81f8b0SPeter Munch       MatrixFree<dim, Number>::AdditionalData::TasksParallelScheme::none;
7998c81f8b0SPeter Munch 
8008c81f8b0SPeter Munch     // create MatrixFree
8018c81f8b0SPeter Munch     matrix_free.reinit(mapping, dof_handler, constraints, quadrature, additional_data);
8028c81f8b0SPeter Munch   }
8038c81f8b0SPeter Munch 
8048c81f8b0SPeter Munch   /**
8058c81f8b0SPeter Munch    * Matrix-vector product.
8068c81f8b0SPeter Munch    */
8078c81f8b0SPeter Munch   void
8088c81f8b0SPeter Munch   vmult(VectorType &dst, const VectorType &src) const override
8098c81f8b0SPeter Munch   {
8108c81f8b0SPeter Munch     if (dof_handler.get_fe().n_components() == 1)
8118c81f8b0SPeter Munch       {
8128c81f8b0SPeter Munch         matrix_free.cell_loop(&OperatorDealii::do_cell_integral_range<1>, this, dst, src, true);
8138c81f8b0SPeter Munch       }
8148c81f8b0SPeter Munch     else
8158c81f8b0SPeter Munch       {
8168c81f8b0SPeter Munch         AssertThrow(dof_handler.get_fe().n_components() == dim, ExcInternalError());
8178c81f8b0SPeter Munch 
8188c81f8b0SPeter Munch         matrix_free.cell_loop(&OperatorDealii::do_cell_integral_range<dim>, this, dst, src, true);
8198c81f8b0SPeter Munch       }
8208c81f8b0SPeter Munch   }
8218c81f8b0SPeter Munch 
8228c81f8b0SPeter Munch   /**
8238c81f8b0SPeter Munch    * Initialize vector.
8248c81f8b0SPeter Munch    */
8258c81f8b0SPeter Munch   void
8268c81f8b0SPeter Munch   initialize_dof_vector(VectorType &vec) const override
8278c81f8b0SPeter Munch   {
8288c81f8b0SPeter Munch     matrix_free.initialize_dof_vector(vec);
8298c81f8b0SPeter Munch   }
8308c81f8b0SPeter Munch 
8318c81f8b0SPeter Munch   /**
8328c81f8b0SPeter Munch    * Compute inverse of diagonal.
8338c81f8b0SPeter Munch    */
8348c81f8b0SPeter Munch   void
8358c81f8b0SPeter Munch   compute_inverse_diagonal(VectorType &diagonal) const override
8368c81f8b0SPeter Munch   {
8378c81f8b0SPeter Munch     this->initialize_dof_vector(diagonal);
8388c81f8b0SPeter Munch 
8398c81f8b0SPeter Munch     if (dof_handler.get_fe().n_components() == 1)
8408c81f8b0SPeter Munch       {
8418c81f8b0SPeter Munch         MatrixFreeTools::compute_diagonal(matrix_free,
8428c81f8b0SPeter Munch                                           diagonal,
8438c81f8b0SPeter Munch                                           &OperatorDealii::do_cell_integral_local<1>,
8448c81f8b0SPeter Munch                                           this);
8458c81f8b0SPeter Munch       }
8468c81f8b0SPeter Munch     else
8478c81f8b0SPeter Munch       {
8488c81f8b0SPeter Munch         AssertThrow(dof_handler.get_fe().n_components() == dim, ExcInternalError());
8498c81f8b0SPeter Munch 
8508c81f8b0SPeter Munch         MatrixFreeTools::compute_diagonal(matrix_free,
8518c81f8b0SPeter Munch                                           diagonal,
8528c81f8b0SPeter Munch                                           &OperatorDealii::do_cell_integral_local<dim>,
8538c81f8b0SPeter Munch                                           this);
8548c81f8b0SPeter Munch       }
8558c81f8b0SPeter Munch 
8568c81f8b0SPeter Munch     for (auto &i : diagonal)
8578c81f8b0SPeter Munch       i = (std::abs(i) > 1.0e-10) ? (1.0 / i) : 1.0;
8588c81f8b0SPeter Munch   }
8598c81f8b0SPeter Munch 
8608c81f8b0SPeter Munch private:
8618c81f8b0SPeter Munch   /**
8628c81f8b0SPeter Munch    * Cell integral without vector access.
8638c81f8b0SPeter Munch    */
8648c81f8b0SPeter Munch   template <int n_components>
8658c81f8b0SPeter Munch   void
8668c81f8b0SPeter Munch   do_cell_integral_local(FEEvaluation<dim, -1, 0, n_components, Number> &phi) const
8678c81f8b0SPeter Munch   {
8688c81f8b0SPeter Munch     if (bp <= BPType::BP2) // mass matrix
8698c81f8b0SPeter Munch       {
8708c81f8b0SPeter Munch         phi.evaluate(EvaluationFlags::values);
8718c81f8b0SPeter Munch         for (const auto q : phi.quadrature_point_indices())
8728c81f8b0SPeter Munch           phi.submit_value(phi.get_value(q), q);
8738c81f8b0SPeter Munch         phi.integrate(EvaluationFlags::values);
8748c81f8b0SPeter Munch       }
8758c81f8b0SPeter Munch     else // Poisson operator
8768c81f8b0SPeter Munch       {
8778c81f8b0SPeter Munch         phi.evaluate(EvaluationFlags::gradients);
8788c81f8b0SPeter Munch         for (const auto q : phi.quadrature_point_indices())
8798c81f8b0SPeter Munch           phi.submit_gradient(phi.get_gradient(q), q);
8808c81f8b0SPeter Munch         phi.integrate(EvaluationFlags::gradients);
8818c81f8b0SPeter Munch       }
8828c81f8b0SPeter Munch   }
8838c81f8b0SPeter Munch 
8848c81f8b0SPeter Munch   /**
8858c81f8b0SPeter Munch    * Cell integral on a range of cells.
8868c81f8b0SPeter Munch    */
8878c81f8b0SPeter Munch   template <int n_components>
8888c81f8b0SPeter Munch   void
8898c81f8b0SPeter Munch   do_cell_integral_range(const MatrixFree<dim, Number>               &matrix_free,
8908c81f8b0SPeter Munch                          VectorType                                  &dst,
8918c81f8b0SPeter Munch                          const VectorType                            &src,
8928c81f8b0SPeter Munch                          const std::pair<unsigned int, unsigned int> &range) const
8938c81f8b0SPeter Munch   {
8948c81f8b0SPeter Munch     FEEvaluation<dim, -1, 0, n_components, Number> phi(matrix_free, range);
8958c81f8b0SPeter Munch 
8968c81f8b0SPeter Munch     for (unsigned cell = range.first; cell < range.second; ++cell)
8978c81f8b0SPeter Munch       {
8988c81f8b0SPeter Munch         phi.reinit(cell);
8998c81f8b0SPeter Munch         phi.read_dof_values(src);            // read source vector
9008c81f8b0SPeter Munch         do_cell_integral_local(phi);         // cell integral
9018c81f8b0SPeter Munch         phi.distribute_local_to_global(dst); // write to destination vector
9028c81f8b0SPeter Munch       }
9038c81f8b0SPeter Munch   }
9048c81f8b0SPeter Munch 
9058c81f8b0SPeter Munch   /**
9068c81f8b0SPeter Munch    * Mapping object passed to the constructor.
9078c81f8b0SPeter Munch    */
9088c81f8b0SPeter Munch   const Mapping<dim> &mapping;
9098c81f8b0SPeter Munch 
9108c81f8b0SPeter Munch   /**
9118c81f8b0SPeter Munch    * DoFHandler object passed to the constructor.
9128c81f8b0SPeter Munch    */
9138c81f8b0SPeter Munch   const DoFHandler<dim> &dof_handler;
9148c81f8b0SPeter Munch 
9158c81f8b0SPeter Munch   /**
9168c81f8b0SPeter Munch    * Constraints object passed to the constructor.
9178c81f8b0SPeter Munch    */
9188c81f8b0SPeter Munch   const AffineConstraints<Number> &constraints;
9198c81f8b0SPeter Munch 
9208c81f8b0SPeter Munch   /**
9218c81f8b0SPeter Munch    * Quadrature rule object passed to the constructor.
9228c81f8b0SPeter Munch    */
9238c81f8b0SPeter Munch   const Quadrature<dim> &quadrature;
9248c81f8b0SPeter Munch 
9258c81f8b0SPeter Munch   /**
9268c81f8b0SPeter Munch    * Selected BP.
9278c81f8b0SPeter Munch    */
9288c81f8b0SPeter Munch   const BPType bp;
9298c81f8b0SPeter Munch 
9308c81f8b0SPeter Munch   /**
9318c81f8b0SPeter Munch    * MatrixFree object.
9328c81f8b0SPeter Munch    */
9338c81f8b0SPeter Munch   MatrixFree<dim, Number> matrix_free;
9348c81f8b0SPeter Munch };
935