www.mooseframework.org
Public Member Functions | Public Attributes | Protected Attributes | List of all members
XFEMInterface Class Referenceabstract

This is the XFEMInterface class. More...

#include <XFEMInterface.h>

Inheritance diagram for XFEMInterface:
[legend]

Public Member Functions

 XFEMInterface (const InputParameters &params)
 Constructor. More...
 
virtual ~XFEMInterface ()
 Destructor. More...
 
void setMesh (MooseMesh *mesh)
 Set the pointer to the master mesh that is modified by XFEM. More...
 
void setDisplacedMesh (MooseMesh *displaced_mesh)
 Set the pointer to the displaced mesh that is modified by XFEM. More...
 
void setMaterialData (const std::vector< MaterialData *> &data)
 Set the pointer to the MaterialData. More...
 
void setBoundaryMaterialData (const std::vector< MaterialData *> &data)
 Set the pointer to the Boundary MaterialData. More...
 
virtual bool update (Real time, const std::vector< std::shared_ptr< NonlinearSystemBase >> &nl, AuxiliarySystem &aux)=0
 Method to update the mesh due to modified cut definitions. More...
 
virtual void initSolution (const std::vector< std::shared_ptr< NonlinearSystemBase >> &nl, AuxiliarySystem &aux)=0
 Initialize the solution on newly created nodes. More...
 
virtual bool getXFEMWeights (MooseArray< Real > &weights, const Elem *elem, QBase *qrule, const MooseArray< Point > &q_points)=0
 Get the factors for the QP weighs for XFEM partial elements. More...
 
virtual bool getXFEMFaceWeights (MooseArray< Real > &weights, const Elem *elem, QBase *qrule, const MooseArray< Point > &q_points, unsigned int side)=0
 Get the factors for the face QP weighs for XFEM partial elements. More...
 
virtual bool updateHeal ()=0
 Potentially update the mesh by healing previous XFEM cuts. More...
 

Public Attributes

const ConsoleStream _console
 An instance of helper class to write streams to the Console objects. More...
 

Protected Attributes

FEProblemBase_fe_problem
 
std::vector< MaterialData * > _material_data
 
std::vector< MaterialData * > _bnd_material_data
 
MooseMesh_moose_mesh
 
MooseMesh_moose_displaced_mesh
 
MeshBase * _mesh
 
MeshBase * _displaced_mesh
 

Detailed Description

This is the XFEMInterface class.

This is an abstract base class that defines interfaces with a class that dynamically modifies the mesh in support of a phantom node approach for XFEM

Definition at line 37 of file XFEMInterface.h.

Constructor & Destructor Documentation

◆ XFEMInterface()

XFEMInterface::XFEMInterface ( const InputParameters params)
inlineexplicit

Constructor.

Definition at line 43 of file XFEMInterface.h.

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  }
T getCheckedPointerParam(const std::string &name, const std::string &error_string="") const
Verifies that the requested parameter exists and is not NULL and returns it to the caller...
Base class for MOOSE-based applications.
Definition: MooseApp.h:73
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
FEProblemBase * _fe_problem
MooseMesh * _moose_displaced_mesh
ConsoleStreamInterface(MooseApp &app)
A class for providing a helper stream object for writting message to all the Output objects...
MeshBase * _mesh
MooseMesh * _moose_mesh
MeshBase * _displaced_mesh

◆ ~XFEMInterface()

virtual XFEMInterface::~XFEMInterface ( )
inlinevirtual

Destructor.

Definition at line 56 of file XFEMInterface.h.

56 {}

Member Function Documentation

◆ getXFEMFaceWeights()

virtual bool XFEMInterface::getXFEMFaceWeights ( MooseArray< Real > &  weights,
const Elem *  elem,
QBase *  qrule,
const MooseArray< Point > &  q_points,
unsigned int  side 
)
pure virtual

Get the factors for the face QP weighs for XFEM partial elements.

Parameters
weightsThe new weights at element face quadrature points
elemThe element for which the weights are adjusted
qruleThe quadrature rule for the face integration
q_pointsThe vector of quadrature points at element face
sideThe side of element for which the weights are adjusted

◆ getXFEMWeights()

virtual bool XFEMInterface::getXFEMWeights ( MooseArray< Real > &  weights,
const Elem *  elem,
QBase *  qrule,
const MooseArray< Point > &  q_points 
)
pure virtual

Get the factors for the QP weighs for XFEM partial elements.

Parameters
weightsThe new weights at element quadrature points
elemThe element for which the weights are adjusted
qruleThe quadrature rule for the volume integration
q_pointsThe vector of quadrature points

◆ initSolution()

