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