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 : #pragma once 11 : 12 : #include "ConsoleStreamInterface.h" 13 : #include "MooseTypes.h" 14 : #include "InputParameters.h" 15 : #include "MooseMesh.h" 16 : 17 : class MooseApp; 18 : class AuxiliarySystem; 19 : class NonlinearSystemBase; 20 : class MaterialData; 21 : class FEProblemBase; 22 : 23 : namespace libMesh 24 : { 25 : class MeshBase; 26 : class QBase; 27 : } // namespace libMesh 28 : 29 : /** 30 : * This is the XFEMInterface class. This is an abstract base 31 : * class that defines interfaces with a class that dynamically 32 : * modifies the mesh in support of a phantom node approach for XFEM 33 : */ 34 : 35 : // ------------------------------------------------------------ 36 : // XFEMInterface class definition 37 : class XFEMInterface : public ConsoleStreamInterface 38 : { 39 : public: 40 : /** 41 : * Constructor 42 : */ 43 : explicit XFEMInterface(const InputParameters & params) 44 : : ConsoleStreamInterface(*params.getCheckedPointerParam<MooseApp *>("_moose_app")), 45 : _fe_problem(params.getCheckedPointerParam<FEProblemBase *>("_fe_problem_base")), 46 : _moose_mesh(nullptr), 47 : _moose_displaced_mesh(nullptr), 48 : _mesh(nullptr), 49 : _displaced_mesh(nullptr) 50 : { 51 : } 52 : 53 : /** 54 : * Destructor 55 : */ 56 : virtual ~XFEMInterface() {} 57 : 58 : /** 59 : * Set the pointer to the master mesh that is modified by XFEM 60 : */ 61 0 : void setMesh(MooseMesh * mesh) 62 : { 63 0 : _moose_mesh = mesh; 64 0 : _mesh = &mesh->getMesh(); 65 0 : } 66 : 67 : /** 68 : * Set the pointer to the displaced mesh that is modified by XFEM 69 : */ 70 0 : void setDisplacedMesh(MooseMesh * displaced_mesh) 71 : { 72 0 : _moose_displaced_mesh = displaced_mesh; 73 0 : _displaced_mesh = &displaced_mesh->getMesh(); 74 0 : } 75 : 76 : /** 77 : * Set the pointer to the MaterialData 78 : */ 79 0 : void setMaterialData(const std::vector<MaterialData *> & data) { _material_data = data; } 80 : 81 : /** 82 : * Set the pointer to the Boundary MaterialData 83 : */ 84 0 : void setBoundaryMaterialData(const std::vector<MaterialData *> & data) 85 : { 86 0 : _bnd_material_data = data; 87 0 : } 88 : 89 : /** 90 : * Method to update the mesh due to modified cut definitions 91 : */ 92 : virtual bool update(Real time, 93 : const std::vector<std::shared_ptr<NonlinearSystemBase>> & nl, 94 : AuxiliarySystem & aux) = 0; 95 : 96 : /** 97 : * Initialize the solution on newly created nodes 98 : */ 99 : virtual void initSolution(const std::vector<std::shared_ptr<NonlinearSystemBase>> & nl, 100 : AuxiliarySystem & aux) = 0; 101 : 102 : /** 103 : * Get the factors for the QP weighs for XFEM partial elements 104 : * @param weights The new weights at element quadrature points 105 : * @param elem The element for which the weights are adjusted 106 : * @param qrule The quadrature rule for the volume integration 107 : * @param q_points The vector of quadrature points 108 : */ 109 : 110 : virtual bool getXFEMWeights(MooseArray<Real> & weights, 111 : const Elem * elem, 112 : libMesh::QBase * qrule, 113 : const MooseArray<Point> & q_points) = 0; 114 : 115 : /** 116 : * Get the factors for the face QP weighs for XFEM partial elements 117 : * @param weights The new weights at element face quadrature points 118 : * @param elem The element for which the weights are adjusted 119 : * @param qrule The quadrature rule for the face integration 120 : * @param q_points The vector of quadrature points at element face 121 : * @param side The side of element for which the weights are adjusted 122 : */ 123 : virtual bool getXFEMFaceWeights(MooseArray<Real> & weights, 124 : const Elem * elem, 125 : libMesh::QBase * qrule, 126 : const MooseArray<Point> & q_points, 127 : unsigned int side) = 0; 128 : 129 : /** 130 : * Potentially update the mesh by healing previous XFEM cuts. 131 : * @return true if the mesh has been updated due to healing 132 : **/ 133 : virtual bool updateHeal() = 0; 134 : 135 : protected: 136 : FEProblemBase * _fe_problem; 137 : std::vector<MaterialData *> _material_data; 138 : std::vector<MaterialData *> _bnd_material_data; 139 : 140 : MooseMesh * _moose_mesh; 141 : MooseMesh * _moose_displaced_mesh; 142 : MeshBase * _mesh; 143 : MeshBase * _displaced_mesh; 144 : };