virtual void XFEMInterface::initSolution ( const std::vector< std::shared_ptr< NonlinearSystemBase >> &  nl,
AuxiliarySystem aux 
)
pure virtual

Initialize the solution on newly created nodes.

◆ setBoundaryMaterialData()

void XFEMInterface::setBoundaryMaterialData ( const std::vector< MaterialData *> &  data)
inline

Set the pointer to the Boundary MaterialData.

Definition at line 84 of file XFEMInterface.h.

85  {
86  _bnd_material_data = data;
87  }
std::vector< MaterialData * > _bnd_material_data

◆ setDisplacedMesh()

void XFEMInterface::setDisplacedMesh ( MooseMesh displaced_mesh)
inline

Set the pointer to the displaced mesh that is modified by XFEM.

Definition at line 70 of file XFEMInterface.h.

71  {
72  _moose_displaced_mesh = displaced_mesh;
73  _displaced_mesh = &displaced_mesh->getMesh();
74  }
MooseMesh * _moose_displaced_mesh
MeshBase & getMesh()
Accessor for the underlying libMesh Mesh object.
Definition: MooseMesh.C:3198
MeshBase * _displaced_mesh

◆ setMaterialData()

void XFEMInterface::setMaterialData ( const std::vector< MaterialData *> &  data)
inline

Set the pointer to the MaterialData.

Definition at line 79 of file XFEMInterface.h.

79 { _material_data = data; }
std::vector< MaterialData * > _material_data

◆ setMesh()

void XFEMInterface::setMesh ( MooseMesh mesh)
inline

Set the pointer to the master mesh that is modified by XFEM.

Definition at line 61 of file XFEMInterface.h.

62  {
63  _moose_mesh = mesh;
64  _mesh = &mesh->getMesh();
65  }
MeshBase & mesh
MeshBase * _mesh
MooseMesh * _moose_mesh

◆ update()

virtual bool XFEMInterface::update ( Real  time,
const std::vector< std::shared_ptr< NonlinearSystemBase >> &  nl,
AuxiliarySystem aux 
)
pure virtual

Method to update the mesh due to modified cut definitions.

◆ updateHeal()

virtual bool XFEMInterface::updateHeal ( )
pure virtual

Potentially update the mesh by healing previous XFEM cuts.

Returns
true if the mesh has been updated due to healing

Member Data Documentation

◆ _bnd_material_data

std::vector<MaterialData *> XFEMInterface::_bnd_material_data
protected

Definition at line 138 of file XFEMInterface.h.

Referenced by setBoundaryMaterialData().

◆ _console

const ConsoleStream ConsoleStreamInterface::_console
inherited

An instance of helper class to write streams to the Console objects.

Definition at line 31 of file ConsoleStreamInterface.h.

