Line data Source code
1 : //* This file is part of the MOOSE framework 2 : //* https://mooseframework.inl.gov 3 : //* 4 : //* All rights reserved, see COPYRIGHT for full restrictions 5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT 6 : //* 7 : //* Licensed under LGPL 2.1, please see LICENSE for details 8 : //* https://www.gnu.org/licenses/lgpl-2.1.html 9 : 10 : #ifdef MFEM_ENABLED 11 : 12 : #include "ScaleIntegrator.h" 13 : 14 : namespace Moose::MFEM 15 : { 16 : 17 : void 18 1 : ScaleIntegrator::SetIntRule(const mfem::IntegrationRule * ir) 19 : { 20 1 : IntRule = ir; 21 1 : _integrator->SetIntRule(ir); 22 1 : } 23 : 24 : void 25 161697 : ScaleIntegrator::AssembleElementMatrix(const mfem::FiniteElement & el, 26 : mfem::ElementTransformation & Trans, 27 : mfem::DenseMatrix & elmat) 28 : { 29 161697 : CheckIntegrator(); 30 161697 : _integrator->AssembleElementMatrix(el, Trans, elmat); 31 161697 : elmat *= _scale; 32 161697 : } 33 : 34 : void 35 1 : ScaleIntegrator::AssembleElementMatrix2(const mfem::FiniteElement & el1, 36 : const mfem::FiniteElement & el2, 37 : mfem::ElementTransformation & Trans, 38 : mfem::DenseMatrix & elmat) 39 : { 40 1 : CheckIntegrator(); 41 1 : _integrator->AssembleElementMatrix2(el1, el2, Trans, elmat); 42 1 : elmat *= _scale; 43 1 : } 44 : 45 : void 46 11 : ScaleIntegrator::AssembleFaceMatrix(const mfem::FiniteElement & el1, 47 : const mfem::FiniteElement & el2, 48 : mfem::FaceElementTransformations & Trans, 49 : mfem::DenseMatrix & elmat) 50 : { 51 11 : CheckIntegrator(); 52 11 : _integrator->AssembleFaceMatrix(el1, el2, Trans, elmat); 53 11 : elmat *= _scale; 54 11 : } 55 : 56 : void 57 1 : ScaleIntegrator::AssemblePA(const mfem::FiniteElementSpace & fes) 58 : { 59 1 : CheckIntegrator(); 60 1 : _integrator->AssemblePA(fes); 61 1 : } 62 : 63 : void 64 17 : ScaleIntegrator::AssembleDiagonalPA(mfem::Vector & diag) 65 : { 66 17 : _integrator->AssembleDiagonalPA(diag); 67 17 : diag *= _scale; 68 17 : } 69 : 70 : void 71 1 : ScaleIntegrator::AssemblePAInteriorFaces(const mfem::FiniteElementSpace & fes) 72 : { 73 1 : _integrator->AssemblePAInteriorFaces(fes); 74 1 : } 75 : 76 : void 77 1 : ScaleIntegrator::AssemblePABoundaryFaces(const mfem::FiniteElementSpace & fes) 78 : { 79 1 : _integrator->AssemblePABoundaryFaces(fes); 80 1 : } 81 : 82 : void 83 1 : ScaleIntegrator::AddMultPA(const mfem::Vector & x, mfem::Vector & y) const 84 : { 85 : // y += Mx*scale 86 1 : mfem::Vector Mx(y.Size()); 87 1 : Mx = 0.0; 88 1 : _integrator->AddMultPA(x, Mx); 89 1 : Mx *= _scale; 90 1 : y += Mx; 91 1 : } 92 : 93 : void 94 1 : ScaleIntegrator::AddMultTransposePA(const mfem::Vector & x, mfem::Vector & y) const 95 : { 96 : // y += M^T x*scale 97 1 : mfem::Vector MTx(y.Size()); 98 1 : MTx = 0.0; 99 1 : _integrator->AddMultTransposePA(x, MTx); 100 1 : MTx *= _scale; 101 1 : y += MTx; 102 1 : } 103 : 104 : void 105 20 : ScaleIntegrator::AssembleEA(const mfem::FiniteElementSpace & fes, 106 : mfem::Vector & emat, 107 : const bool add) 108 : { 109 20 : CheckIntegrator(); 110 20 : if (add) 111 : { 112 20 : mfem::Vector emat_scale(emat.Size()); 113 20 : _integrator->AssembleEA(fes, emat_scale, false); 114 20 : emat_scale *= _scale; 115 20 : emat += emat_scale; 116 20 : } 117 : else 118 : { 119 0 : _integrator->AssembleEA(fes, emat, add); 120 0 : emat *= _scale; 121 : } 122 20 : } 123 : 124 219 : ScaleIntegrator::~ScaleIntegrator() 125 : { 126 111 : if (_own_integrator) 127 94 : delete _integrator; 128 219 : } 129 : 130 : } // namespace Moose::MFEM 131 : 132 : #endif