xref: /libCEED/examples/deal.II/bps-kokkos.h (revision a44bca27e5d331f553a82497b0c34f31883da506)
1*d416dc2bSJeremy L Thompson // ---------------------------------------------------------------------
2*d416dc2bSJeremy L Thompson //
3*d416dc2bSJeremy L Thompson // Copyright (C) 2023 by the deal.II authors
4*d416dc2bSJeremy L Thompson //
5*d416dc2bSJeremy L Thompson // This file is part of the deal.II library.
6*d416dc2bSJeremy L Thompson //
7*d416dc2bSJeremy L Thompson // The deal.II library is free software; you can use it, redistribute
8*d416dc2bSJeremy L Thompson // it, and/or modify it under the terms of the GNU Lesser General
9*d416dc2bSJeremy L Thompson // Public License as published by the Free Software Foundation; either
10*d416dc2bSJeremy L Thompson // version 2.1 of the License, or (at your option) any later version.
11*d416dc2bSJeremy L Thompson // The full text of the license can be found in the file LICENSE.md at
12*d416dc2bSJeremy L Thompson // the top level directory of deal.II.
13*d416dc2bSJeremy L Thompson //
14*d416dc2bSJeremy L Thompson //  Authors: Peter Munch, Martin Kronbichler
15*d416dc2bSJeremy L Thompson //
16*d416dc2bSJeremy L Thompson // ---------------------------------------------------------------------
17*d416dc2bSJeremy L Thompson 
18*d416dc2bSJeremy L Thompson #pragma once
19*d416dc2bSJeremy L Thompson #ifndef bps_kokkos_h
20*d416dc2bSJeremy L Thompson #  define bps_kokkos_h
21*d416dc2bSJeremy L Thompson 
22*d416dc2bSJeremy L Thompson // deal.II includes
23*d416dc2bSJeremy L Thompson #  include <deal.II/dofs/dof_tools.h>
24*d416dc2bSJeremy L Thompson 
25*d416dc2bSJeremy L Thompson #  include <deal.II/fe/mapping.h>
26*d416dc2bSJeremy L Thompson 
27*d416dc2bSJeremy L Thompson #  include <deal.II/lac/la_parallel_vector.h>
28*d416dc2bSJeremy L Thompson 
29*d416dc2bSJeremy L Thompson #  include <deal.II/matrix_free/fe_evaluation.h>
30*d416dc2bSJeremy L Thompson #  include <deal.II/matrix_free/matrix_free.h>
31*d416dc2bSJeremy L Thompson #  include <deal.II/matrix_free/shape_info.h>
32*d416dc2bSJeremy L Thompson #  include <deal.II/matrix_free/tools.h>
33*d416dc2bSJeremy L Thompson 
34*d416dc2bSJeremy L Thompson // local includes
35*d416dc2bSJeremy L Thompson #  include "bps.h"
36*d416dc2bSJeremy L Thompson 
37*d416dc2bSJeremy L Thompson using namespace dealii;
38*d416dc2bSJeremy L Thompson 
39*d416dc2bSJeremy L Thompson 
40*d416dc2bSJeremy L Thompson 
41*d416dc2bSJeremy L Thompson template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename Number>
42*d416dc2bSJeremy L Thompson class OperatorDealiiMassQuad
43*d416dc2bSJeremy L Thompson {
44*d416dc2bSJeremy L Thompson public:
45*d416dc2bSJeremy L Thompson   DEAL_II_HOST_DEVICE void
operator()46*d416dc2bSJeremy L Thompson   operator()(Portable::FEEvaluation<dim, fe_degree, n_q_points_1d, n_components, Number> *fe_eval,
47*d416dc2bSJeremy L Thompson              const int q_point) const
48*d416dc2bSJeremy L Thompson   {
49*d416dc2bSJeremy L Thompson     fe_eval->submit_value(fe_eval->get_value(q_point), q_point);
50*d416dc2bSJeremy L Thompson   }
51*d416dc2bSJeremy L Thompson };
52*d416dc2bSJeremy L Thompson 
53*d416dc2bSJeremy L Thompson 
54*d416dc2bSJeremy L Thompson 
55*d416dc2bSJeremy L Thompson template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename Number>
56*d416dc2bSJeremy L Thompson class OperatorDealiiLaplaceQuad
57*d416dc2bSJeremy L Thompson {
58*d416dc2bSJeremy L Thompson public:
59*d416dc2bSJeremy L Thompson   DEAL_II_HOST_DEVICE void
operator()60*d416dc2bSJeremy L Thompson   operator()(Portable::FEEvaluation<dim, fe_degree, n_q_points_1d, n_components, Number> *fe_eval,
61*d416dc2bSJeremy L Thompson              const int q_point) const
62*d416dc2bSJeremy L Thompson   {
63*d416dc2bSJeremy L Thompson     fe_eval->submit_gradient(fe_eval->get_gradient(q_point), q_point);
64*d416dc2bSJeremy L Thompson   }
65*d416dc2bSJeremy L Thompson };
66*d416dc2bSJeremy L Thompson 
67*d416dc2bSJeremy L Thompson 
68*d416dc2bSJeremy L Thompson 
69*d416dc2bSJeremy L Thompson template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename Number>
70*d416dc2bSJeremy L Thompson class OperatorDealiiMassLocal
71*d416dc2bSJeremy L Thompson {
72*d416dc2bSJeremy L Thompson public:
73*d416dc2bSJeremy L Thompson   DEAL_II_HOST_DEVICE void
operator()74*d416dc2bSJeremy L Thompson   operator()(const typename Portable::MatrixFree<dim, Number>::Data *data,
75*d416dc2bSJeremy L Thompson              const Portable::DeviceVector<Number>                   &src,
76*d416dc2bSJeremy L Thompson              Portable::DeviceVector<Number>                         &dst) const
77*d416dc2bSJeremy L Thompson   {
78*d416dc2bSJeremy L Thompson     Portable::FEEvaluation<dim, fe_degree, n_q_points_1d, n_components, Number> fe_eval(data);
79*d416dc2bSJeremy L Thompson     fe_eval.read_dof_values(src);
80*d416dc2bSJeremy L Thompson     fe_eval.evaluate(EvaluationFlags::values);
81*d416dc2bSJeremy L Thompson     fe_eval.apply_for_each_quad_point(
82*d416dc2bSJeremy L Thompson       OperatorDealiiMassQuad<dim, fe_degree, n_q_points_1d, n_components, Number>());
83*d416dc2bSJeremy L Thompson     fe_eval.integrate(EvaluationFlags::values);
84*d416dc2bSJeremy L Thompson     fe_eval.distribute_local_to_global(dst);
85*d416dc2bSJeremy L Thompson   }
86*d416dc2bSJeremy L Thompson 
87*d416dc2bSJeremy L Thompson   static const unsigned int n_local_dofs = Utilities::pow(fe_degree + 1, dim) * n_components;
88*d416dc2bSJeremy L Thompson   static const unsigned int n_q_points   = Utilities::pow(n_q_points_1d, dim);
89*d416dc2bSJeremy L Thompson };
90*d416dc2bSJeremy L Thompson 
91*d416dc2bSJeremy L Thompson 
92*d416dc2bSJeremy L Thompson 
93*d416dc2bSJeremy L Thompson template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename Number>
94*d416dc2bSJeremy L Thompson class OperatorDealiiLaplaceLocal
95*d416dc2bSJeremy L Thompson {
96*d416dc2bSJeremy L Thompson public:
97*d416dc2bSJeremy L Thompson   DEAL_II_HOST_DEVICE void
operator()98*d416dc2bSJeremy L Thompson   operator()(const typename Portable::MatrixFree<dim, Number>::Data *data,
99*d416dc2bSJeremy L Thompson              const Portable::DeviceVector<Number>                   &src,
100*d416dc2bSJeremy L Thompson              Portable::DeviceVector<Number>                         &dst) const
101*d416dc2bSJeremy L Thompson   {
102*d416dc2bSJeremy L Thompson     Portable::FEEvaluation<dim, fe_degree, n_q_points_1d, n_components, Number> fe_eval(data);
103*d416dc2bSJeremy L Thompson     fe_eval.read_dof_values(src);
104*d416dc2bSJeremy L Thompson     fe_eval.evaluate(EvaluationFlags::gradients);
105*d416dc2bSJeremy L Thompson     fe_eval.apply_for_each_quad_point(
106*d416dc2bSJeremy L Thompson       OperatorDealiiLaplaceQuad<dim, fe_degree, n_q_points_1d, n_components, Number>());
107*d416dc2bSJeremy L Thompson     fe_eval.integrate(EvaluationFlags::gradients);
108*d416dc2bSJeremy L Thompson     fe_eval.distribute_local_to_global(dst);
109*d416dc2bSJeremy L Thompson   }
110*d416dc2bSJeremy L Thompson 
111*d416dc2bSJeremy L Thompson   static const unsigned int n_local_dofs = Utilities::pow(fe_degree + 1, dim) * n_components;
112*d416dc2bSJeremy L Thompson   static const unsigned int n_q_points   = Utilities::pow(n_q_points_1d, dim);
113*d416dc2bSJeremy L Thompson };
114*d416dc2bSJeremy L Thompson 
115*d416dc2bSJeremy L Thompson 
116*d416dc2bSJeremy L Thompson 
117*d416dc2bSJeremy L Thompson /**
118*d416dc2bSJeremy L Thompson  * Operator GPU implementation using deal.II.
119*d416dc2bSJeremy L Thompson  */
120*d416dc2bSJeremy L Thompson template <int dim, typename Number>
121*d416dc2bSJeremy L Thompson class OperatorDealii : public OperatorBase<Number, MemorySpace::Default>
122*d416dc2bSJeremy L Thompson {
123*d416dc2bSJeremy L Thompson public:
124*d416dc2bSJeremy L Thompson   using VectorType = typename OperatorBase<Number, MemorySpace::Default>::VectorType;
125*d416dc2bSJeremy L Thompson 
126*d416dc2bSJeremy L Thompson   /**
127*d416dc2bSJeremy L Thompson    * Constructor.
128*d416dc2bSJeremy L Thompson    */
OperatorDealii(const Mapping<dim> & mapping,const DoFHandler<dim> & dof_handler,const AffineConstraints<Number> & constraints,const Quadrature<dim> & quadrature,const BPType & bp)129*d416dc2bSJeremy L Thompson   OperatorDealii(const Mapping<dim>              &mapping,
130*d416dc2bSJeremy L Thompson                  const DoFHandler<dim>           &dof_handler,
131*d416dc2bSJeremy L Thompson                  const AffineConstraints<Number> &constraints,
132*d416dc2bSJeremy L Thompson                  const Quadrature<dim>           &quadrature,
133*d416dc2bSJeremy L Thompson                  const BPType                    &bp)
134*d416dc2bSJeremy L Thompson     : mapping(mapping)
135*d416dc2bSJeremy L Thompson     , dof_handler(dof_handler)
136*d416dc2bSJeremy L Thompson     , constraints(constraints)
137*d416dc2bSJeremy L Thompson     , quadrature(quadrature)
138*d416dc2bSJeremy L Thompson     , bp(bp)
139*d416dc2bSJeremy L Thompson   {
140*d416dc2bSJeremy L Thompson     reinit();
141*d416dc2bSJeremy L Thompson   }
142*d416dc2bSJeremy L Thompson 
143*d416dc2bSJeremy L Thompson   /**
144*d416dc2bSJeremy L Thompson    * Destructor.
145*d416dc2bSJeremy L Thompson    */
146*d416dc2bSJeremy L Thompson   ~OperatorDealii() = default;
147*d416dc2bSJeremy L Thompson 
148*d416dc2bSJeremy L Thompson   /**
149*d416dc2bSJeremy L Thompson    * Initialized internal data structures, particularly, MatrixFree.
150*d416dc2bSJeremy L Thompson    */
151*d416dc2bSJeremy L Thompson   void
reinit()152*d416dc2bSJeremy L Thompson   reinit() override
153*d416dc2bSJeremy L Thompson   {
154*d416dc2bSJeremy L Thompson     // configure MatrixFree
155*d416dc2bSJeremy L Thompson     typename Portable::MatrixFree<dim, Number>::AdditionalData additional_data;
156*d416dc2bSJeremy L Thompson 
157*d416dc2bSJeremy L Thompson     if (bp <= BPType::BP2) // mass matrix
158*d416dc2bSJeremy L Thompson       additional_data.mapping_update_flags = update_JxW_values | update_values;
159*d416dc2bSJeremy L Thompson     else
160*d416dc2bSJeremy L Thompson       additional_data.mapping_update_flags = update_JxW_values | update_gradients;
161*d416dc2bSJeremy L Thompson 
162*d416dc2bSJeremy L Thompson     // create MatrixFree
163*d416dc2bSJeremy L Thompson     AssertThrow(quadrature.is_tensor_product(), ExcNotImplemented());
164*d416dc2bSJeremy L Thompson     matrix_free.reinit(
165*d416dc2bSJeremy L Thompson       mapping, dof_handler, constraints, quadrature.get_tensor_basis()[0], additional_data);
166*d416dc2bSJeremy L Thompson   }
167*d416dc2bSJeremy L Thompson 
168*d416dc2bSJeremy L Thompson   /**
169*d416dc2bSJeremy L Thompson    * Matrix-vector product.
170*d416dc2bSJeremy L Thompson    */
171*d416dc2bSJeremy L Thompson   void
vmult(VectorType & dst,const VectorType & src)172*d416dc2bSJeremy L Thompson   vmult(VectorType &dst, const VectorType &src) const override
173*d416dc2bSJeremy L Thompson   {
174*d416dc2bSJeremy L Thompson     dst = 0.0;
175*d416dc2bSJeremy L Thompson 
176*d416dc2bSJeremy L Thompson     const unsigned int n_components  = dof_handler.get_fe().n_components();
177*d416dc2bSJeremy L Thompson     const unsigned int fe_degree     = dof_handler.get_fe().tensor_degree();
178*d416dc2bSJeremy L Thompson     const unsigned int n_q_points_1d = quadrature.get_tensor_basis()[0].size();
179*d416dc2bSJeremy L Thompson 
180*d416dc2bSJeremy L Thompson     if (n_components == 1 && fe_degree == 1 && n_q_points_1d == 2)
181*d416dc2bSJeremy L Thompson       this->vmult_internal<1, 1, 2>(dst, src);
182*d416dc2bSJeremy L Thompson     else if (n_components == 1 && fe_degree == 2 && n_q_points_1d == 3)
183*d416dc2bSJeremy L Thompson       this->vmult_internal<1, 2, 3>(dst, src);
184*d416dc2bSJeremy L Thompson     else if (n_components == dim && fe_degree == 1 && n_q_points_1d == 2)
185*d416dc2bSJeremy L Thompson       this->vmult_internal<dim, 1, 2>(dst, src);
186*d416dc2bSJeremy L Thompson     else if (n_components == dim && fe_degree == 2 && n_q_points_1d == 3)
187*d416dc2bSJeremy L Thompson       this->vmult_internal<dim, 2, 3>(dst, src);
188*d416dc2bSJeremy L Thompson     else if (n_components == 1 && fe_degree == 1 && n_q_points_1d == 3)
189*d416dc2bSJeremy L Thompson       this->vmult_internal<1, 1, 3>(dst, src);
190*d416dc2bSJeremy L Thompson     else if (n_components == 1 && fe_degree == 2 && n_q_points_1d == 4)
191*d416dc2bSJeremy L Thompson       this->vmult_internal<1, 2, 4>(dst, src);
192*d416dc2bSJeremy L Thompson     else if (n_components == dim && fe_degree == 1 && n_q_points_1d == 3)
193*d416dc2bSJeremy L Thompson       this->vmult_internal<dim, 1, 3>(dst, src);
194*d416dc2bSJeremy L Thompson     else if (n_components == dim && fe_degree == 2 && n_q_points_1d == 4)
195*d416dc2bSJeremy L Thompson       this->vmult_internal<dim, 2, 4>(dst, src);
196*d416dc2bSJeremy L Thompson     else
197*d416dc2bSJeremy L Thompson       AssertThrow(false, ExcInternalError());
198*d416dc2bSJeremy L Thompson 
199*d416dc2bSJeremy L Thompson     matrix_free.copy_constrained_values(src, dst);
200*d416dc2bSJeremy L Thompson   }
201*d416dc2bSJeremy L Thompson 
202*d416dc2bSJeremy L Thompson   /**
203*d416dc2bSJeremy L Thompson    * Initialize vector.
204*d416dc2bSJeremy L Thompson    */
205*d416dc2bSJeremy L Thompson   void
initialize_dof_vector(VectorType & vec)206*d416dc2bSJeremy L Thompson   initialize_dof_vector(VectorType &vec) const override
207*d416dc2bSJeremy L Thompson   {
208*d416dc2bSJeremy L Thompson     matrix_free.initialize_dof_vector(vec);
209*d416dc2bSJeremy L Thompson   }
210*d416dc2bSJeremy L Thompson 
211*d416dc2bSJeremy L Thompson   /**
212*d416dc2bSJeremy L Thompson    * Compute inverse of diagonal.
213*d416dc2bSJeremy L Thompson    */
214*d416dc2bSJeremy L Thompson   void
compute_inverse_diagonal(VectorType & diagonal)215*d416dc2bSJeremy L Thompson   compute_inverse_diagonal(VectorType &diagonal) const override
216*d416dc2bSJeremy L Thompson   {
217*d416dc2bSJeremy L Thompson     this->initialize_dof_vector(diagonal);
218*d416dc2bSJeremy L Thompson 
219*d416dc2bSJeremy L Thompson     const unsigned int n_components  = dof_handler.get_fe().n_components();
220*d416dc2bSJeremy L Thompson     const unsigned int fe_degree     = dof_handler.get_fe().tensor_degree();
221*d416dc2bSJeremy L Thompson     const unsigned int n_q_points_1d = quadrature.get_tensor_basis()[0].size();
222*d416dc2bSJeremy L Thompson 
223*d416dc2bSJeremy L Thompson     if (n_components == 1 && fe_degree == 1 && n_q_points_1d == 2)
224*d416dc2bSJeremy L Thompson       this->compute_inverse_diagonal_internal<1, 1, 2>(diagonal);
225*d416dc2bSJeremy L Thompson     else if (n_components == 1 && fe_degree == 2 && n_q_points_1d == 3)
226*d416dc2bSJeremy L Thompson       this->compute_inverse_diagonal_internal<1, 2, 3>(diagonal);
227*d416dc2bSJeremy L Thompson     else if (n_components == dim && fe_degree == 1 && n_q_points_1d == 2)
228*d416dc2bSJeremy L Thompson       this->compute_inverse_diagonal_internal<dim, 1, 2>(diagonal);
229*d416dc2bSJeremy L Thompson     else if (n_components == dim && fe_degree == 2 && n_q_points_1d == 3)
230*d416dc2bSJeremy L Thompson       this->compute_inverse_diagonal_internal<dim, 2, 3>(diagonal);
231*d416dc2bSJeremy L Thompson     else if (n_components == 1 && fe_degree == 1 && n_q_points_1d == 3)
232*d416dc2bSJeremy L Thompson       this->compute_inverse_diagonal_internal<1, 1, 3>(diagonal);
233*d416dc2bSJeremy L Thompson     else if (n_components == 1 && fe_degree == 2 && n_q_points_1d == 4)
234*d416dc2bSJeremy L Thompson       this->compute_inverse_diagonal_internal<1, 2, 4>(diagonal);
235*d416dc2bSJeremy L Thompson     else if (n_components == dim && fe_degree == 1 && n_q_points_1d == 3)
236*d416dc2bSJeremy L Thompson       this->compute_inverse_diagonal_internal<dim, 1, 3>(diagonal);
237*d416dc2bSJeremy L Thompson     else if (n_components == dim && fe_degree == 2 && n_q_points_1d == 4)
238*d416dc2bSJeremy L Thompson       this->compute_inverse_diagonal_internal<dim, 2, 4>(diagonal);
239*d416dc2bSJeremy L Thompson     else
240*d416dc2bSJeremy L Thompson       AssertThrow(false, ExcInternalError());
241*d416dc2bSJeremy L Thompson   }
242*d416dc2bSJeremy L Thompson 
243*d416dc2bSJeremy L Thompson private:
244*d416dc2bSJeremy L Thompson   /**
245*d416dc2bSJeremy L Thompson    * Templated vmult function.
246*d416dc2bSJeremy L Thompson    */
247*d416dc2bSJeremy L Thompson   template <int n_components, int fe_degree, int n_q_points_1d>
248*d416dc2bSJeremy L Thompson   void
vmult_internal(VectorType & dst,const VectorType & src)249*d416dc2bSJeremy L Thompson   vmult_internal(VectorType &dst, const VectorType &src) const
250*d416dc2bSJeremy L Thompson   {
251*d416dc2bSJeremy L Thompson     if (bp <= BPType::BP2) // mass matrix
252*d416dc2bSJeremy L Thompson       {
253*d416dc2bSJeremy L Thompson         OperatorDealiiMassLocal<dim, fe_degree, n_q_points_1d, n_components, Number> mass_operator;
254*d416dc2bSJeremy L Thompson         matrix_free.cell_loop(mass_operator, src, dst);
255*d416dc2bSJeremy L Thompson       }
256*d416dc2bSJeremy L Thompson     else
257*d416dc2bSJeremy L Thompson       {
258*d416dc2bSJeremy L Thompson         OperatorDealiiLaplaceLocal<dim, fe_degree, n_q_points_1d, n_components, Number>
259*d416dc2bSJeremy L Thompson           local_operator;
260*d416dc2bSJeremy L Thompson         matrix_free.cell_loop(local_operator, src, dst);
261*d416dc2bSJeremy L Thompson       }
262*d416dc2bSJeremy L Thompson   }
263*d416dc2bSJeremy L Thompson 
264*d416dc2bSJeremy L Thompson   /**
265*d416dc2bSJeremy L Thompson    * Templated compute_inverse_diagonal function.
266*d416dc2bSJeremy L Thompson    */
267*d416dc2bSJeremy L Thompson   template <int n_components, int fe_degree, int n_q_points_1d>
268*d416dc2bSJeremy L Thompson   void
compute_inverse_diagonal_internal(VectorType & diagonal)269*d416dc2bSJeremy L Thompson   compute_inverse_diagonal_internal(VectorType &diagonal) const
270*d416dc2bSJeremy L Thompson   {
271*d416dc2bSJeremy L Thompson     if (bp <= BPType::BP2) // mass matrix
272*d416dc2bSJeremy L Thompson       {
273*d416dc2bSJeremy L Thompson         OperatorDealiiMassQuad<dim, fe_degree, n_q_points_1d, n_components, Number> op_quad;
274*d416dc2bSJeremy L Thompson 
275*d416dc2bSJeremy L Thompson         MatrixFreeTools::compute_diagonal<dim, fe_degree, n_q_points_1d, n_components, Number>(
276*d416dc2bSJeremy L Thompson           matrix_free, diagonal, op_quad, EvaluationFlags::values, EvaluationFlags::values);
277*d416dc2bSJeremy L Thompson       }
278*d416dc2bSJeremy L Thompson     else
279*d416dc2bSJeremy L Thompson       {
280*d416dc2bSJeremy L Thompson         OperatorDealiiLaplaceQuad<dim, fe_degree, n_q_points_1d, n_components, Number> op_quad;
281*d416dc2bSJeremy L Thompson 
282*d416dc2bSJeremy L Thompson         MatrixFreeTools::compute_diagonal<dim, fe_degree, n_q_points_1d, n_components, Number>(
283*d416dc2bSJeremy L Thompson           matrix_free, diagonal, op_quad, EvaluationFlags::gradients, EvaluationFlags::gradients);
284*d416dc2bSJeremy L Thompson       }
285*d416dc2bSJeremy L Thompson 
286*d416dc2bSJeremy L Thompson 
287*d416dc2bSJeremy L Thompson     Number *diagonal_ptr = diagonal.get_values();
288*d416dc2bSJeremy L Thompson 
289*d416dc2bSJeremy L Thompson     Kokkos::parallel_for(
290*d416dc2bSJeremy L Thompson       "lethe::invert_vector",
291*d416dc2bSJeremy L Thompson       Kokkos::RangePolicy<MemorySpace::Default::kokkos_space::execution_space>(
292*d416dc2bSJeremy L Thompson         0, diagonal.locally_owned_size()),
293*d416dc2bSJeremy L Thompson       KOKKOS_LAMBDA(int i) { diagonal_ptr[i] = 1.0 / diagonal_ptr[i]; });
294*d416dc2bSJeremy L Thompson   }
295*d416dc2bSJeremy L Thompson 
296*d416dc2bSJeremy L Thompson   /**
297*d416dc2bSJeremy L Thompson    * Mapping object passed to the constructor.
298*d416dc2bSJeremy L Thompson    */
299*d416dc2bSJeremy L Thompson   const Mapping<dim> &mapping;
300*d416dc2bSJeremy L Thompson 
301*d416dc2bSJeremy L Thompson   /**
302*d416dc2bSJeremy L Thompson    * DoFHandler object passed to the constructor.
303*d416dc2bSJeremy L Thompson    */
304*d416dc2bSJeremy L Thompson   const DoFHandler<dim> &dof_handler;
305*d416dc2bSJeremy L Thompson 
306*d416dc2bSJeremy L Thompson   /**
307*d416dc2bSJeremy L Thompson    * Constraints object passed to the constructor.
308*d416dc2bSJeremy L Thompson    */
309*d416dc2bSJeremy L Thompson   const AffineConstraints<Number> &constraints;
310*d416dc2bSJeremy L Thompson 
311*d416dc2bSJeremy L Thompson   /**
312*d416dc2bSJeremy L Thompson    * Quadrature rule object passed to the constructor.
313*d416dc2bSJeremy L Thompson    */
314*d416dc2bSJeremy L Thompson   const Quadrature<dim> &quadrature;
315*d416dc2bSJeremy L Thompson 
316*d416dc2bSJeremy L Thompson   /**
317*d416dc2bSJeremy L Thompson    * Selected BP.
318*d416dc2bSJeremy L Thompson    */
319*d416dc2bSJeremy L Thompson   const BPType bp;
320*d416dc2bSJeremy L Thompson 
321*d416dc2bSJeremy L Thompson   /**
322*d416dc2bSJeremy L Thompson    * MatrixFree object.
323*d416dc2bSJeremy L Thompson    */
324*d416dc2bSJeremy L Thompson   Portable::MatrixFree<dim, Number> matrix_free;
325*d416dc2bSJeremy L Thompson };
326*d416dc2bSJeremy L Thompson 
327*d416dc2bSJeremy L Thompson #endif
328