Referenced by IterationAdaptiveDT::acceptStep(), MeshOnlyAction::act(), SetupDebugAction::act(), MaterialOutputAction::act(), Adaptivity::adaptMesh(), FEProblemBase::adaptMesh(), PerfGraph::addToExecutionList(), SimplePredictor::apply(), SystemBase::applyScalingFactors(), MultiApp::backup(), FEProblemBase::backupMultiApps(), CoarsenedPiecewiseLinear::buildCoarsenedGrid(), MeshDiagnosticsGenerator::checkElementOverlap(), MeshDiagnosticsGenerator::checkElementTypes(), MeshDiagnosticsGenerator::checkElementVolumes(), FEProblemBase::checkExceptionAndStopSolve(), MeshDiagnosticsGenerator::checkLocalJacobians(), MeshDiagnosticsGenerator::checkNonConformalMesh(), MeshDiagnosticsGenerator::checkNonConformalMeshFromAdaptivity(), MeshDiagnosticsGenerator::checkNonPlanarSides(), FEProblemBase::checkProblemIntegrity(), ReferenceResidualProblem::checkRelativeConvergence(), MeshDiagnosticsGenerator::checkSidesetsOrientation(), IterationAdaptiveDT::computeAdaptiveDT(), Transient::computeConstrainedDT(), FixedPointSolve::computeCustomConvergencePostprocessor(), NonlinearSystemBase::computeDamping(), IterationAdaptiveDT::computeDT(), IterationAdaptiveDT::computeFailedDT(), IterationAdaptiveDT::computeInitialDT(), IterationAdaptiveDT::computeInterpolationDT(), NonlinearSystemBase::computeScaling(), Problem::console(), IterationAdaptiveDT::constrainStep(), TimeStepper::constrainStep(), MultiApp::createApp(), FEProblemBase::execMultiApps(), FEProblemBase::execMultiAppTransfers(), MessageFromInput::execute(), Steady::execute(), Eigenvalue::execute(), ActionWarehouse::executeActionsWithAction(), ActionWarehouse::executeAllActions(), ElementQualityChecker::finalize(), FEProblemBase::finishMultiAppStep(), MeshRepairGenerator::fixOverlappingNodes(), CoarsenBlockGenerator::generate(), MeshGenerator::generateInternal(), VariableCondensationPreconditioner::getDofToCondense(), InversePowerMethod::init(), NonlinearEigen::init(), FEProblemBase::initialAdaptMesh(), EigenExecutionerBase::inversePowerIteration(), FEProblemBase::joinAndFinalize(), Transient::keepGoing(), IterationAdaptiveDT::limitDTByFunction(), IterationAdaptiveDT::limitDTToPostprocessorValue(), FEProblemBase::logAdd(), EigenExecutionerBase::makeBXConsistent(), Console::meshChanged(), MooseBaseErrorInterface::mooseDeprecated(), MooseBaseErrorInterface::mooseInfo(), MooseBaseErrorInterface::mooseWarning(), MooseBaseErrorInterface::mooseWarningNonPrefixed(), ReferenceResidualProblem::nonlinearConvergenceSetup(), ReporterDebugOutput::output(), PerfGraphOutput::output(), MaterialPropertyDebugOutput::output(), DOFMapOutput::output(), VariableResidualNormsDebugOutput::output(), Console::output(), ControlOutput::outputActiveObjects(), ControlOutput::outputChangedControls(), ControlOutput::outputControls(), Console::outputInput(), Console::outputPostprocessors(), PseudoTimestep::outputPseudoTimestep(), Console::outputReporters(), Console::outputScalarVariables(), Console::outputSystemInformation(), FEProblemBase::possiblyRebuildGeomSearchPatches(), EigenExecutionerBase::postExecute(), AB2PredictorCorrector::postSolve(), ActionWarehouse::printActionDependencySets(), SolutionInvalidity::printDebug(), EigenExecutionerBase::printEigenvalue(), SecantSolve::printFixedPointConvergenceHistory(), SteffensenSolve::printFixedPointConvergenceHistory(), PicardSolve::printFixedPointConvergenceHistory(), FixedPointSolve::printFixedPointConvergenceReason(), PerfGraphLivePrint::printLiveMessage(), MaterialPropertyDebugOutput::printMaterialMap(), PerfGraphLivePrint::printStats(), AutomaticMortarGeneration::projectPrimaryNodesSinglePair(), AutomaticMortarGeneration::projectSecondaryNodesSinglePair(), CoarsenBlockGenerator::recursiveCoarsen(), SolutionTimeAdaptiveDT::rejectStep(), MultiApp::restore(), FEProblemBase::restoreMultiApps(), SimplePredictor::shouldApply(), Checkpoint::shouldOutput(), SubProblem::showFunctorRequestors(), SubProblem::showFunctors(), FullSolveMultiApp::showStatusMessage(), FEProblemSolve::solve(), FixedPointSolve::solve(), NonlinearSystem::solve(), EigenProblem::solve(), LStableDirk2::solve(), LStableDirk3::solve(), ImplicitMidpoint::solve(), ExplicitTVDRK2::solve(), LStableDirk4::solve(), AStableDirk4::solve(), ExplicitRK2::solve(), TransientMultiApp::solveStep(), FixedPointSolve::solveStep(), PerfGraphLivePrint::start(), AB2PredictorCorrector::step(), NonlinearEigen::takeStep(), Transient::takeStep(), Console::writeTimestepInformation(), Console::writeVariableNorms(), and FEProblemBase::~FEProblemBase().

◆ _displaced_mesh

MeshBase* XFEMInterface::_displaced_mesh
protected

Definition at line 143 of file XFEMInterface.h.

Referenced by setDisplacedMesh().

◆ _fe_problem

FEProblemBase* XFEMInterface::_fe_problem
protected

Definition at line 136 of file XFEMInterface.h.

◆ _material_data

std::vector<MaterialData *> XFEMInterface::_material_data
protected

Definition at line 137 of file XFEMInterface.h.

Referenced by setMaterialData().

◆ _mesh

MeshBase* XFEMInterface::_mesh
protected

Definition at line 142 of file XFEMInterface.h.

Referenced by setMesh().

◆ _moose_displaced_mesh

MooseMesh* XFEMInterface::_moose_displaced_mesh
protected

Definition at line 141 of file XFEMInterface.h.

Referenced by setDisplacedMesh().

◆ _moose_mesh

MooseMesh* XFEMInterface::_moose_mesh
protected

Definition at line 140 of file XFEMInterface.h.

Referenced by setMesh().


The documentation for this class was generated from the following file: