xref: /libCEED/examples/deal.II/bps.h (revision 704703efc9822a8339f73fbc4da0d2a9396e454a)
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2023 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 //  Authors: Peter Munch, Martin Kronbichler
15 //
16 // ---------------------------------------------------------------------
17 
18 // deal.II includes
19 #include <deal.II/dofs/dof_tools.h>
20 
21 #include <deal.II/fe/mapping.h>
22 
23 #include <deal.II/lac/la_parallel_vector.h>
24 
25 #include <deal.II/matrix_free/fe_evaluation.h>
26 #include <deal.II/matrix_free/matrix_free.h>
27 #include <deal.II/matrix_free/shape_info.h>
28 #include <deal.II/matrix_free/tools.h>
29 
30 // libCEED includes
31 #include <ceed.h>
32 #include <ceed/backend.h>
33 
34 // QFunction source
35 #include "bps-qfunctions.h"
36 
37 using namespace dealii;
38 
39 /**
40  * BP types. For more details, see https://ceed.exascaleproject.org/bps/.
41  */
42 enum class BPType : unsigned int
43 {
44   BP1,
45   BP2,
46   BP3,
47   BP4,
48   BP5,
49   BP6
50 };
51 
52 
53 
54 /**
55  * Struct storing relevant information regarding each BP.
56  */
57 struct BPInfo
58 {
59   BPInfo(const BPType type, const int dim, const int fe_degree)
60     : type(type)
61     , dim(dim)
62     , fe_degree(fe_degree)
63   {
64     if (type == BPType::BP1)
65       type_string = "BP1";
66     else if (type == BPType::BP2)
67       type_string = "BP2";
68     else if (type == BPType::BP3)
69       type_string = "BP3";
70     else if (type == BPType::BP4)
71       type_string = "BP4";
72     else if (type == BPType::BP5)
73       type_string = "BP5";
74     else if (type == BPType::BP6)
75       type_string = "BP6";
76 
77     this->n_q_points_1d = (type <= BPType::BP4) ? (fe_degree + 2) : (fe_degree + 1);
78 
79     this->n_components =
80       (type == BPType::BP1 || type == BPType::BP3 || type == BPType::BP5) ? 1 : dim;
81   }
82 
83 
84   BPType       type;
85   std::string  type_string;
86   unsigned int dim;
87   unsigned int fe_degree;
88   unsigned int n_q_points_1d;
89   unsigned int n_components;
90 };
91 
92 
93 
94 /**
95  * Base class of operators.
96  */
97 template <typename Number>
98 class OperatorBase
99 {
100 public:
101   /**
102    * deal.II vector type
103    */
104   using VectorType = LinearAlgebra::distributed::Vector<Number>;
105 
106   /**
107    * Initialize vector.
108    */
109   virtual void
110   reinit() = 0;
111 
112   /**
113    * Perform matrix-vector product
114    */
115   virtual void
116   vmult(VectorType &dst, const VectorType &src) const = 0;
117 
118   /**
119    * Initialize vector.
120    */
121   virtual void
122   initialize_dof_vector(VectorType &vec) const = 0;
123 
124   /**
125    * Compute inverse of diagonal.
126    */
127   virtual void
128   compute_inverse_diagonal(VectorType &diagonal) const = 0;
129 };
130 
131 
132 /**
133  * Operator implementation using libCEED.
134  */
135 template <int dim, typename Number>
136 class OperatorCeed : public OperatorBase<Number>
137 {
138 public:
139   using VectorType = typename OperatorBase<Number>::VectorType;
140 
141   /**
142    * Constructor.
143    */
144   OperatorCeed(const Mapping<dim>              &mapping,
145                const DoFHandler<dim>           &dof_handler,
146                const AffineConstraints<Number> &constraints,
147                const Quadrature<dim>           &quadrature,
148                const BPType                    &bp,
149                const std::string               &resource)
150     : mapping(mapping)
151     , dof_handler(dof_handler)
152     , constraints(constraints)
153     , quadrature(quadrature)
154     , bp(bp)
155     , resource(resource)
156   {
157     reinit();
158   }
159 
160   /**
161    * Destructor.
162    */
163   ~OperatorCeed()
164   {
165     CeedVectorDestroy(&src_ceed);
166     CeedVectorDestroy(&dst_ceed);
167     CeedOperatorDestroy(&op_apply);
168     CeedDestroy(&ceed);
169   }
170 
171   /**
172    * Initialized internal data structures, particularly, libCEED.
173    */
174   void
175   reinit() override
176   {
177     CeedVector           q_data;
178     CeedBasis            sol_basis;
179     CeedElemRestriction  sol_restriction;
180     CeedElemRestriction  q_data_restriction;
181     BuildContext         build_ctx_data;
182     CeedQFunctionContext build_ctx;
183     CeedQFunction        qf_apply;
184 
185     const auto &tria = dof_handler.get_triangulation();
186     const auto &fe   = dof_handler.get_fe();
187 
188     const auto n_components = fe.n_components();
189 
190     if (bp == BPType::BP1 || bp == BPType::BP3 || bp == BPType::BP5)
191       {
192         AssertThrow(n_components == 1, ExcInternalError());
193       }
194     else
195       {
196         AssertThrow(n_components == dim, ExcInternalError());
197       }
198 
199     // 1) create CEED instance -> "MatrixFree"
200     const char *ceed_spec = resource.c_str();
201     CeedInit(ceed_spec, &ceed);
202 
203     // 2) create shape functions -> "ShapeInfo"
204     const unsigned int fe_degree  = fe.tensor_degree();
205     const unsigned int n_q_points = quadrature.get_tensor_basis()[0].size();
206     {
207       const dealii::internal::MatrixFreeFunctions::ShapeInfo<double> shape_info(quadrature, fe, 0);
208       const auto             &shape_data = shape_info.get_shape_data();
209       std::vector<CeedScalar> q_ref_1d;
210       for (const auto q : shape_data.quadrature.get_points())
211         q_ref_1d.push_back(q(0));
212 
213       // transpose bases for compatibility with restriction
214       std::vector<CeedScalar> interp_1d(shape_data.shape_values.size());
215       std::vector<CeedScalar> grad_1d(shape_data.shape_gradients.size());
216       for (unsigned int i = 0; i < n_q_points; ++i)
217         for (unsigned int j = 0; j < fe_degree + 1; ++j)
218           {
219             interp_1d[j + i * (fe_degree + 1)] = shape_data.shape_values[j * n_q_points + i];
220             grad_1d[j + i * (fe_degree + 1)]   = shape_data.shape_gradients[j * n_q_points + i];
221           }
222 
223       CeedBasisCreateTensorH1(ceed,
224                               dim,
225                               n_components,
226                               fe_degree + 1,
227                               n_q_points,
228                               interp_1d.data(),
229                               grad_1d.data(),
230                               q_ref_1d.data(),
231                               quadrature.get_tensor_basis()[0].get_weights().data(),
232                               &sol_basis);
233     }
234 
235     // 3) create restriction matrix -> DoFInfo
236     unsigned int n_local_active_cells = 0;
237 
238     for (const auto &cell : dof_handler.active_cell_iterators())
239       if (cell->is_locally_owned())
240         n_local_active_cells++;
241 
242     partitioner =
243       std::make_shared<Utilities::MPI::Partitioner>(dof_handler.locally_owned_dofs(),
244                                                     DoFTools::extract_locally_active_dofs(
245                                                       dof_handler),
246                                                     dof_handler.get_communicator());
247 
248     std::vector<CeedInt> indices;
249     indices.reserve(n_local_active_cells * fe.n_dofs_per_cell() / n_components);
250 
251     const auto dof_mapping = FETools::lexicographic_to_hierarchic_numbering<dim>(fe_degree);
252 
253     std::vector<types::global_dof_index> local_indices(fe.n_dofs_per_cell());
254 
255     for (const auto &cell : dof_handler.active_cell_iterators())
256       if (cell->is_locally_owned())
257         {
258           cell->get_dof_indices(local_indices);
259 
260           for (const auto i : dof_mapping)
261             indices.emplace_back(
262               partitioner->global_to_local(local_indices[fe.component_to_system_index(0, i)]));
263         }
264 
265     CeedElemRestrictionCreate(ceed,
266                               n_local_active_cells,
267                               fe.n_dofs_per_cell() / n_components,
268                               n_components,
269                               1,
270                               this->extended_local_size(),
271                               CEED_MEM_HOST,
272                               CEED_COPY_VALUES,
273                               indices.data(),
274                               &sol_restriction);
275 
276     // 4) create mapping -> MappingInfo
277     const unsigned int n_components_metric = (bp <= BPType::BP2) ? 1 : (dim * (dim + 1) / 2);
278 
279     this->weights = compute_metric_data(ceed, mapping, tria, quadrature, bp);
280 
281     strides = {{1,
282                 static_cast<int>(quadrature.size()),
283                 static_cast<int>(quadrature.size() * n_components_metric)}};
284     CeedVectorCreate(ceed, weights.size(), &q_data);
285     CeedVectorSetArray(q_data, CEED_MEM_HOST, CEED_USE_POINTER, weights.data());
286     CeedElemRestrictionCreateStrided(ceed,
287                                      n_local_active_cells,
288                                      quadrature.size(),
289                                      n_components_metric,
290                                      weights.size(),
291                                      strides.data(),
292                                      &q_data_restriction);
293 
294     build_ctx_data.dim       = dim;
295     build_ctx_data.space_dim = dim;
296 
297     CeedQFunctionContextCreate(ceed, &build_ctx);
298     CeedQFunctionContextSetData(
299       build_ctx, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(build_ctx_data), &build_ctx_data);
300 
301     // 5) create q operation
302     if (bp == BPType::BP1)
303       CeedQFunctionCreateInterior(ceed, 1, f_apply_mass, f_apply_mass_loc, &qf_apply);
304     else if (bp == BPType::BP2)
305       CeedQFunctionCreateInterior(ceed, 1, f_apply_mass_vec, f_apply_mass_vec_loc, &qf_apply);
306     else if (bp == BPType::BP3 || bp == BPType::BP5)
307       CeedQFunctionCreateInterior(ceed, 1, f_apply_poisson, f_apply_poisson_loc, &qf_apply);
308     else if (bp == BPType::BP4 || bp == BPType::BP6)
309       CeedQFunctionCreateInterior(ceed, 1, f_apply_poisson_vec, f_apply_poisson_vec_loc, &qf_apply);
310     else
311       AssertThrow(false, ExcInternalError());
312 
313     if (bp <= BPType::BP2)
314       CeedQFunctionAddInput(qf_apply, "u", n_components, CEED_EVAL_INTERP);
315     else
316       CeedQFunctionAddInput(qf_apply, "u", dim * n_components, CEED_EVAL_GRAD);
317 
318     CeedQFunctionAddInput(qf_apply, "qdata", n_components_metric, CEED_EVAL_NONE);
319 
320     if (bp <= BPType::BP2)
321       CeedQFunctionAddOutput(qf_apply, "v", n_components, CEED_EVAL_INTERP);
322     else
323       CeedQFunctionAddOutput(qf_apply, "v", dim * n_components, CEED_EVAL_GRAD);
324 
325     CeedQFunctionSetContext(qf_apply, build_ctx);
326 
327     // 6) put everything together
328     CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_apply);
329 
330     CeedOperatorSetField(op_apply, "u", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
331     CeedOperatorSetField(op_apply, "qdata", q_data_restriction, CEED_BASIS_NONE, q_data);
332     CeedOperatorSetField(op_apply, "v", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE);
333 
334     // 7) libCEED vectors
335     CeedElemRestrictionCreateVector(sol_restriction, &src_ceed, NULL);
336     CeedElemRestrictionCreateVector(sol_restriction, &dst_ceed, NULL);
337 
338     // 8) cleanup
339     CeedVectorDestroy(&q_data);
340     CeedElemRestrictionDestroy(&q_data_restriction);
341     CeedElemRestrictionDestroy(&sol_restriction);
342     CeedBasisDestroy(&sol_basis);
343     CeedQFunctionContextDestroy(&build_ctx);
344     CeedQFunctionDestroy(&qf_apply);
345   }
346 
347   /**
348    * Perform matrix-vector product.
349    */
350   void
351   vmult(VectorType &dst, const VectorType &src) const override
352   {
353     // communicate: update ghost values
354     src.update_ghost_values();
355 
356     {
357       // pass memory buffers to libCEED
358       VectorTypeCeed x(src_ceed);
359       VectorTypeCeed y(dst_ceed);
360       x.import_array(src, CEED_MEM_HOST);
361       y.import_array(dst, CEED_MEM_HOST);
362 
363       // apply operator
364       CeedOperatorApply(op_apply, x(), y(), CEED_REQUEST_IMMEDIATE);
365 
366       // pull arrays back to deal.II
367       x.sync_array();
368       y.sync_array();
369     }
370     // communicate: compress
371     src.zero_out_ghost_values();
372     dst.compress(VectorOperation::add);
373 
374     // apply constraints: we assume homogeneous DBC
375     constraints.set_zero(dst);
376   }
377 
378   /**
379    * Initialized vector.
380    */
381   void
382   initialize_dof_vector(VectorType &vec) const override
383   {
384     vec.reinit(partitioner);
385   }
386 
387   /**
388    * Compute inverse of diagonal.
389    */
390   void
391   compute_inverse_diagonal(VectorType &diagonal) const override
392   {
393     this->initialize_dof_vector(diagonal);
394 
395     {
396       // pass memory buffer to libCEED
397       VectorTypeCeed y(dst_ceed);
398       y.import_array(diagonal, CEED_MEM_HOST);
399 
400       CeedOperatorLinearAssembleDiagonal(op_apply, y(), CEED_REQUEST_IMMEDIATE);
401 
402       // pull array back to deal.II
403       y.sync_array();
404     }
405 
406     diagonal.compress(VectorOperation::add);
407 
408     for (auto &i : diagonal)
409       i = (std::abs(i) > 1.0e-10) ? (1.0 / i) : 1.0;
410   }
411 
412 private:
413   /**
414    * Wrapper around a deal.II vector to create a libCEED vector view.
415    */
416   class VectorTypeCeed
417   {
418   public:
419     /**
420      * Constructor.
421      */
422     VectorTypeCeed(const CeedVector &vec_orig)
423     {
424       vec_ceed = NULL;
425       CeedVectorReferenceCopy(vec_orig, &vec_ceed);
426     }
427 
428     /**
429      * Return libCEED vector view.
430      */
431     CeedVector &
432     operator()()
433     {
434       return vec_ceed;
435     }
436 
437     /**
438      * Set deal.II memory in libCEED vector.
439      */
440     void
441     import_array(const VectorType &vec, const CeedMemType space)
442     {
443       mem_space = space;
444       CeedVectorSetArray(vec_ceed, mem_space, CEED_USE_POINTER, vec.get_values());
445     }
446 
447     /**
448      * Sync memory from device to host.
449      */
450     void
451     sync_array()
452     {
453       CeedVectorSyncArray(vec_ceed, mem_space);
454     }
455 
456     /**
457      * Destructor: destroy vector view.
458      */
459     ~VectorTypeCeed()
460     {
461       bool has_array;
462       CeedVectorHasBorrowedArrayOfType(vec_ceed, mem_space, &has_array);
463       if (has_array)
464         {
465           CeedScalar *ptr;
466           CeedVectorTakeArray(vec_ceed, mem_space, &ptr);
467         }
468       CeedVectorDestroy(&vec_ceed);
469     }
470 
471   private:
472     /**
473      * libCEED vector view.
474      */
475     CeedMemType mem_space;
476     CeedVector  vec_ceed;
477   };
478 
479   /**
480    * Number of locally active DoFs.
481    */
482   unsigned int
483   extended_local_size() const
484   {
485     return partitioner->locally_owned_size() + partitioner->n_ghost_indices();
486   }
487 
488   /**
489    * Compute metric data: Jacobian, ...
490    */
491   static std::vector<double>
492   compute_metric_data(const Ceed               &ceed,
493                       const Mapping<dim>       &mapping,
494                       const Triangulation<dim> &tria,
495                       const Quadrature<dim>    &quadrature,
496                       const BPType              bp)
497   {
498     std::vector<double> weights;
499 
500     if (false)
501       {
502         FE_Nothing<dim> dummy_fe;
503         FEValues<dim>   fe_values(mapping, dummy_fe, quadrature, update_JxW_values);
504 
505         for (const auto &cell : tria.active_cell_iterators())
506           if (cell->is_locally_owned())
507             {
508               fe_values.reinit(cell);
509 
510               for (const auto q : fe_values.quadrature_point_indices())
511                 weights.emplace_back(fe_values.JxW(q));
512             }
513 
514         return weights;
515       }
516 
517     CeedBasis            geo_basis;
518     CeedVector           q_data;
519     CeedElemRestriction  q_data_restriction;
520     CeedVector           node_coords;
521     CeedElemRestriction  geo_restriction;
522     CeedQFunctionContext build_ctx;
523     CeedQFunction        qf_build;
524     CeedOperator         op_build;
525 
526     const unsigned int n_q_points = quadrature.get_tensor_basis()[0].size();
527 
528     const unsigned int n_components = (bp <= BPType::BP2) ? 1 : (dim * (dim + 1) / 2);
529 
530     const auto mapping_q = dynamic_cast<const MappingQ<dim> *>(&mapping);
531 
532     AssertThrow(mapping_q, ExcMessage("Wrong mapping!"));
533 
534     const unsigned int fe_degree = mapping_q->get_degree();
535 
536     FE_Q<dim> geo_fe(fe_degree);
537 
538     {
539       const dealii::internal::MatrixFreeFunctions::ShapeInfo<double> shape_info(quadrature,
540                                                                                 geo_fe,
541                                                                                 0);
542       const auto             &shape_data = shape_info.get_shape_data();
543       std::vector<CeedScalar> q_ref_1d;
544       for (const auto q : shape_data.quadrature.get_points())
545         q_ref_1d.push_back(q(0));
546 
547       // transpose bases for compatibility with restriction
548       std::vector<CeedScalar> interp_1d(shape_data.shape_values.size());
549       std::vector<CeedScalar> grad_1d(shape_data.shape_gradients.size());
550       for (unsigned int i = 0; i < n_q_points; ++i)
551         for (unsigned int j = 0; j < fe_degree + 1; ++j)
552           {
553             interp_1d[j + i * (fe_degree + 1)] = shape_data.shape_values[j * n_q_points + i];
554             grad_1d[j + i * (fe_degree + 1)]   = shape_data.shape_gradients[j * n_q_points + i];
555           }
556 
557       CeedBasisCreateTensorH1(ceed,
558                               dim,
559                               dim,
560                               fe_degree + 1,
561                               n_q_points,
562                               interp_1d.data(),
563                               grad_1d.data(),
564                               q_ref_1d.data(),
565                               quadrature.get_tensor_basis()[0].get_weights().data(),
566                               &geo_basis);
567     }
568 
569     unsigned int n_local_active_cells = 0;
570 
571     for (const auto &cell : tria.active_cell_iterators())
572       if (cell->is_locally_owned())
573         n_local_active_cells++;
574 
575     std::vector<double>  geo_support_points;
576     std::vector<CeedInt> geo_indices;
577 
578     DoFHandler<dim> geo_dof_handler(tria);
579     geo_dof_handler.distribute_dofs(geo_fe);
580 
581     const auto geo_partitioner =
582       std::make_shared<Utilities::MPI::Partitioner>(geo_dof_handler.locally_owned_dofs(),
583                                                     DoFTools::extract_locally_active_dofs(
584                                                       geo_dof_handler),
585                                                     geo_dof_handler.get_communicator());
586 
587     geo_indices.reserve(n_local_active_cells * geo_fe.n_dofs_per_cell());
588 
589     const auto dof_mapping = FETools::lexicographic_to_hierarchic_numbering<dim>(fe_degree);
590 
591     FEValues<dim> fe_values(mapping,
592                             geo_fe,
593                             geo_fe.get_unit_support_points(),
594                             update_quadrature_points);
595 
596     std::vector<types::global_dof_index> local_indices(geo_fe.n_dofs_per_cell());
597 
598     const unsigned int n_points =
599       geo_partitioner->locally_owned_size() + geo_partitioner->n_ghost_indices();
600 
601     geo_support_points.resize(dim * n_points);
602 
603     for (const auto &cell : geo_dof_handler.active_cell_iterators())
604       if (cell->is_locally_owned())
605         {
606           fe_values.reinit(cell);
607           cell->get_dof_indices(local_indices);
608 
609           for (const auto i : dof_mapping)
610             {
611               const auto index = geo_partitioner->global_to_local(local_indices[i]);
612               geo_indices.emplace_back(index * dim);
613 
614               const auto point = fe_values.quadrature_point(i);
615 
616               for (unsigned int d = 0; d < dim; ++d)
617                 geo_support_points[index * dim + d] = point[d];
618             }
619         }
620 
621     weights.resize(n_local_active_cells * quadrature.size() * n_components);
622 
623     CeedInt strides[3] = {1,
624                           static_cast<int>(quadrature.size()),
625                           static_cast<int>(quadrature.size() * n_components)};
626 
627     CeedVectorCreate(ceed, weights.size(), &q_data);
628     CeedVectorSetArray(q_data, CEED_MEM_HOST, CEED_USE_POINTER, weights.data());
629     CeedElemRestrictionCreateStrided(ceed,
630                                      n_local_active_cells,
631                                      quadrature.size(),
632                                      n_components,
633                                      weights.size(),
634                                      strides,
635                                      &q_data_restriction);
636 
637     CeedVectorCreate(ceed, geo_support_points.size(), &node_coords);
638     CeedVectorSetArray(node_coords, CEED_MEM_HOST, CEED_USE_POINTER, geo_support_points.data());
639 
640     CeedElemRestrictionCreate(ceed,
641                               n_local_active_cells,
642                               geo_fe.n_dofs_per_cell(),
643                               dim,
644                               1,
645                               geo_support_points.size(),
646                               CEED_MEM_HOST,
647                               CEED_COPY_VALUES,
648                               geo_indices.data(),
649                               &geo_restriction);
650 
651     BuildContext build_ctx_data;
652     build_ctx_data.dim       = dim;
653     build_ctx_data.space_dim = dim;
654 
655     CeedQFunctionContextCreate(ceed, &build_ctx);
656     CeedQFunctionContextSetData(
657       build_ctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(build_ctx_data), &build_ctx_data);
658 
659     // 5) create q operation
660     if (bp <= BPType::BP2)
661       CeedQFunctionCreateInterior(ceed, 1, f_build_mass, f_build_mass_loc, &qf_build);
662     else
663       CeedQFunctionCreateInterior(ceed, 1, f_build_poisson, f_build_poisson_loc, &qf_build);
664 
665     CeedQFunctionAddInput(qf_build, "geo", dim * dim, CEED_EVAL_GRAD);
666     CeedQFunctionAddInput(qf_build, "weights", 1, CEED_EVAL_WEIGHT);
667     CeedQFunctionAddOutput(qf_build, "qdata", n_components, CEED_EVAL_NONE);
668     CeedQFunctionSetContext(qf_build, build_ctx);
669 
670     // 6) put everything together
671     CeedOperatorCreate(ceed, qf_build, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_build);
672     CeedOperatorSetField(op_build, "geo", geo_restriction, geo_basis, CEED_VECTOR_ACTIVE);
673     CeedOperatorSetField(
674       op_build, "weights", CEED_ELEMRESTRICTION_NONE, geo_basis, CEED_VECTOR_NONE);
675     CeedOperatorSetField(
676       op_build, "qdata", q_data_restriction, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
677 
678     CeedOperatorApply(op_build, node_coords, q_data, CEED_REQUEST_IMMEDIATE);
679 
680     CeedVectorDestroy(&node_coords);
681     CeedVectorSyncArray(q_data, CEED_MEM_HOST);
682     CeedVectorDestroy(&q_data);
683     CeedElemRestrictionDestroy(&geo_restriction);
684     CeedElemRestrictionDestroy(&q_data_restriction);
685     CeedBasisDestroy(&geo_basis);
686     CeedQFunctionContextDestroy(&build_ctx);
687     CeedQFunctionDestroy(&qf_build);
688     CeedOperatorDestroy(&op_build);
689 
690     return weights;
691   }
692 
693   /**
694    * Mapping object passed to the constructor.
695    */
696   const Mapping<dim> &mapping;
697 
698   /**
699    * DoFHandler object passed to the constructor.
700    */
701   const DoFHandler<dim> &dof_handler;
702 
703   /**
704    * Constraints object passed to the constructor.
705    */
706   const AffineConstraints<Number> &constraints;
707 
708   /**
709    * Quadrature rule object passed to the constructor.
710    */
711   const Quadrature<dim> &quadrature;
712 
713   /**
714    * Selected BP.
715    */
716   const BPType bp;
717 
718   /**
719    * Resource name.
720    */
721   const std::string resource;
722 
723   /**
724    * Partitioner for distributed vectors.
725    */
726   std::shared_ptr<Utilities::MPI::Partitioner> partitioner;
727 
728   /**
729    * libCEED data structures.
730    */
731   Ceed                   ceed;
732   std::vector<double>    weights;
733   std::array<CeedInt, 3> strides;
734   CeedVector             src_ceed;
735   CeedVector             dst_ceed;
736   CeedOperator           op_apply;
737 
738   /**
739    * Temporal (tempral) vectors.
740    *
741    * @note Only needed for multiple components.
742    */
743   mutable VectorType src_tmp;
744   mutable VectorType dst_tmp;
745 };
746 
747 
748 
749 template <int dim, typename Number>
750 class OperatorDealii : public OperatorBase<Number>
751 {
752 public:
753   using VectorType = typename OperatorBase<Number>::VectorType;
754 
755   /**
756    * Constructor.
757    */
758   OperatorDealii(const Mapping<dim>              &mapping,
759                  const DoFHandler<dim>           &dof_handler,
760                  const AffineConstraints<Number> &constraints,
761                  const Quadrature<dim>           &quadrature,
762                  const BPType                    &bp)
763     : mapping(mapping)
764     , dof_handler(dof_handler)
765     , constraints(constraints)
766     , quadrature(quadrature)
767     , bp(bp)
768   {
769     reinit();
770   }
771 
772   /**
773    * Destructor.
774    */
775   ~OperatorDealii() = default;
776 
777   /**
778    * Initialized internal data structures, particularly, MatrixFree.
779    */
780   void
781   reinit() override
782   {
783     // configure MatrixFree
784     typename MatrixFree<dim, Number>::AdditionalData additional_data;
785     additional_data.tasks_parallel_scheme =
786       MatrixFree<dim, Number>::AdditionalData::TasksParallelScheme::none;
787 
788     // create MatrixFree
789     matrix_free.reinit(mapping, dof_handler, constraints, quadrature, additional_data);
790   }
791 
792   /**
793    * Matrix-vector product.
794    */
795   void
796   vmult(VectorType &dst, const VectorType &src) const override
797   {
798     if (dof_handler.get_fe().n_components() == 1)
799       {
800         matrix_free.cell_loop(&OperatorDealii::do_cell_integral_range<1>, this, dst, src, true);
801       }
802     else
803       {
804         AssertThrow(dof_handler.get_fe().n_components() == dim, ExcInternalError());
805 
806         matrix_free.cell_loop(&OperatorDealii::do_cell_integral_range<dim>, this, dst, src, true);
807       }
808   }
809 
810   /**
811    * Initialize vector.
812    */
813   void
814   initialize_dof_vector(VectorType &vec) const override
815   {
816     matrix_free.initialize_dof_vector(vec);
817   }
818 
819   /**
820    * Compute inverse of diagonal.
821    */
822   void
823   compute_inverse_diagonal(VectorType &diagonal) const override
824   {
825     this->initialize_dof_vector(diagonal);
826 
827     if (dof_handler.get_fe().n_components() == 1)
828       {
829         MatrixFreeTools::compute_diagonal(matrix_free,
830                                           diagonal,
831                                           &OperatorDealii::do_cell_integral_local<1>,
832                                           this);
833       }
834     else
835       {
836         AssertThrow(dof_handler.get_fe().n_components() == dim, ExcInternalError());
837 
838         MatrixFreeTools::compute_diagonal(matrix_free,
839                                           diagonal,
840                                           &OperatorDealii::do_cell_integral_local<dim>,
841                                           this);
842       }
843 
844     for (auto &i : diagonal)
845       i = (std::abs(i) > 1.0e-10) ? (1.0 / i) : 1.0;
846   }
847 
848 private:
849   /**
850    * Cell integral without vector access.
851    */
852   template <int n_components>
853   void
854   do_cell_integral_local(FEEvaluation<dim, -1, 0, n_components, Number> &phi) const
855   {
856     if (bp <= BPType::BP2) // mass matrix
857       {
858         phi.evaluate(EvaluationFlags::values);
859         for (const auto q : phi.quadrature_point_indices())
860           phi.submit_value(phi.get_value(q), q);
861         phi.integrate(EvaluationFlags::values);
862       }
863     else // Poisson operator
864       {
865         phi.evaluate(EvaluationFlags::gradients);
866         for (const auto q : phi.quadrature_point_indices())
867           phi.submit_gradient(phi.get_gradient(q), q);
868         phi.integrate(EvaluationFlags::gradients);
869       }
870   }
871 
872   /**
873    * Cell integral on a range of cells.
874    */
875   template <int n_components>
876   void
877   do_cell_integral_range(const MatrixFree<dim, Number>               &matrix_free,
878                          VectorType                                  &dst,
879                          const VectorType                            &src,
880                          const std::pair<unsigned int, unsigned int> &range) const
881   {
882     FEEvaluation<dim, -1, 0, n_components, Number> phi(matrix_free, range);
883 
884     for (unsigned cell = range.first; cell < range.second; ++cell)
885       {
886         phi.reinit(cell);
887         phi.read_dof_values(src);            // read source vector
888         do_cell_integral_local(phi);         // cell integral
889         phi.distribute_local_to_global(dst); // write to destination vector
890       }
891   }
892 
893   /**
894    * Mapping object passed to the constructor.
895    */
896   const Mapping<dim> &mapping;
897 
898   /**
899    * DoFHandler object passed to the constructor.
900    */
901   const DoFHandler<dim> &dof_handler;
902 
903   /**
904    * Constraints object passed to the constructor.
905    */
906   const AffineConstraints<Number> &constraints;
907 
908   /**
909    * Quadrature rule object passed to the constructor.
910    */
911   const Quadrature<dim> &quadrature;
912 
913   /**
914    * Selected BP.
915    */
916   const BPType bp;
917 
918   /**
919    * MatrixFree object.
920    */
921   MatrixFree<dim, Number> matrix_free;
922 };
923