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_ceed_h 20*d416dc2bSJeremy L Thompson # define bps_ceed_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 // libCEED includes 38*d416dc2bSJeremy L Thompson # include <ceed.h> 39*d416dc2bSJeremy L Thompson # include <ceed/backend.h> 40*d416dc2bSJeremy L Thompson 41*d416dc2bSJeremy L Thompson // QFunction source 42*d416dc2bSJeremy L Thompson # include "bps-qfunctions.h" 43*d416dc2bSJeremy L Thompson 44*d416dc2bSJeremy L Thompson using namespace dealii; 45*d416dc2bSJeremy L Thompson 46*d416dc2bSJeremy L Thompson 47*d416dc2bSJeremy L Thompson /** 48*d416dc2bSJeremy L Thompson * Operator implementation using libCEED. 49*d416dc2bSJeremy L Thompson */ 50*d416dc2bSJeremy L Thompson template <int dim, typename Number, typename MemorySpace = MemorySpace::Host> 51*d416dc2bSJeremy L Thompson class OperatorCeed : public OperatorBase<Number, MemorySpace> 52*d416dc2bSJeremy L Thompson { 53*d416dc2bSJeremy L Thompson public: 54*d416dc2bSJeremy L Thompson using VectorType = typename OperatorBase<Number, MemorySpace>::VectorType; 55*d416dc2bSJeremy L Thompson 56*d416dc2bSJeremy L Thompson /** 57*d416dc2bSJeremy L Thompson * Constructor. 58*d416dc2bSJeremy L Thompson */ OperatorCeed(const Mapping<dim> & mapping,const DoFHandler<dim> & dof_handler,const AffineConstraints<Number> & constraints,const Quadrature<dim> & quadrature,const BPType & bp,const std::string & resource)59*d416dc2bSJeremy L Thompson OperatorCeed(const Mapping<dim> &mapping, 60*d416dc2bSJeremy L Thompson const DoFHandler<dim> &dof_handler, 61*d416dc2bSJeremy L Thompson const AffineConstraints<Number> &constraints, 62*d416dc2bSJeremy L Thompson const Quadrature<dim> &quadrature, 63*d416dc2bSJeremy L Thompson const BPType &bp, 64*d416dc2bSJeremy L Thompson const std::string &resource) 65*d416dc2bSJeremy L Thompson : mapping(mapping) 66*d416dc2bSJeremy L Thompson , dof_handler(dof_handler) 67*d416dc2bSJeremy L Thompson , constraints(constraints) 68*d416dc2bSJeremy L Thompson , quadrature(quadrature) 69*d416dc2bSJeremy L Thompson , bp(bp) 70*d416dc2bSJeremy L Thompson , resource(resource) 71*d416dc2bSJeremy L Thompson { 72*d416dc2bSJeremy L Thompson reinit(); 73*d416dc2bSJeremy L Thompson } 74*d416dc2bSJeremy L Thompson 75*d416dc2bSJeremy L Thompson /** 76*d416dc2bSJeremy L Thompson * Destructor. 77*d416dc2bSJeremy L Thompson */ ~OperatorCeed()78*d416dc2bSJeremy L Thompson ~OperatorCeed() 79*d416dc2bSJeremy L Thompson { 80*d416dc2bSJeremy L Thompson CeedVectorDestroy(&src_ceed); 81*d416dc2bSJeremy L Thompson CeedVectorDestroy(&dst_ceed); 82*d416dc2bSJeremy L Thompson CeedOperatorDestroy(&op_apply); 83*d416dc2bSJeremy L Thompson CeedDestroy(&ceed); 84*d416dc2bSJeremy L Thompson } 85*d416dc2bSJeremy L Thompson 86*d416dc2bSJeremy L Thompson /** 87*d416dc2bSJeremy L Thompson * Initialized internal data structures, particularly, libCEED. 88*d416dc2bSJeremy L Thompson */ 89*d416dc2bSJeremy L Thompson void reinit()90*d416dc2bSJeremy L Thompson reinit() override 91*d416dc2bSJeremy L Thompson { 92*d416dc2bSJeremy L Thompson CeedVector metric_data; 93*d416dc2bSJeremy L Thompson CeedBasis sol_basis; 94*d416dc2bSJeremy L Thompson CeedElemRestriction sol_restriction; 95*d416dc2bSJeremy L Thompson CeedElemRestriction metric_data_restriction; 96*d416dc2bSJeremy L Thompson BuildContext build_ctx_data; 97*d416dc2bSJeremy L Thompson CeedQFunctionContext build_ctx; 98*d416dc2bSJeremy L Thompson CeedQFunction qf_apply; 99*d416dc2bSJeremy L Thompson 100*d416dc2bSJeremy L Thompson const auto &tria = dof_handler.get_triangulation(); 101*d416dc2bSJeremy L Thompson const auto &fe = dof_handler.get_fe(); 102*d416dc2bSJeremy L Thompson 103*d416dc2bSJeremy L Thompson const auto n_components = fe.n_components(); 104*d416dc2bSJeremy L Thompson 105*d416dc2bSJeremy L Thompson if (bp == BPType::BP1 || bp == BPType::BP3 || bp == BPType::BP5) 106*d416dc2bSJeremy L Thompson { 107*d416dc2bSJeremy L Thompson AssertThrow(n_components == 1, ExcInternalError()); 108*d416dc2bSJeremy L Thompson } 109*d416dc2bSJeremy L Thompson else 110*d416dc2bSJeremy L Thompson { 111*d416dc2bSJeremy L Thompson AssertThrow(n_components == dim, ExcInternalError()); 112*d416dc2bSJeremy L Thompson } 113*d416dc2bSJeremy L Thompson 114*d416dc2bSJeremy L Thompson // 1) create CEED instance -> "MatrixFree" 115*d416dc2bSJeremy L Thompson const char *ceed_spec = resource.c_str(); 116*d416dc2bSJeremy L Thompson CeedInit(ceed_spec, &ceed); 117*d416dc2bSJeremy L Thompson 118*d416dc2bSJeremy L Thompson // 2) create shape functions -> "ShapeInfo" 119*d416dc2bSJeremy L Thompson const unsigned int fe_degree = fe.tensor_degree(); 120*d416dc2bSJeremy L Thompson const unsigned int n_q_points = quadrature.get_tensor_basis()[0].size(); 121*d416dc2bSJeremy L Thompson { 122*d416dc2bSJeremy L Thompson const dealii::internal::MatrixFreeFunctions::ShapeInfo<double> shape_info(quadrature, fe, 0); 123*d416dc2bSJeremy L Thompson const auto &shape_data = shape_info.get_shape_data(); 124*d416dc2bSJeremy L Thompson std::vector<CeedScalar> q_ref_1d; 125*d416dc2bSJeremy L Thompson for (const auto q : shape_data.quadrature.get_points()) 126*d416dc2bSJeremy L Thompson q_ref_1d.push_back(q(0)); 127*d416dc2bSJeremy L Thompson 128*d416dc2bSJeremy L Thompson // transpose bases for compatibility with restriction 129*d416dc2bSJeremy L Thompson std::vector<CeedScalar> interp_1d(shape_data.shape_values.size()); 130*d416dc2bSJeremy L Thompson std::vector<CeedScalar> grad_1d(shape_data.shape_gradients.size()); 131*d416dc2bSJeremy L Thompson for (unsigned int i = 0; i < n_q_points; ++i) 132*d416dc2bSJeremy L Thompson for (unsigned int j = 0; j < fe_degree + 1; ++j) 133*d416dc2bSJeremy L Thompson { 134*d416dc2bSJeremy L Thompson interp_1d[j + i * (fe_degree + 1)] = shape_data.shape_values[j * n_q_points + i]; 135*d416dc2bSJeremy L Thompson grad_1d[j + i * (fe_degree + 1)] = shape_data.shape_gradients[j * n_q_points + i]; 136*d416dc2bSJeremy L Thompson } 137*d416dc2bSJeremy L Thompson 138*d416dc2bSJeremy L Thompson CeedBasisCreateTensorH1(ceed, 139*d416dc2bSJeremy L Thompson dim, 140*d416dc2bSJeremy L Thompson n_components, 141*d416dc2bSJeremy L Thompson fe_degree + 1, 142*d416dc2bSJeremy L Thompson n_q_points, 143*d416dc2bSJeremy L Thompson interp_1d.data(), 144*d416dc2bSJeremy L Thompson grad_1d.data(), 145*d416dc2bSJeremy L Thompson q_ref_1d.data(), 146*d416dc2bSJeremy L Thompson quadrature.get_tensor_basis()[0].get_weights().data(), 147*d416dc2bSJeremy L Thompson &sol_basis); 148*d416dc2bSJeremy L Thompson } 149*d416dc2bSJeremy L Thompson 150*d416dc2bSJeremy L Thompson // 3) create restriction matrix -> DoFInfo 151*d416dc2bSJeremy L Thompson unsigned int n_local_active_cells = 0; 152*d416dc2bSJeremy L Thompson 153*d416dc2bSJeremy L Thompson for (const auto &cell : dof_handler.active_cell_iterators()) 154*d416dc2bSJeremy L Thompson if (cell->is_locally_owned()) 155*d416dc2bSJeremy L Thompson n_local_active_cells++; 156*d416dc2bSJeremy L Thompson 157*d416dc2bSJeremy L Thompson partitioner = 158*d416dc2bSJeremy L Thompson std::make_shared<Utilities::MPI::Partitioner>(dof_handler.locally_owned_dofs(), 159*d416dc2bSJeremy L Thompson DoFTools::extract_locally_active_dofs( 160*d416dc2bSJeremy L Thompson dof_handler), 161*d416dc2bSJeremy L Thompson dof_handler.get_communicator()); 162*d416dc2bSJeremy L Thompson 163*d416dc2bSJeremy L Thompson std::vector<CeedInt> indices; 164*d416dc2bSJeremy L Thompson indices.reserve(n_local_active_cells * fe.n_dofs_per_cell() / n_components); 165*d416dc2bSJeremy L Thompson 166*d416dc2bSJeremy L Thompson const auto dof_mapping = FETools::lexicographic_to_hierarchic_numbering<dim>(fe_degree); 167*d416dc2bSJeremy L Thompson 168*d416dc2bSJeremy L Thompson std::vector<types::global_dof_index> local_indices(fe.n_dofs_per_cell()); 169*d416dc2bSJeremy L Thompson 170*d416dc2bSJeremy L Thompson for (const auto &cell : dof_handler.active_cell_iterators()) 171*d416dc2bSJeremy L Thompson if (cell->is_locally_owned()) 172*d416dc2bSJeremy L Thompson { 173*d416dc2bSJeremy L Thompson cell->get_dof_indices(local_indices); 174*d416dc2bSJeremy L Thompson 175*d416dc2bSJeremy L Thompson for (const auto i : dof_mapping) 176*d416dc2bSJeremy L Thompson indices.emplace_back( 177*d416dc2bSJeremy L Thompson partitioner->global_to_local(local_indices[fe.component_to_system_index(0, i)])); 178*d416dc2bSJeremy L Thompson } 179*d416dc2bSJeremy L Thompson 180*d416dc2bSJeremy L Thompson CeedElemRestrictionCreate(ceed, 181*d416dc2bSJeremy L Thompson n_local_active_cells, 182*d416dc2bSJeremy L Thompson fe.n_dofs_per_cell() / n_components, 183*d416dc2bSJeremy L Thompson n_components, 184*d416dc2bSJeremy L Thompson 1, 185*d416dc2bSJeremy L Thompson this->extended_local_size(), 186*d416dc2bSJeremy L Thompson CEED_MEM_HOST, 187*d416dc2bSJeremy L Thompson CEED_COPY_VALUES, 188*d416dc2bSJeremy L Thompson indices.data(), 189*d416dc2bSJeremy L Thompson &sol_restriction); 190*d416dc2bSJeremy L Thompson 191*d416dc2bSJeremy L Thompson // 4) create mapping -> MappingInfo 192*d416dc2bSJeremy L Thompson const unsigned int n_components_metric = (bp <= BPType::BP2) ? 1 : (dim * (dim + 1) / 2); 193*d416dc2bSJeremy L Thompson 194*d416dc2bSJeremy L Thompson metric_data_raw = compute_metric_data(ceed, mapping, tria, quadrature, bp); 195*d416dc2bSJeremy L Thompson 196*d416dc2bSJeremy L Thompson strides = {{1, 197*d416dc2bSJeremy L Thompson static_cast<int>(quadrature.size()), 198*d416dc2bSJeremy L Thompson static_cast<int>(quadrature.size() * n_components_metric)}}; 199*d416dc2bSJeremy L Thompson CeedVectorCreate(ceed, metric_data_raw.size(), &metric_data); 200*d416dc2bSJeremy L Thompson CeedVectorSetArray(metric_data, CEED_MEM_HOST, CEED_USE_POINTER, metric_data_raw.data()); 201*d416dc2bSJeremy L Thompson CeedElemRestrictionCreateStrided(ceed, 202*d416dc2bSJeremy L Thompson n_local_active_cells, 203*d416dc2bSJeremy L Thompson quadrature.size(), 204*d416dc2bSJeremy L Thompson n_components_metric, 205*d416dc2bSJeremy L Thompson metric_data_raw.size(), 206*d416dc2bSJeremy L Thompson strides.data(), 207*d416dc2bSJeremy L Thompson &metric_data_restriction); 208*d416dc2bSJeremy L Thompson 209*d416dc2bSJeremy L Thompson build_ctx_data.dim = dim; 210*d416dc2bSJeremy L Thompson build_ctx_data.space_dim = dim; 211*d416dc2bSJeremy L Thompson 212*d416dc2bSJeremy L Thompson CeedQFunctionContextCreate(ceed, &build_ctx); 213*d416dc2bSJeremy L Thompson CeedQFunctionContextSetData( 214*d416dc2bSJeremy L Thompson build_ctx, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(build_ctx_data), &build_ctx_data); 215*d416dc2bSJeremy L Thompson 216*d416dc2bSJeremy L Thompson // 5) create q operation 217*d416dc2bSJeremy L Thompson if (bp == BPType::BP1) 218*d416dc2bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, f_apply_mass, f_apply_mass_loc, &qf_apply); 219*d416dc2bSJeremy L Thompson else if (bp == BPType::BP2) 220*d416dc2bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, f_apply_mass_vec, f_apply_mass_vec_loc, &qf_apply); 221*d416dc2bSJeremy L Thompson else if (bp == BPType::BP3 || bp == BPType::BP5) 222*d416dc2bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, f_apply_poisson, f_apply_poisson_loc, &qf_apply); 223*d416dc2bSJeremy L Thompson else if (bp == BPType::BP4 || bp == BPType::BP6) 224*d416dc2bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, f_apply_poisson_vec, f_apply_poisson_vec_loc, &qf_apply); 225*d416dc2bSJeremy L Thompson else 226*d416dc2bSJeremy L Thompson AssertThrow(false, ExcInternalError()); 227*d416dc2bSJeremy L Thompson 228*d416dc2bSJeremy L Thompson if (bp <= BPType::BP2) 229*d416dc2bSJeremy L Thompson CeedQFunctionAddInput(qf_apply, "u", n_components, CEED_EVAL_INTERP); 230*d416dc2bSJeremy L Thompson else 231*d416dc2bSJeremy L Thompson CeedQFunctionAddInput(qf_apply, "u", dim * n_components, CEED_EVAL_GRAD); 232*d416dc2bSJeremy L Thompson 233*d416dc2bSJeremy L Thompson CeedQFunctionAddInput(qf_apply, "metric data", n_components_metric, CEED_EVAL_NONE); 234*d416dc2bSJeremy L Thompson 235*d416dc2bSJeremy L Thompson if (bp <= BPType::BP2) 236*d416dc2bSJeremy L Thompson CeedQFunctionAddOutput(qf_apply, "v", n_components, CEED_EVAL_INTERP); 237*d416dc2bSJeremy L Thompson else 238*d416dc2bSJeremy L Thompson CeedQFunctionAddOutput(qf_apply, "v", dim * n_components, CEED_EVAL_GRAD); 239*d416dc2bSJeremy L Thompson 240*d416dc2bSJeremy L Thompson CeedQFunctionSetContext(qf_apply, build_ctx); 241*d416dc2bSJeremy L Thompson 242*d416dc2bSJeremy L Thompson // 6) put everything together 243*d416dc2bSJeremy L Thompson CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_apply); 244*d416dc2bSJeremy L Thompson 245*d416dc2bSJeremy L Thompson CeedOperatorSetField(op_apply, "u", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE); 246*d416dc2bSJeremy L Thompson CeedOperatorSetField( 247*d416dc2bSJeremy L Thompson op_apply, "metric data", metric_data_restriction, CEED_BASIS_NONE, metric_data); 248*d416dc2bSJeremy L Thompson CeedOperatorSetField(op_apply, "v", sol_restriction, sol_basis, CEED_VECTOR_ACTIVE); 249*d416dc2bSJeremy L Thompson 250*d416dc2bSJeremy L Thompson // 7) libCEED vectors 251*d416dc2bSJeremy L Thompson CeedElemRestrictionCreateVector(sol_restriction, &src_ceed, NULL); 252*d416dc2bSJeremy L Thompson CeedElemRestrictionCreateVector(sol_restriction, &dst_ceed, NULL); 253*d416dc2bSJeremy L Thompson 254*d416dc2bSJeremy L Thompson // 8) cleanup 255*d416dc2bSJeremy L Thompson CeedVectorDestroy(&metric_data); 256*d416dc2bSJeremy L Thompson CeedElemRestrictionDestroy(&metric_data_restriction); 257*d416dc2bSJeremy L Thompson CeedElemRestrictionDestroy(&sol_restriction); 258*d416dc2bSJeremy L Thompson CeedBasisDestroy(&sol_basis); 259*d416dc2bSJeremy L Thompson CeedQFunctionContextDestroy(&build_ctx); 260*d416dc2bSJeremy L Thompson CeedQFunctionDestroy(&qf_apply); 261*d416dc2bSJeremy L Thompson } 262*d416dc2bSJeremy L Thompson 263*d416dc2bSJeremy L Thompson /** 264*d416dc2bSJeremy L Thompson * Perform matrix-vector product. 265*d416dc2bSJeremy L Thompson */ 266*d416dc2bSJeremy L Thompson void vmult(VectorType & dst,const VectorType & src)267*d416dc2bSJeremy L Thompson vmult(VectorType &dst, const VectorType &src) const override 268*d416dc2bSJeremy L Thompson { 269*d416dc2bSJeremy L Thompson // communicate: update ghost values 270*d416dc2bSJeremy L Thompson src.update_ghost_values(); 271*d416dc2bSJeremy L Thompson 272*d416dc2bSJeremy L Thompson // pass memory buffers to libCEED 273*d416dc2bSJeremy L Thompson VectorTypeCeed x(src_ceed); 274*d416dc2bSJeremy L Thompson VectorTypeCeed y(dst_ceed); 275*d416dc2bSJeremy L Thompson x.import_array(src, CEED_MEM_HOST); 276*d416dc2bSJeremy L Thompson y.import_array(dst, CEED_MEM_HOST); 277*d416dc2bSJeremy L Thompson 278*d416dc2bSJeremy L Thompson // apply operator 279*d416dc2bSJeremy L Thompson CeedOperatorApply(op_apply, x(), y(), CEED_REQUEST_IMMEDIATE); 280*d416dc2bSJeremy L Thompson 281*d416dc2bSJeremy L Thompson // pull arrays back to deal.II 282*d416dc2bSJeremy L Thompson x.take_array(); 283*d416dc2bSJeremy L Thompson y.take_array(); 284*d416dc2bSJeremy L Thompson 285*d416dc2bSJeremy L Thompson // communicate: compress 286*d416dc2bSJeremy L Thompson src.zero_out_ghost_values(); 287*d416dc2bSJeremy L Thompson dst.compress(VectorOperation::add); 288*d416dc2bSJeremy L Thompson 289*d416dc2bSJeremy L Thompson // apply constraints: we assume homogeneous DBC 290*d416dc2bSJeremy L Thompson constraints.set_zero(dst); 291*d416dc2bSJeremy L Thompson } 292*d416dc2bSJeremy L Thompson 293*d416dc2bSJeremy L Thompson /** 294*d416dc2bSJeremy L Thompson * Initialized vector. 295*d416dc2bSJeremy L Thompson */ 296*d416dc2bSJeremy L Thompson void initialize_dof_vector(VectorType & vec)297*d416dc2bSJeremy L Thompson initialize_dof_vector(VectorType &vec) const override 298*d416dc2bSJeremy L Thompson { 299*d416dc2bSJeremy L Thompson vec.reinit(partitioner); 300*d416dc2bSJeremy L Thompson } 301*d416dc2bSJeremy L Thompson 302*d416dc2bSJeremy L Thompson /** 303*d416dc2bSJeremy L Thompson * Compute inverse of diagonal. 304*d416dc2bSJeremy L Thompson */ 305*d416dc2bSJeremy L Thompson void compute_inverse_diagonal(VectorType & diagonal)306*d416dc2bSJeremy L Thompson compute_inverse_diagonal(VectorType &diagonal) const override 307*d416dc2bSJeremy L Thompson { 308*d416dc2bSJeremy L Thompson this->initialize_dof_vector(diagonal); 309*d416dc2bSJeremy L Thompson 310*d416dc2bSJeremy L Thompson // pass memory buffer to libCEED 311*d416dc2bSJeremy L Thompson VectorTypeCeed y(dst_ceed); 312*d416dc2bSJeremy L Thompson y.import_array(diagonal, CEED_MEM_HOST); 313*d416dc2bSJeremy L Thompson 314*d416dc2bSJeremy L Thompson CeedOperatorLinearAssembleDiagonal(op_apply, y(), CEED_REQUEST_IMMEDIATE); 315*d416dc2bSJeremy L Thompson 316*d416dc2bSJeremy L Thompson // pull array back to deal.II 317*d416dc2bSJeremy L Thompson y.take_array(); 318*d416dc2bSJeremy L Thompson 319*d416dc2bSJeremy L Thompson diagonal.compress(VectorOperation::add); 320*d416dc2bSJeremy L Thompson 321*d416dc2bSJeremy L Thompson // apply constraints: we assume homogeneous DBC 322*d416dc2bSJeremy L Thompson constraints.set_zero(diagonal); 323*d416dc2bSJeremy L Thompson 324*d416dc2bSJeremy L Thompson for (auto &i : diagonal) 325*d416dc2bSJeremy L Thompson i = (std::abs(i) > 1.0e-10) ? (1.0 / i) : 1.0; 326*d416dc2bSJeremy L Thompson } 327*d416dc2bSJeremy L Thompson 328*d416dc2bSJeremy L Thompson private: 329*d416dc2bSJeremy L Thompson /** 330*d416dc2bSJeremy L Thompson * Wrapper around a deal.II vector to create a libCEED vector view. 331*d416dc2bSJeremy L Thompson */ 332*d416dc2bSJeremy L Thompson class VectorTypeCeed 333*d416dc2bSJeremy L Thompson { 334*d416dc2bSJeremy L Thompson public: 335*d416dc2bSJeremy L Thompson /** 336*d416dc2bSJeremy L Thompson * Constructor. 337*d416dc2bSJeremy L Thompson */ VectorTypeCeed(const CeedVector & vec_orig)338*d416dc2bSJeremy L Thompson VectorTypeCeed(const CeedVector &vec_orig) 339*d416dc2bSJeremy L Thompson { 340*d416dc2bSJeremy L Thompson vec_ceed = NULL; 341*d416dc2bSJeremy L Thompson CeedVectorReferenceCopy(vec_orig, &vec_ceed); 342*d416dc2bSJeremy L Thompson } 343*d416dc2bSJeremy L Thompson 344*d416dc2bSJeremy L Thompson /** 345*d416dc2bSJeremy L Thompson * Return libCEED vector view. 346*d416dc2bSJeremy L Thompson */ 347*d416dc2bSJeremy L Thompson CeedVector & operator()348*d416dc2bSJeremy L Thompson operator()() 349*d416dc2bSJeremy L Thompson { 350*d416dc2bSJeremy L Thompson return vec_ceed; 351*d416dc2bSJeremy L Thompson } 352*d416dc2bSJeremy L Thompson 353*d416dc2bSJeremy L Thompson /** 354*d416dc2bSJeremy L Thompson * Set deal.II memory in libCEED vector. 355*d416dc2bSJeremy L Thompson */ 356*d416dc2bSJeremy L Thompson void import_array(const VectorType & vec,const CeedMemType space)357*d416dc2bSJeremy L Thompson import_array(const VectorType &vec, const CeedMemType space) 358*d416dc2bSJeremy L Thompson { 359*d416dc2bSJeremy L Thompson mem_space = space; 360*d416dc2bSJeremy L Thompson CeedVectorSetArray(vec_ceed, mem_space, CEED_USE_POINTER, vec.get_values()); 361*d416dc2bSJeremy L Thompson } 362*d416dc2bSJeremy L Thompson 363*d416dc2bSJeremy L Thompson /** 364*d416dc2bSJeremy L Thompson * Sync memory from device to host. 365*d416dc2bSJeremy L Thompson */ 366*d416dc2bSJeremy L Thompson void sync_array()367*d416dc2bSJeremy L Thompson sync_array() 368*d416dc2bSJeremy L Thompson { 369*d416dc2bSJeremy L Thompson CeedVectorSyncArray(vec_ceed, mem_space); 370*d416dc2bSJeremy L Thompson } 371*d416dc2bSJeremy L Thompson 372*d416dc2bSJeremy L Thompson /** 373*d416dc2bSJeremy L Thompson * Take previously set deal.II array from libCEED vector 374*d416dc2bSJeremy L Thompson */ 375*d416dc2bSJeremy L Thompson void take_array()376*d416dc2bSJeremy L Thompson take_array() 377*d416dc2bSJeremy L Thompson { 378*d416dc2bSJeremy L Thompson CeedScalar *ptr; 379*d416dc2bSJeremy L Thompson CeedVectorTakeArray(vec_ceed, mem_space, &ptr); 380*d416dc2bSJeremy L Thompson } 381*d416dc2bSJeremy L Thompson 382*d416dc2bSJeremy L Thompson /** 383*d416dc2bSJeremy L Thompson * Destructor: destroy vector view. 384*d416dc2bSJeremy L Thompson */ ~VectorTypeCeed()385*d416dc2bSJeremy L Thompson ~VectorTypeCeed() 386*d416dc2bSJeremy L Thompson { 387*d416dc2bSJeremy L Thompson bool has_array; 388*d416dc2bSJeremy L Thompson CeedVectorHasBorrowedArrayOfType(vec_ceed, mem_space, &has_array); 389*d416dc2bSJeremy L Thompson if (has_array) 390*d416dc2bSJeremy L Thompson { 391*d416dc2bSJeremy L Thompson CeedScalar *ptr; 392*d416dc2bSJeremy L Thompson CeedVectorTakeArray(vec_ceed, mem_space, &ptr); 393*d416dc2bSJeremy L Thompson } 394*d416dc2bSJeremy L Thompson CeedVectorDestroy(&vec_ceed); 395*d416dc2bSJeremy L Thompson } 396*d416dc2bSJeremy L Thompson 397*d416dc2bSJeremy L Thompson private: 398*d416dc2bSJeremy L Thompson /** 399*d416dc2bSJeremy L Thompson * libCEED vector view. 400*d416dc2bSJeremy L Thompson */ 401*d416dc2bSJeremy L Thompson CeedMemType mem_space; 402*d416dc2bSJeremy L Thompson CeedVector vec_ceed; 403*d416dc2bSJeremy L Thompson }; 404*d416dc2bSJeremy L Thompson 405*d416dc2bSJeremy L Thompson /** 406*d416dc2bSJeremy L Thompson * Number of locally active DoFs. 407*d416dc2bSJeremy L Thompson */ 408*d416dc2bSJeremy L Thompson unsigned int extended_local_size()409*d416dc2bSJeremy L Thompson extended_local_size() const 410*d416dc2bSJeremy L Thompson { 411*d416dc2bSJeremy L Thompson return partitioner->locally_owned_size() + partitioner->n_ghost_indices(); 412*d416dc2bSJeremy L Thompson } 413*d416dc2bSJeremy L Thompson 414*d416dc2bSJeremy L Thompson /** 415*d416dc2bSJeremy L Thompson * Compute metric data: Jacobian, ... 416*d416dc2bSJeremy L Thompson */ 417*d416dc2bSJeremy L Thompson static std::vector<double> compute_metric_data(const Ceed & ceed,const Mapping<dim> & mapping,const Triangulation<dim> & tria,const Quadrature<dim> & quadrature,const BPType bp)418*d416dc2bSJeremy L Thompson compute_metric_data(const Ceed &ceed, 419*d416dc2bSJeremy L Thompson const Mapping<dim> &mapping, 420*d416dc2bSJeremy L Thompson const Triangulation<dim> &tria, 421*d416dc2bSJeremy L Thompson const Quadrature<dim> &quadrature, 422*d416dc2bSJeremy L Thompson const BPType bp) 423*d416dc2bSJeremy L Thompson { 424*d416dc2bSJeremy L Thompson std::vector<double> metric_data_raw; 425*d416dc2bSJeremy L Thompson 426*d416dc2bSJeremy L Thompson CeedBasis geo_basis; 427*d416dc2bSJeremy L Thompson CeedVector metric_data; 428*d416dc2bSJeremy L Thompson CeedElemRestriction metric_data_restriction; 429*d416dc2bSJeremy L Thompson CeedVector node_coords; 430*d416dc2bSJeremy L Thompson CeedElemRestriction geo_restriction; 431*d416dc2bSJeremy L Thompson CeedQFunctionContext build_ctx; 432*d416dc2bSJeremy L Thompson CeedQFunction qf_build; 433*d416dc2bSJeremy L Thompson CeedOperator op_build; 434*d416dc2bSJeremy L Thompson 435*d416dc2bSJeremy L Thompson const unsigned int n_q_points = quadrature.get_tensor_basis()[0].size(); 436*d416dc2bSJeremy L Thompson 437*d416dc2bSJeremy L Thompson const unsigned int n_components_metric = (bp <= BPType::BP2) ? 1 : (dim * (dim + 1) / 2); 438*d416dc2bSJeremy L Thompson 439*d416dc2bSJeremy L Thompson const auto mapping_q = dynamic_cast<const MappingQ<dim> *>(&mapping); 440*d416dc2bSJeremy L Thompson 441*d416dc2bSJeremy L Thompson AssertThrow(mapping_q, ExcMessage("Wrong mapping!")); 442*d416dc2bSJeremy L Thompson 443*d416dc2bSJeremy L Thompson const unsigned int fe_degree = mapping_q->get_degree(); 444*d416dc2bSJeremy L Thompson 445*d416dc2bSJeremy L Thompson FE_Q<dim> geo_fe(fe_degree); 446*d416dc2bSJeremy L Thompson 447*d416dc2bSJeremy L Thompson { 448*d416dc2bSJeremy L Thompson const dealii::internal::MatrixFreeFunctions::ShapeInfo<double> shape_info(quadrature, 449*d416dc2bSJeremy L Thompson geo_fe, 450*d416dc2bSJeremy L Thompson 0); 451*d416dc2bSJeremy L Thompson const auto &shape_data = shape_info.get_shape_data(); 452*d416dc2bSJeremy L Thompson std::vector<CeedScalar> q_ref_1d; 453*d416dc2bSJeremy L Thompson for (const auto q : shape_data.quadrature.get_points()) 454*d416dc2bSJeremy L Thompson q_ref_1d.push_back(q(0)); 455*d416dc2bSJeremy L Thompson 456*d416dc2bSJeremy L Thompson // transpose bases for compatibility with restriction 457*d416dc2bSJeremy L Thompson std::vector<CeedScalar> interp_1d(shape_data.shape_values.size()); 458*d416dc2bSJeremy L Thompson std::vector<CeedScalar> grad_1d(shape_data.shape_gradients.size()); 459*d416dc2bSJeremy L Thompson for (unsigned int i = 0; i < n_q_points; ++i) 460*d416dc2bSJeremy L Thompson for (unsigned int j = 0; j < fe_degree + 1; ++j) 461*d416dc2bSJeremy L Thompson { 462*d416dc2bSJeremy L Thompson interp_1d[j + i * (fe_degree + 1)] = shape_data.shape_values[j * n_q_points + i]; 463*d416dc2bSJeremy L Thompson grad_1d[j + i * (fe_degree + 1)] = shape_data.shape_gradients[j * n_q_points + i]; 464*d416dc2bSJeremy L Thompson } 465*d416dc2bSJeremy L Thompson 466*d416dc2bSJeremy L Thompson CeedBasisCreateTensorH1(ceed, 467*d416dc2bSJeremy L Thompson dim, 468*d416dc2bSJeremy L Thompson dim, 469*d416dc2bSJeremy L Thompson fe_degree + 1, 470*d416dc2bSJeremy L Thompson n_q_points, 471*d416dc2bSJeremy L Thompson interp_1d.data(), 472*d416dc2bSJeremy L Thompson grad_1d.data(), 473*d416dc2bSJeremy L Thompson q_ref_1d.data(), 474*d416dc2bSJeremy L Thompson quadrature.get_tensor_basis()[0].get_weights().data(), 475*d416dc2bSJeremy L Thompson &geo_basis); 476*d416dc2bSJeremy L Thompson } 477*d416dc2bSJeremy L Thompson 478*d416dc2bSJeremy L Thompson unsigned int n_local_active_cells = 0; 479*d416dc2bSJeremy L Thompson 480*d416dc2bSJeremy L Thompson for (const auto &cell : tria.active_cell_iterators()) 481*d416dc2bSJeremy L Thompson if (cell->is_locally_owned()) 482*d416dc2bSJeremy L Thompson n_local_active_cells++; 483*d416dc2bSJeremy L Thompson 484*d416dc2bSJeremy L Thompson std::vector<double> geo_support_points; 485*d416dc2bSJeremy L Thompson std::vector<CeedInt> geo_indices; 486*d416dc2bSJeremy L Thompson 487*d416dc2bSJeremy L Thompson DoFHandler<dim> geo_dof_handler(tria); 488*d416dc2bSJeremy L Thompson geo_dof_handler.distribute_dofs(geo_fe); 489*d416dc2bSJeremy L Thompson 490*d416dc2bSJeremy L Thompson const auto geo_partitioner = 491*d416dc2bSJeremy L Thompson std::make_shared<Utilities::MPI::Partitioner>(geo_dof_handler.locally_owned_dofs(), 492*d416dc2bSJeremy L Thompson DoFTools::extract_locally_active_dofs( 493*d416dc2bSJeremy L Thompson geo_dof_handler), 494*d416dc2bSJeremy L Thompson geo_dof_handler.get_communicator()); 495*d416dc2bSJeremy L Thompson 496*d416dc2bSJeremy L Thompson geo_indices.reserve(n_local_active_cells * geo_fe.n_dofs_per_cell()); 497*d416dc2bSJeremy L Thompson 498*d416dc2bSJeremy L Thompson const auto dof_mapping = FETools::lexicographic_to_hierarchic_numbering<dim>(fe_degree); 499*d416dc2bSJeremy L Thompson 500*d416dc2bSJeremy L Thompson FEValues<dim> fe_values(mapping, 501*d416dc2bSJeremy L Thompson geo_fe, 502*d416dc2bSJeremy L Thompson geo_fe.get_unit_support_points(), 503*d416dc2bSJeremy L Thompson update_quadrature_points); 504*d416dc2bSJeremy L Thompson 505*d416dc2bSJeremy L Thompson std::vector<types::global_dof_index> local_indices(geo_fe.n_dofs_per_cell()); 506*d416dc2bSJeremy L Thompson 507*d416dc2bSJeremy L Thompson const unsigned int n_points = 508*d416dc2bSJeremy L Thompson geo_partitioner->locally_owned_size() + geo_partitioner->n_ghost_indices(); 509*d416dc2bSJeremy L Thompson 510*d416dc2bSJeremy L Thompson geo_support_points.resize(dim * n_points); 511*d416dc2bSJeremy L Thompson 512*d416dc2bSJeremy L Thompson for (const auto &cell : geo_dof_handler.active_cell_iterators()) 513*d416dc2bSJeremy L Thompson if (cell->is_locally_owned()) 514*d416dc2bSJeremy L Thompson { 515*d416dc2bSJeremy L Thompson fe_values.reinit(cell); 516*d416dc2bSJeremy L Thompson cell->get_dof_indices(local_indices); 517*d416dc2bSJeremy L Thompson 518*d416dc2bSJeremy L Thompson for (const auto i : dof_mapping) 519*d416dc2bSJeremy L Thompson { 520*d416dc2bSJeremy L Thompson const auto index = geo_partitioner->global_to_local(local_indices[i]); 521*d416dc2bSJeremy L Thompson geo_indices.emplace_back(index * dim); 522*d416dc2bSJeremy L Thompson 523*d416dc2bSJeremy L Thompson const auto point = fe_values.quadrature_point(i); 524*d416dc2bSJeremy L Thompson 525*d416dc2bSJeremy L Thompson for (unsigned int d = 0; d < dim; ++d) 526*d416dc2bSJeremy L Thompson geo_support_points[index * dim + d] = point[d]; 527*d416dc2bSJeremy L Thompson } 528*d416dc2bSJeremy L Thompson } 529*d416dc2bSJeremy L Thompson 530*d416dc2bSJeremy L Thompson metric_data_raw.resize(n_local_active_cells * quadrature.size() * n_components_metric); 531*d416dc2bSJeremy L Thompson 532*d416dc2bSJeremy L Thompson CeedInt strides[3] = {1, 533*d416dc2bSJeremy L Thompson static_cast<int>(quadrature.size()), 534*d416dc2bSJeremy L Thompson static_cast<int>(quadrature.size() * n_components_metric)}; 535*d416dc2bSJeremy L Thompson 536*d416dc2bSJeremy L Thompson CeedVectorCreate(ceed, metric_data_raw.size(), &metric_data); 537*d416dc2bSJeremy L Thompson CeedVectorSetArray(metric_data, CEED_MEM_HOST, CEED_USE_POINTER, metric_data_raw.data()); 538*d416dc2bSJeremy L Thompson CeedElemRestrictionCreateStrided(ceed, 539*d416dc2bSJeremy L Thompson n_local_active_cells, 540*d416dc2bSJeremy L Thompson quadrature.size(), 541*d416dc2bSJeremy L Thompson n_components_metric, 542*d416dc2bSJeremy L Thompson metric_data_raw.size(), 543*d416dc2bSJeremy L Thompson strides, 544*d416dc2bSJeremy L Thompson &metric_data_restriction); 545*d416dc2bSJeremy L Thompson 546*d416dc2bSJeremy L Thompson CeedVectorCreate(ceed, geo_support_points.size(), &node_coords); 547*d416dc2bSJeremy L Thompson CeedVectorSetArray(node_coords, CEED_MEM_HOST, CEED_USE_POINTER, geo_support_points.data()); 548*d416dc2bSJeremy L Thompson 549*d416dc2bSJeremy L Thompson CeedElemRestrictionCreate(ceed, 550*d416dc2bSJeremy L Thompson n_local_active_cells, 551*d416dc2bSJeremy L Thompson geo_fe.n_dofs_per_cell(), 552*d416dc2bSJeremy L Thompson dim, 553*d416dc2bSJeremy L Thompson 1, 554*d416dc2bSJeremy L Thompson geo_support_points.size(), 555*d416dc2bSJeremy L Thompson CEED_MEM_HOST, 556*d416dc2bSJeremy L Thompson CEED_COPY_VALUES, 557*d416dc2bSJeremy L Thompson geo_indices.data(), 558*d416dc2bSJeremy L Thompson &geo_restriction); 559*d416dc2bSJeremy L Thompson 560*d416dc2bSJeremy L Thompson BuildContext build_ctx_data; 561*d416dc2bSJeremy L Thompson build_ctx_data.dim = dim; 562*d416dc2bSJeremy L Thompson build_ctx_data.space_dim = dim; 563*d416dc2bSJeremy L Thompson 564*d416dc2bSJeremy L Thompson CeedQFunctionContextCreate(ceed, &build_ctx); 565*d416dc2bSJeremy L Thompson CeedQFunctionContextSetData( 566*d416dc2bSJeremy L Thompson build_ctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(build_ctx_data), &build_ctx_data); 567*d416dc2bSJeremy L Thompson 568*d416dc2bSJeremy L Thompson // 5) create q operation 569*d416dc2bSJeremy L Thompson if (bp <= BPType::BP2) 570*d416dc2bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, f_build_mass, f_build_mass_loc, &qf_build); 571*d416dc2bSJeremy L Thompson else 572*d416dc2bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, f_build_poisson, f_build_poisson_loc, &qf_build); 573*d416dc2bSJeremy L Thompson 574*d416dc2bSJeremy L Thompson CeedQFunctionAddInput(qf_build, "geo", dim * dim, CEED_EVAL_GRAD); 575*d416dc2bSJeremy L Thompson CeedQFunctionAddInput(qf_build, "weight", 1, CEED_EVAL_WEIGHT); 576*d416dc2bSJeremy L Thompson CeedQFunctionAddOutput(qf_build, "metric data", n_components_metric, CEED_EVAL_NONE); 577*d416dc2bSJeremy L Thompson CeedQFunctionSetContext(qf_build, build_ctx); 578*d416dc2bSJeremy L Thompson 579*d416dc2bSJeremy L Thompson // 6) put everything together 580*d416dc2bSJeremy L Thompson CeedOperatorCreate(ceed, qf_build, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_build); 581*d416dc2bSJeremy L Thompson CeedOperatorSetField(op_build, "geo", geo_restriction, geo_basis, CEED_VECTOR_ACTIVE); 582*d416dc2bSJeremy L Thompson CeedOperatorSetField( 583*d416dc2bSJeremy L Thompson op_build, "weight", CEED_ELEMRESTRICTION_NONE, geo_basis, CEED_VECTOR_NONE); 584*d416dc2bSJeremy L Thompson CeedOperatorSetField( 585*d416dc2bSJeremy L Thompson op_build, "metric data", metric_data_restriction, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE); 586*d416dc2bSJeremy L Thompson 587*d416dc2bSJeremy L Thompson CeedOperatorApply(op_build, node_coords, metric_data, CEED_REQUEST_IMMEDIATE); 588*d416dc2bSJeremy L Thompson 589*d416dc2bSJeremy L Thompson CeedVectorDestroy(&node_coords); 590*d416dc2bSJeremy L Thompson CeedVectorSyncArray(metric_data, CEED_MEM_HOST); 591*d416dc2bSJeremy L Thompson CeedVectorDestroy(&metric_data); 592*d416dc2bSJeremy L Thompson CeedElemRestrictionDestroy(&geo_restriction); 593*d416dc2bSJeremy L Thompson CeedElemRestrictionDestroy(&metric_data_restriction); 594*d416dc2bSJeremy L Thompson CeedBasisDestroy(&geo_basis); 595*d416dc2bSJeremy L Thompson CeedQFunctionContextDestroy(&build_ctx); 596*d416dc2bSJeremy L Thompson CeedQFunctionDestroy(&qf_build); 597*d416dc2bSJeremy L Thompson CeedOperatorDestroy(&op_build); 598*d416dc2bSJeremy L Thompson 599*d416dc2bSJeremy L Thompson return metric_data_raw; 600*d416dc2bSJeremy L Thompson } 601*d416dc2bSJeremy L Thompson 602*d416dc2bSJeremy L Thompson /** 603*d416dc2bSJeremy L Thompson * Mapping object passed to the constructor. 604*d416dc2bSJeremy L Thompson */ 605*d416dc2bSJeremy L Thompson const Mapping<dim> &mapping; 606*d416dc2bSJeremy L Thompson 607*d416dc2bSJeremy L Thompson /** 608*d416dc2bSJeremy L Thompson * DoFHandler object passed to the constructor. 609*d416dc2bSJeremy L Thompson */ 610*d416dc2bSJeremy L Thompson const DoFHandler<dim> &dof_handler; 611*d416dc2bSJeremy L Thompson 612*d416dc2bSJeremy L Thompson /** 613*d416dc2bSJeremy L Thompson * Constraints object passed to the constructor. 614*d416dc2bSJeremy L Thompson */ 615*d416dc2bSJeremy L Thompson const AffineConstraints<Number> &constraints; 616*d416dc2bSJeremy L Thompson 617*d416dc2bSJeremy L Thompson /** 618*d416dc2bSJeremy L Thompson * Quadrature rule object passed to the constructor. 619*d416dc2bSJeremy L Thompson */ 620*d416dc2bSJeremy L Thompson const Quadrature<dim> &quadrature; 621*d416dc2bSJeremy L Thompson 622*d416dc2bSJeremy L Thompson /** 623*d416dc2bSJeremy L Thompson * Selected BP. 624*d416dc2bSJeremy L Thompson */ 625*d416dc2bSJeremy L Thompson const BPType bp; 626*d416dc2bSJeremy L Thompson 627*d416dc2bSJeremy L Thompson /** 628*d416dc2bSJeremy L Thompson * Resource name. 629*d416dc2bSJeremy L Thompson */ 630*d416dc2bSJeremy L Thompson const std::string resource; 631*d416dc2bSJeremy L Thompson 632*d416dc2bSJeremy L Thompson /** 633*d416dc2bSJeremy L Thompson * Partitioner for distributed vectors. 634*d416dc2bSJeremy L Thompson */ 635*d416dc2bSJeremy L Thompson std::shared_ptr<Utilities::MPI::Partitioner> partitioner; 636*d416dc2bSJeremy L Thompson 637*d416dc2bSJeremy L Thompson /** 638*d416dc2bSJeremy L Thompson * libCEED data structures. 639*d416dc2bSJeremy L Thompson */ 640*d416dc2bSJeremy L Thompson Ceed ceed; 641*d416dc2bSJeremy L Thompson std::vector<double> metric_data_raw; 642*d416dc2bSJeremy L Thompson std::array<CeedInt, 3> strides; 643*d416dc2bSJeremy L Thompson CeedVector src_ceed; 644*d416dc2bSJeremy L Thompson CeedVector dst_ceed; 645*d416dc2bSJeremy L Thompson CeedOperator op_apply; 646*d416dc2bSJeremy L Thompson }; 647*d416dc2bSJeremy L Thompson 648*d416dc2bSJeremy L Thompson #endif 649