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