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