xref: /libCEED/examples/deal.II/bps.h (revision 39472f1800a7201e795b8e67cf33eaec17608dee)
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->metric_data = 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, metric_data.size(), &q_data);
285     CeedVectorSetArray(q_data, CEED_MEM_HOST, CEED_USE_POINTER, metric_data.data());
286     CeedElemRestrictionCreateStrided(ceed,
287                                      n_local_active_cells,
288                                      quadrature.size(),
289                                      n_components_metric,
290                                      metric_data.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     // pass memory buffers to libCEED
357     VectorTypeCeed x(src_ceed);
358     VectorTypeCeed y(dst_ceed);
359     x.import_array(src, CEED_MEM_HOST);
360     y.import_array(dst, CEED_MEM_HOST);
361 
362     // apply operator
363     CeedOperatorApply(op_apply, x(), y(), CEED_REQUEST_IMMEDIATE);
364 
365     // pull arrays back to deal.II
366     x.take_array();
367     y.take_array();
368 
369     // communicate: compress
370     src.zero_out_ghost_values();
371     dst.compress(VectorOperation::add);
372 
373     // apply constraints: we assume homogeneous DBC
374     constraints.set_zero(dst);
375   }
376 
377   /**
378    * Initialized vector.
379    */
380   void
381   initialize_dof_vector(VectorType &vec) const override
382   {
383     vec.reinit(partitioner);
384   }
385 
386   /**
387    * Compute inverse of diagonal.
388    */
389   void
390   compute_inverse_diagonal(VectorType &diagonal) const override
391   {
392     this->initialize_dof_vector(diagonal);
393 
394     // pass memory buffer to libCEED
395     VectorTypeCeed y(dst_ceed);
396     y.import_array(diagonal, CEED_MEM_HOST);
397 
398     CeedOperatorLinearAssembleDiagonal(op_apply, y(), CEED_REQUEST_IMMEDIATE);
399 
400     // pull array back to deal.II
401     y.take_array();
402 
403     diagonal.compress(VectorOperation::add);
404 
405     // apply constraints: we assume homogeneous DBC
406     constraints.set_zero(diagonal);
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      * Take previously set deal.II array from libCEED vector
458      */
459     void
460     take_array()
461     {
462       CeedScalar *ptr;
463       CeedVectorTakeArray(vec_ceed, mem_space, &ptr);
464     }
465 
466     /**
467      * Destructor: destroy vector view.
468      */
469     ~VectorTypeCeed()
470     {
471       bool has_array;
472       CeedVectorHasBorrowedArrayOfType(vec_ceed, mem_space, &has_array);
473       if (has_array)
474         {
475           CeedScalar *ptr;
476           CeedVectorTakeArray(vec_ceed, mem_space, &ptr);
477         }
478       CeedVectorDestroy(&vec_ceed);
479     }
480 
481   private:
482     /**
483      * libCEED vector view.
484      */
485     CeedMemType mem_space;
486     CeedVector  vec_ceed;
487   };
488 
489   /**
490    * Number of locally active DoFs.
491    */
492   unsigned int
493   extended_local_size() const
494   {
495     return partitioner->locally_owned_size() + partitioner->n_ghost_indices();
496   }
497 
498   /**
499    * Compute metric data: Jacobian, ...
500    */
501   static std::vector<double>
502   compute_metric_data(const Ceed               &ceed,
503                       const Mapping<dim>       &mapping,
504                       const Triangulation<dim> &tria,
505                       const Quadrature<dim>    &quadrature,
506                       const BPType              bp)
507   {
508     std::vector<double> metric_data;
509 
510     CeedBasis            geo_basis;
511     CeedVector           q_data;
512     CeedElemRestriction  q_data_restriction;
513     CeedVector           node_coords;
514     CeedElemRestriction  geo_restriction;
515     CeedQFunctionContext build_ctx;
516     CeedQFunction        qf_build;
517     CeedOperator         op_build;
518 
519     const unsigned int n_q_points = quadrature.get_tensor_basis()[0].size();
520 
521     const unsigned int n_components_metric = (bp <= BPType::BP2) ? 1 : (dim * (dim + 1) / 2);
522 
523     const auto mapping_q = dynamic_cast<const MappingQ<dim> *>(&mapping);
524 
525     AssertThrow(mapping_q, ExcMessage("Wrong mapping!"));
526 
527     const unsigned int fe_degree = mapping_q->get_degree();
528 
529     FE_Q<dim> geo_fe(fe_degree);
530 
531     {
532       const dealii::internal::MatrixFreeFunctions::ShapeInfo<double> shape_info(quadrature,
533                                                                                 geo_fe,
534                                                                                 0);
535       const auto             &shape_data = shape_info.get_shape_data();
536       std::vector<CeedScalar> q_ref_1d;
537       for (const auto q : shape_data.quadrature.get_points())
538         q_ref_1d.push_back(q(0));
539 
540       // transpose bases for compatibility with restriction
541       std::vector<CeedScalar> interp_1d(shape_data.shape_values.size());
542       std::vector<CeedScalar> grad_1d(shape_data.shape_gradients.size());
543       for (unsigned int i = 0; i < n_q_points; ++i)
544         for (unsigned int j = 0; j < fe_degree + 1; ++j)
545           {
546             interp_1d[j + i * (fe_degree + 1)] = shape_data.shape_values[j * n_q_points + i];
547             grad_1d[j + i * (fe_degree + 1)]   = shape_data.shape_gradients[j * n_q_points + i];
548           }
549 
550       CeedBasisCreateTensorH1(ceed,
551                               dim,
552                               dim,
553                               fe_degree + 1,
554                               n_q_points,
555                               interp_1d.data(),
556                               grad_1d.data(),
557                               q_ref_1d.data(),
558                               quadrature.get_tensor_basis()[0].get_weights().data(),
559                               &geo_basis);
560     }
561 
562     unsigned int n_local_active_cells = 0;
563 
564     for (const auto &cell : tria.active_cell_iterators())
565       if (cell->is_locally_owned())
566         n_local_active_cells++;
567 
568     std::vector<double>  geo_support_points;
569     std::vector<CeedInt> geo_indices;
570 
571     DoFHandler<dim> geo_dof_handler(tria);
572     geo_dof_handler.distribute_dofs(geo_fe);
573 
574     const auto geo_partitioner =
575       std::make_shared<Utilities::MPI::Partitioner>(geo_dof_handler.locally_owned_dofs(),
576                                                     DoFTools::extract_locally_active_dofs(
577                                                       geo_dof_handler),
578                                                     geo_dof_handler.get_communicator());
579 
580     geo_indices.reserve(n_local_active_cells * geo_fe.n_dofs_per_cell());
581 
582     const auto dof_mapping = FETools::lexicographic_to_hierarchic_numbering<dim>(fe_degree);
583 
584     FEValues<dim> fe_values(mapping,
585                             geo_fe,
586                             geo_fe.get_unit_support_points(),
587                             update_quadrature_points);
588 
589     std::vector<types::global_dof_index> local_indices(geo_fe.n_dofs_per_cell());
590 
591     const unsigned int n_points =
592       geo_partitioner->locally_owned_size() + geo_partitioner->n_ghost_indices();
593 
594     geo_support_points.resize(dim * n_points);
595 
596     for (const auto &cell : geo_dof_handler.active_cell_iterators())
597       if (cell->is_locally_owned())
598         {
599           fe_values.reinit(cell);
600           cell->get_dof_indices(local_indices);
601 
602           for (const auto i : dof_mapping)
603             {
604               const auto index = geo_partitioner->global_to_local(local_indices[i]);
605               geo_indices.emplace_back(index * dim);
606 
607               const auto point = fe_values.quadrature_point(i);
608 
609               for (unsigned int d = 0; d < dim; ++d)
610                 geo_support_points[index * dim + d] = point[d];
611             }
612         }
613 
614     metric_data.resize(n_local_active_cells * quadrature.size() * n_components_metric);
615 
616     CeedInt strides[3] = {1,
617                           static_cast<int>(quadrature.size()),
618                           static_cast<int>(quadrature.size() * n_components_metric)};
619 
620     CeedVectorCreate(ceed, metric_data.size(), &q_data);
621     CeedVectorSetArray(q_data, CEED_MEM_HOST, CEED_USE_POINTER, metric_data.data());
622     CeedElemRestrictionCreateStrided(ceed,
623                                      n_local_active_cells,
624                                      quadrature.size(),
625                                      n_components_metric,
626                                      metric_data.size(),
627                                      strides,
628                                      &q_data_restriction);
629 
630     CeedVectorCreate(ceed, geo_support_points.size(), &node_coords);
631     CeedVectorSetArray(node_coords, CEED_MEM_HOST, CEED_USE_POINTER, geo_support_points.data());
632 
633     CeedElemRestrictionCreate(ceed,
634                               n_local_active_cells,
635                               geo_fe.n_dofs_per_cell(),
636                               dim,
637                               1,
638                               geo_support_points.size(),
639                               CEED_MEM_HOST,
640                               CEED_COPY_VALUES,
641                               geo_indices.data(),
642                               &geo_restriction);
643 
644     BuildContext build_ctx_data;
645     build_ctx_data.dim       = dim;
646     build_ctx_data.space_dim = dim;
647 
648     CeedQFunctionContextCreate(ceed, &build_ctx);
649     CeedQFunctionContextSetData(
650       build_ctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(build_ctx_data), &build_ctx_data);
651 
652     // 5) create q operation
653     if (bp <= BPType::BP2)
654       CeedQFunctionCreateInterior(ceed, 1, f_build_mass, f_build_mass_loc, &qf_build);
655     else
656       CeedQFunctionCreateInterior(ceed, 1, f_build_poisson, f_build_poisson_loc, &qf_build);
657 
658     CeedQFunctionAddInput(qf_build, "geo", dim * dim, CEED_EVAL_GRAD);
659     CeedQFunctionAddInput(qf_build, "metric_data", 1, CEED_EVAL_WEIGHT);
660     CeedQFunctionAddOutput(qf_build, "qdata", n_components_metric, CEED_EVAL_NONE);
661     CeedQFunctionSetContext(qf_build, build_ctx);
662 
663     // 6) put everything together
664     CeedOperatorCreate(ceed, qf_build, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_build);
665     CeedOperatorSetField(op_build, "geo", geo_restriction, geo_basis, CEED_VECTOR_ACTIVE);
666     CeedOperatorSetField(
667       op_build, "metric_data", CEED_ELEMRESTRICTION_NONE, geo_basis, CEED_VECTOR_NONE);
668     CeedOperatorSetField(
669       op_build, "qdata", q_data_restriction, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
670 
671     CeedOperatorApply(op_build, node_coords, q_data, CEED_REQUEST_IMMEDIATE);
672 
673     CeedVectorDestroy(&node_coords);
674     CeedVectorSyncArray(q_data, CEED_MEM_HOST);
675     CeedVectorDestroy(&q_data);
676     CeedElemRestrictionDestroy(&geo_restriction);
677     CeedElemRestrictionDestroy(&q_data_restriction);
678     CeedBasisDestroy(&geo_basis);
679     CeedQFunctionContextDestroy(&build_ctx);
680     CeedQFunctionDestroy(&qf_build);
681     CeedOperatorDestroy(&op_build);
682 
683     return metric_data;
684   }
685 
686   /**
687    * Mapping object passed to the constructor.
688    */
689   const Mapping<dim> &mapping;
690 
691   /**
692    * DoFHandler object passed to the constructor.
693    */
694   const DoFHandler<dim> &dof_handler;
695 
696   /**
697    * Constraints object passed to the constructor.
698    */
699   const AffineConstraints<Number> &constraints;
700 
701   /**
702    * Quadrature rule object passed to the constructor.
703    */
704   const Quadrature<dim> &quadrature;
705 
706   /**
707    * Selected BP.
708    */
709   const BPType bp;
710 
711   /**
712    * Resource name.
713    */
714   const std::string resource;
715 
716   /**
717    * Partitioner for distributed vectors.
718    */
719   std::shared_ptr<Utilities::MPI::Partitioner> partitioner;
720 
721   /**
722    * libCEED data structures.
723    */
724   Ceed                   ceed;
725   std::vector<double>    metric_data;
726   std::array<CeedInt, 3> strides;
727   CeedVector             src_ceed;
728   CeedVector             dst_ceed;
729   CeedOperator           op_apply;
730 };
731 
732 
733 
734 template <int dim, typename Number>
735 class OperatorDealii : public OperatorBase<Number>
736 {
737 public:
738   using VectorType = typename OperatorBase<Number>::VectorType;
739 
740   /**
741    * Constructor.
742    */
743   OperatorDealii(const Mapping<dim>              &mapping,
744                  const DoFHandler<dim>           &dof_handler,
745                  const AffineConstraints<Number> &constraints,
746                  const Quadrature<dim>           &quadrature,
747                  const BPType                    &bp)
748     : mapping(mapping)
749     , dof_handler(dof_handler)
750     , constraints(constraints)
751     , quadrature(quadrature)
752     , bp(bp)
753   {
754     reinit();
755   }
756 
757   /**
758    * Destructor.
759    */
760   ~OperatorDealii() = default;
761 
762   /**
763    * Initialized internal data structures, particularly, MatrixFree.
764    */
765   void
766   reinit() override
767   {
768     // configure MatrixFree
769     typename MatrixFree<dim, Number>::AdditionalData additional_data;
770     additional_data.tasks_parallel_scheme =
771       MatrixFree<dim, Number>::AdditionalData::TasksParallelScheme::none;
772 
773     // create MatrixFree
774     matrix_free.reinit(mapping, dof_handler, constraints, quadrature, additional_data);
775   }
776 
777   /**
778    * Matrix-vector product.
779    */
780   void
781   vmult(VectorType &dst, const VectorType &src) const override
782   {
783     if (dof_handler.get_fe().n_components() == 1)
784       {
785         matrix_free.cell_loop(&OperatorDealii::do_cell_integral_range<1>, this, dst, src, true);
786       }
787     else
788       {
789         AssertThrow(dof_handler.get_fe().n_components() == dim, ExcInternalError());
790 
791         matrix_free.cell_loop(&OperatorDealii::do_cell_integral_range<dim>, this, dst, src, true);
792       }
793   }
794 
795   /**
796    * Initialize vector.
797    */
798   void
799   initialize_dof_vector(VectorType &vec) const override
800   {
801     matrix_free.initialize_dof_vector(vec);
802   }
803 
804   /**
805    * Compute inverse of diagonal.
806    */
807   void
808   compute_inverse_diagonal(VectorType &diagonal) const override
809   {
810     this->initialize_dof_vector(diagonal);
811 
812     if (dof_handler.get_fe().n_components() == 1)
813       {
814         MatrixFreeTools::compute_diagonal(matrix_free,
815                                           diagonal,
816                                           &OperatorDealii::do_cell_integral_local<1>,
817                                           this);
818       }
819     else
820       {
821         AssertThrow(dof_handler.get_fe().n_components() == dim, ExcInternalError());
822 
823         MatrixFreeTools::compute_diagonal(matrix_free,
824                                           diagonal,
825                                           &OperatorDealii::do_cell_integral_local<dim>,
826                                           this);
827       }
828 
829     for (auto &i : diagonal)
830       i = (std::abs(i) > 1.0e-10) ? (1.0 / i) : 1.0;
831   }
832 
833 private:
834   /**
835    * Cell integral without vector access.
836    */
837   template <int n_components>
838   void
839   do_cell_integral_local(FEEvaluation<dim, -1, 0, n_components, Number> &phi) const
840   {
841     if (bp <= BPType::BP2) // mass matrix
842       {
843         phi.evaluate(EvaluationFlags::values);
844         for (const auto q : phi.quadrature_point_indices())
845           phi.submit_value(phi.get_value(q), q);
846         phi.integrate(EvaluationFlags::values);
847       }
848     else // Poisson operator
849       {
850         phi.evaluate(EvaluationFlags::gradients);
851         for (const auto q : phi.quadrature_point_indices())
852           phi.submit_gradient(phi.get_gradient(q), q);
853         phi.integrate(EvaluationFlags::gradients);
854       }
855   }
856 
857   /**
858    * Cell integral on a range of cells.
859    */
860   template <int n_components>
861   void
862   do_cell_integral_range(const MatrixFree<dim, Number>               &matrix_free,
863                          VectorType                                  &dst,
864                          const VectorType                            &src,
865                          const std::pair<unsigned int, unsigned int> &range) const
866   {
867     FEEvaluation<dim, -1, 0, n_components, Number> phi(matrix_free, range);
868 
869     for (unsigned cell = range.first; cell < range.second; ++cell)
870       {
871         phi.reinit(cell);
872         phi.read_dof_values(src);            // read source vector
873         do_cell_integral_local(phi);         // cell integral
874         phi.distribute_local_to_global(dst); // write to destination vector
875       }
876   }
877 
878   /**
879    * Mapping object passed to the constructor.
880    */
881   const Mapping<dim> &mapping;
882 
883   /**
884    * DoFHandler object passed to the constructor.
885    */
886   const DoFHandler<dim> &dof_handler;
887 
888   /**
889    * Constraints object passed to the constructor.
890    */
891   const AffineConstraints<Number> &constraints;
892 
893   /**
894    * Quadrature rule object passed to the constructor.
895    */
896   const Quadrature<dim> &quadrature;
897 
898   /**
899    * Selected BP.
900    */
901   const BPType bp;
902 
903   /**
904    * MatrixFree object.
905    */
906   MatrixFree<dim, Number> matrix_free;
907 };
908