1*d416dc2bSJeremy L Thompson // --------------------------------------------------------------------- 2*d416dc2bSJeremy L Thompson // 3*d416dc2bSJeremy L Thompson // Copyright (C) 2023 by the deal.II authors 4*d416dc2bSJeremy L Thompson // 5*d416dc2bSJeremy L Thompson // This file is part of the deal.II library. 6*d416dc2bSJeremy L Thompson // 7*d416dc2bSJeremy L Thompson // The deal.II library is free software; you can use it, redistribute 8*d416dc2bSJeremy L Thompson // it, and/or modify it under the terms of the GNU Lesser General 9*d416dc2bSJeremy L Thompson // Public License as published by the Free Software Foundation; either 10*d416dc2bSJeremy L Thompson // version 2.1 of the License, or (at your option) any later version. 11*d416dc2bSJeremy L Thompson // The full text of the license can be found in the file LICENSE.md at 12*d416dc2bSJeremy L Thompson // the top level directory of deal.II. 13*d416dc2bSJeremy L Thompson // 14*d416dc2bSJeremy L Thompson // Authors: Peter Munch, Martin Kronbichler 15*d416dc2bSJeremy L Thompson // 16*d416dc2bSJeremy L Thompson // --------------------------------------------------------------------- 17*d416dc2bSJeremy L Thompson 18*d416dc2bSJeremy L Thompson #pragma once 19*d416dc2bSJeremy L Thompson #ifndef bps_cpu_h 20*d416dc2bSJeremy L Thompson # define bps_cpu_h 21*d416dc2bSJeremy L Thompson 22*d416dc2bSJeremy L Thompson // deal.II includes 23*d416dc2bSJeremy L Thompson # include <deal.II/dofs/dof_tools.h> 24*d416dc2bSJeremy L Thompson 25*d416dc2bSJeremy L Thompson # include <deal.II/fe/mapping.h> 26*d416dc2bSJeremy L Thompson 27*d416dc2bSJeremy L Thompson # include <deal.II/lac/la_parallel_vector.h> 28*d416dc2bSJeremy L Thompson 29*d416dc2bSJeremy L Thompson # include <deal.II/matrix_free/fe_evaluation.h> 30*d416dc2bSJeremy L Thompson # include <deal.II/matrix_free/matrix_free.h> 31*d416dc2bSJeremy L Thompson # include <deal.II/matrix_free/shape_info.h> 32*d416dc2bSJeremy L Thompson # include <deal.II/matrix_free/tools.h> 33*d416dc2bSJeremy L Thompson 34*d416dc2bSJeremy L Thompson // local includes 35*d416dc2bSJeremy L Thompson # include "bps.h" 36*d416dc2bSJeremy L Thompson 37*d416dc2bSJeremy L Thompson using namespace dealii; 38*d416dc2bSJeremy L Thompson 39*d416dc2bSJeremy L Thompson 40*d416dc2bSJeremy L Thompson 41*d416dc2bSJeremy L Thompson /** 42*d416dc2bSJeremy L Thompson * Operator CPU implementation using deal.II. 43*d416dc2bSJeremy L Thompson */ 44*d416dc2bSJeremy L Thompson template <int dim, typename Number> 45*d416dc2bSJeremy L Thompson class OperatorDealii : public OperatorBase<Number, MemorySpace::Host> 46*d416dc2bSJeremy L Thompson { 47*d416dc2bSJeremy L Thompson public: 48*d416dc2bSJeremy L Thompson using VectorType = typename OperatorBase<Number, MemorySpace::Host>::VectorType; 49*d416dc2bSJeremy L Thompson 50*d416dc2bSJeremy L Thompson /** 51*d416dc2bSJeremy L Thompson * Constructor. 52*d416dc2bSJeremy L Thompson */ OperatorDealii(const Mapping<dim> & mapping,const DoFHandler<dim> & dof_handler,const AffineConstraints<Number> & constraints,const Quadrature<dim> & quadrature,const BPType & bp)53*d416dc2bSJeremy L Thompson OperatorDealii(const Mapping<dim> &mapping, 54*d416dc2bSJeremy L Thompson const DoFHandler<dim> &dof_handler, 55*d416dc2bSJeremy L Thompson const AffineConstraints<Number> &constraints, 56*d416dc2bSJeremy L Thompson const Quadrature<dim> &quadrature, 57*d416dc2bSJeremy L Thompson const BPType &bp) 58*d416dc2bSJeremy L Thompson : mapping(mapping) 59*d416dc2bSJeremy L Thompson , dof_handler(dof_handler) 60*d416dc2bSJeremy L Thompson , constraints(constraints) 61*d416dc2bSJeremy L Thompson , quadrature(quadrature) 62*d416dc2bSJeremy L Thompson , bp(bp) 63*d416dc2bSJeremy L Thompson { 64*d416dc2bSJeremy L Thompson reinit(); 65*d416dc2bSJeremy L Thompson } 66*d416dc2bSJeremy L Thompson 67*d416dc2bSJeremy L Thompson /** 68*d416dc2bSJeremy L Thompson * Destructor. 69*d416dc2bSJeremy L Thompson */ 70*d416dc2bSJeremy L Thompson ~OperatorDealii() = default; 71*d416dc2bSJeremy L Thompson 72*d416dc2bSJeremy L Thompson /** 73*d416dc2bSJeremy L Thompson * Initialized internal data structures, particularly, MatrixFree. 74*d416dc2bSJeremy L Thompson */ 75*d416dc2bSJeremy L Thompson void reinit()76*d416dc2bSJeremy L Thompson reinit() override 77*d416dc2bSJeremy L Thompson { 78*d416dc2bSJeremy L Thompson // configure MatrixFree 79*d416dc2bSJeremy L Thompson typename MatrixFree<dim, Number>::AdditionalData additional_data; 80*d416dc2bSJeremy L Thompson additional_data.tasks_parallel_scheme = 81*d416dc2bSJeremy L Thompson MatrixFree<dim, Number>::AdditionalData::TasksParallelScheme::none; 82*d416dc2bSJeremy L Thompson 83*d416dc2bSJeremy L Thompson // create MatrixFree 84*d416dc2bSJeremy L Thompson matrix_free.reinit(mapping, dof_handler, constraints, quadrature, additional_data); 85*d416dc2bSJeremy L Thompson } 86*d416dc2bSJeremy L Thompson 87*d416dc2bSJeremy L Thompson /** 88*d416dc2bSJeremy L Thompson * Matrix-vector product. 89*d416dc2bSJeremy L Thompson */ 90*d416dc2bSJeremy L Thompson void vmult(VectorType & dst,const VectorType & src)91*d416dc2bSJeremy L Thompson vmult(VectorType &dst, const VectorType &src) const override 92*d416dc2bSJeremy L Thompson { 93*d416dc2bSJeremy L Thompson if (dof_handler.get_fe().n_components() == 1) 94*d416dc2bSJeremy L Thompson { 95*d416dc2bSJeremy L Thompson matrix_free.cell_loop(&OperatorDealii::do_cell_integral_range<1>, this, dst, src, true); 96*d416dc2bSJeremy L Thompson } 97*d416dc2bSJeremy L Thompson else 98*d416dc2bSJeremy L Thompson { 99*d416dc2bSJeremy L Thompson AssertThrow(dof_handler.get_fe().n_components() == dim, ExcInternalError()); 100*d416dc2bSJeremy L Thompson 101*d416dc2bSJeremy L Thompson matrix_free.cell_loop(&OperatorDealii::do_cell_integral_range<dim>, this, dst, src, true); 102*d416dc2bSJeremy L Thompson } 103*d416dc2bSJeremy L Thompson } 104*d416dc2bSJeremy L Thompson 105*d416dc2bSJeremy L Thompson /** 106*d416dc2bSJeremy L Thompson * Initialize vector. 107*d416dc2bSJeremy L Thompson */ 108*d416dc2bSJeremy L Thompson void initialize_dof_vector(VectorType & vec)109*d416dc2bSJeremy L Thompson initialize_dof_vector(VectorType &vec) const override 110*d416dc2bSJeremy L Thompson { 111*d416dc2bSJeremy L Thompson matrix_free.initialize_dof_vector(vec); 112*d416dc2bSJeremy L Thompson } 113*d416dc2bSJeremy L Thompson 114*d416dc2bSJeremy L Thompson /** 115*d416dc2bSJeremy L Thompson * Compute inverse of diagonal. 116*d416dc2bSJeremy L Thompson */ 117*d416dc2bSJeremy L Thompson void compute_inverse_diagonal(VectorType & diagonal)118*d416dc2bSJeremy L Thompson compute_inverse_diagonal(VectorType &diagonal) const override 119*d416dc2bSJeremy L Thompson { 120*d416dc2bSJeremy L Thompson this->initialize_dof_vector(diagonal); 121*d416dc2bSJeremy L Thompson 122*d416dc2bSJeremy L Thompson if (dof_handler.get_fe().n_components() == 1) 123*d416dc2bSJeremy L Thompson { 124*d416dc2bSJeremy L Thompson MatrixFreeTools::compute_diagonal(matrix_free, 125*d416dc2bSJeremy L Thompson diagonal, 126*d416dc2bSJeremy L Thompson &OperatorDealii::do_cell_integral_local<1>, 127*d416dc2bSJeremy L Thompson this); 128*d416dc2bSJeremy L Thompson } 129*d416dc2bSJeremy L Thompson else 130*d416dc2bSJeremy L Thompson { 131*d416dc2bSJeremy L Thompson AssertThrow(dof_handler.get_fe().n_components() == dim, ExcInternalError()); 132*d416dc2bSJeremy L Thompson 133*d416dc2bSJeremy L Thompson MatrixFreeTools::compute_diagonal(matrix_free, 134*d416dc2bSJeremy L Thompson diagonal, 135*d416dc2bSJeremy L Thompson &OperatorDealii::do_cell_integral_local<dim>, 136*d416dc2bSJeremy L Thompson this); 137*d416dc2bSJeremy L Thompson } 138*d416dc2bSJeremy L Thompson 139*d416dc2bSJeremy L Thompson for (auto &i : diagonal) 140*d416dc2bSJeremy L Thompson i = (std::abs(i) > 1.0e-10) ? (1.0 / i) : 1.0; 141*d416dc2bSJeremy L Thompson } 142*d416dc2bSJeremy L Thompson 143*d416dc2bSJeremy L Thompson private: 144*d416dc2bSJeremy L Thompson /** 145*d416dc2bSJeremy L Thompson * Cell integral without vector access. 146*d416dc2bSJeremy L Thompson */ 147*d416dc2bSJeremy L Thompson template <int n_components> 148*d416dc2bSJeremy L Thompson void do_cell_integral_local(FEEvaluation<dim,-1,0,n_components,Number> & phi)149*d416dc2bSJeremy L Thompson do_cell_integral_local(FEEvaluation<dim, -1, 0, n_components, Number> &phi) const 150*d416dc2bSJeremy L Thompson { 151*d416dc2bSJeremy L Thompson if (bp <= BPType::BP2) // mass matrix 152*d416dc2bSJeremy L Thompson { 153*d416dc2bSJeremy L Thompson phi.evaluate(EvaluationFlags::values); 154*d416dc2bSJeremy L Thompson for (const auto q : phi.quadrature_point_indices()) 155*d416dc2bSJeremy L Thompson phi.submit_value(phi.get_value(q), q); 156*d416dc2bSJeremy L Thompson phi.integrate(EvaluationFlags::values); 157*d416dc2bSJeremy L Thompson } 158*d416dc2bSJeremy L Thompson else // Poisson operator 159*d416dc2bSJeremy L Thompson { 160*d416dc2bSJeremy L Thompson phi.evaluate(EvaluationFlags::gradients); 161*d416dc2bSJeremy L Thompson for (const auto q : phi.quadrature_point_indices()) 162*d416dc2bSJeremy L Thompson phi.submit_gradient(phi.get_gradient(q), q); 163*d416dc2bSJeremy L Thompson phi.integrate(EvaluationFlags::gradients); 164*d416dc2bSJeremy L Thompson } 165*d416dc2bSJeremy L Thompson } 166*d416dc2bSJeremy L Thompson 167*d416dc2bSJeremy L Thompson /** 168*d416dc2bSJeremy L Thompson * Cell integral on a range of cells. 169*d416dc2bSJeremy L Thompson */ 170*d416dc2bSJeremy L Thompson template <int n_components> 171*d416dc2bSJeremy L Thompson void do_cell_integral_range(const MatrixFree<dim,Number> & matrix_free,VectorType & dst,const VectorType & src,const std::pair<unsigned int,unsigned int> & range)172*d416dc2bSJeremy L Thompson do_cell_integral_range(const MatrixFree<dim, Number> &matrix_free, 173*d416dc2bSJeremy L Thompson VectorType &dst, 174*d416dc2bSJeremy L Thompson const VectorType &src, 175*d416dc2bSJeremy L Thompson const std::pair<unsigned int, unsigned int> &range) const 176*d416dc2bSJeremy L Thompson { 177*d416dc2bSJeremy L Thompson FEEvaluation<dim, -1, 0, n_components, Number> phi(matrix_free, range); 178*d416dc2bSJeremy L Thompson 179*d416dc2bSJeremy L Thompson for (unsigned cell = range.first; cell < range.second; ++cell) 180*d416dc2bSJeremy L Thompson { 181*d416dc2bSJeremy L Thompson phi.reinit(cell); 182*d416dc2bSJeremy L Thompson phi.read_dof_values(src); // read source vector 183*d416dc2bSJeremy L Thompson do_cell_integral_local(phi); // cell integral 184*d416dc2bSJeremy L Thompson phi.distribute_local_to_global(dst); // write to destination vector 185*d416dc2bSJeremy L Thompson } 186*d416dc2bSJeremy L Thompson } 187*d416dc2bSJeremy L Thompson 188*d416dc2bSJeremy L Thompson /** 189*d416dc2bSJeremy L Thompson * Mapping object passed to the constructor. 190*d416dc2bSJeremy L Thompson */ 191*d416dc2bSJeremy L Thompson const Mapping<dim> &mapping; 192*d416dc2bSJeremy L Thompson 193*d416dc2bSJeremy L Thompson /** 194*d416dc2bSJeremy L Thompson * DoFHandler object passed to the constructor. 195*d416dc2bSJeremy L Thompson */ 196*d416dc2bSJeremy L Thompson const DoFHandler<dim> &dof_handler; 197*d416dc2bSJeremy L Thompson 198*d416dc2bSJeremy L Thompson /** 199*d416dc2bSJeremy L Thompson * Constraints object passed to the constructor. 200*d416dc2bSJeremy L Thompson */ 201*d416dc2bSJeremy L Thompson const AffineConstraints<Number> &constraints; 202*d416dc2bSJeremy L Thompson 203*d416dc2bSJeremy L Thompson /** 204*d416dc2bSJeremy L Thompson * Quadrature rule object passed to the constructor. 205*d416dc2bSJeremy L Thompson */ 206*d416dc2bSJeremy L Thompson const Quadrature<dim> &quadrature; 207*d416dc2bSJeremy L Thompson 208*d416dc2bSJeremy L Thompson /** 209*d416dc2bSJeremy L Thompson * Selected BP. 210*d416dc2bSJeremy L Thompson */ 211*d416dc2bSJeremy L Thompson const BPType bp; 212*d416dc2bSJeremy L Thompson 213*d416dc2bSJeremy L Thompson /** 214*d416dc2bSJeremy L Thompson * MatrixFree object. 215*d416dc2bSJeremy L Thompson */ 216*d416dc2bSJeremy L Thompson MatrixFree<dim, Number> matrix_free; 217*d416dc2bSJeremy L Thompson }; 218*d416dc2bSJeremy L Thompson 219*d416dc2bSJeremy L Thompson #endif 220