OGS
ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >

Definition at line 15 of file HMPhaseFieldProcess.h.

#include <HMPhaseFieldProcess.h>

Inheritance diagram for ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >:
[legend]

Public Member Functions

 HMPhaseFieldProcess (std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, HMPhaseFieldProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
ODESystem interface
bool isLinear () const override
Public Member Functions inherited from ProcessLib::Process
 Process (std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
void preTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id)
 Preprocessing before starting assembly for new timestep.
void postTimestep (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep.
void postNonLinearSolver (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id)
void preIteration (const unsigned iter, GlobalVector const &x) final
void computeSecondaryVariable (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
 compute secondary variables for the coupled equations or for output.
NumLib::IterationResult postIteration (GlobalVector const &x) final
void initialize (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
void setInitialConditions (std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
void updateDeactivatedSubdomains (double const time, const int process_id)
virtual bool isMonolithicSchemeUsed () const
virtual void extrapolateIntegrationPointValuesToNodes (const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &)
void preAssemble (const double t, double const dt, GlobalVector const &x) final
void assemble (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) final
void assembleWithJacobian (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) final
void preOutput (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * getKnownSolutions (double const t, GlobalVector const &x, int const process_id) const final
MeshLib::MeshgetMesh () const
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables () const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables (const int process_id) const
std::vector< std::size_t > const & getActiveElementIDs () const
SecondaryVariableCollection const & getSecondaryVariables () const
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const & getIntegrationPointWriters () const
virtual Eigen::Vector3d getFlux (std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const
virtual void solveReactionEquation (std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
bool requiresNormalization () const override
Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
virtual std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes)
virtual ~SubmeshAssemblySupport ()=default

Private Types

using LocalAssemblerInterface = HMPhaseFieldLocalAssemblerInterface

Private Member Functions

void constructDofTable () override
void initializeBoundaryConditions (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
void initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
 Process specific initialization called by initialize().
void assembleConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
void assembleWithJacobianConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) override
void preTimestepConcreteProcess (std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
void postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &, const double t, const double dt, int const process_id) override
void postNonLinearSolverConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id) override
void updateConstraints (GlobalVector &lower, GlobalVector &upper, int const process_id) override
NumLib::LocalToGlobalIndexMap const & getDOFTable (const int process_id) const override

Private Attributes

HMPhaseFieldProcessData< DisplacementDim > _process_data
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map_single_component
MeshLib::PropertyVector< double > * _nodal_forces = nullptr
GlobalSparsityPattern _sparsity_pattern_with_single_component
std::unique_ptr< GlobalVector_x_previous_timestep

Additional Inherited Members

Public Attributes inherited from ProcessLib::Process
std::string const name
Static Public Attributes inherited from ProcessLib::Process
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name = "constant_one"
Protected Member Functions inherited from ProcessLib::Process
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables (int const number_of_processes) const
NumLib::ExtrapolatorgetExtrapolator () const
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
void constructMonolithicProcessDofTable ()
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_process_id)
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
void setReleaseNodalForces (GlobalVector const *r_neq, int const process_id) override
Protected Attributes inherited from ProcessLib::Process
MeshLib::Mesh_mesh
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map
SecondaryVariableCollection _secondary_variables
CellAverageData cell_average_data_
std::unique_ptr< ProcessLib::AbstractJacobianAssembler_jacobian_assembler
VectorMatrixAssembler _global_assembler
const bool _use_monolithic_scheme
unsigned const _integration_order
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
GlobalSparsityPattern _sparsity_pattern
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
std::vector< BoundaryConditionCollection_boundary_conditions

Member Typedef Documentation

◆ LocalAssemblerInterface

template<int DisplacementDim>
using ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::LocalAssemblerInterface = HMPhaseFieldLocalAssemblerInterface
private

Definition at line 41 of file HMPhaseFieldProcess.h.

Constructor & Destructor Documentation

◆ HMPhaseFieldProcess()

template<int DisplacementDim>
ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::HMPhaseFieldProcess ( std::string name,
MeshLib::Mesh & mesh,
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && jacobian_assembler,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
unsigned const integration_order,
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > && process_variables,
HMPhaseFieldProcessData< DisplacementDim > && process_data,
SecondaryVariableCollection && secondary_variables,
bool const use_monolithic_scheme )

Definition at line 19 of file HMPhaseFieldProcess.cpp.

34{
35 // For numerical Jacobian assembler
36 if (this->_jacobian_assembler->isPerturbationEnabled())
37 {
39 "The numericial Jacobian assembler is not supported for the "
40 "HMPhaseField process.");
41 }
42
44 {
46 "Monolithic scheme is not implemented for the HMPhaseField "
47 "process.");
48 }
49
52}
#define OGS_FATAL(...)
Definition Error.h:19
HMPhaseFieldProcessData< DisplacementDim > _process_data
MeshLib::PropertyVector< double > * _nodal_forces
std::string const name
Definition Process.h:361
Process(std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
Definition Process.cpp:37
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:375
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)

References ProcessLib::Process::Process(), ProcessLib::Process::_jacobian_assembler, _nodal_forces, _process_data, MeshLib::getOrCreateMeshProperty(), ProcessLib::Process::name, MeshLib::Node, and OGS_FATAL.

Member Function Documentation

◆ assembleConcreteProcess()

template<int DisplacementDim>
void ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::assembleConcreteProcess ( const double t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 182 of file HMPhaseFieldProcess.cpp.

187{
188 OGS_FATAL(
189 "HMPhaseFieldLocalAssembler: assembly for the Picard non-linear solver "
190 "is not implemented.");
191}

References OGS_FATAL.

◆ assembleWithJacobianConcreteProcess()

template<int DisplacementDim>
void ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::assembleWithJacobianConcreteProcess ( const double t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
GlobalVector & b,
GlobalMatrix & Jac )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 194 of file HMPhaseFieldProcess.cpp.

198{
200
201 // For the staggered scheme
202 if (process_id == _process_data._phasefield_process_id)
203 {
204 DBUG(
205 "Assemble the Jacobian equations of phase field in "
206 "HMPhaseFieldProcess for the staggered scheme.");
207 }
208 else if (process_id == _process_data._hydro_process_id)
209 {
210 DBUG(
211 "Assemble the Jacobian equations of pressure in "
212 "HMPhaseFieldProcess for the staggered scheme.");
213 }
214 else
215 {
216 DBUG(
217 "Assemble the Jacobian equations of deformation in "
218 "HMPhaseFieldProcess for the staggered scheme.");
219 }
220
223 dof_tables.emplace_back(_local_to_global_index_map.get());
224
225 // Call global assembler for each local assembly item.
229 process_id, &b, &Jac);
230
231 if (process_id == _process_data._mechanics_related_process_id)
232 {
233 b.copyValues(*_nodal_forces);
235 _nodal_forces->begin(), [](double val) { return -val; });
236 }
237}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
VectorMatrixAssembler _global_assembler
Definition Process.h:376
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:367
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, _local_to_global_index_map_single_component, _nodal_forces, _process_data, ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), MathLib::EigenVector::copyValues(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::Process::getActiveElementIDs().

◆ constructDofTable()

template<int DisplacementDim>
void ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::constructDofTable ( )
overrideprivatevirtual

This function is for general cases, in which all equations of the coupled processes have the same number of unknowns. For the general cases with the staggered scheme, all equations of the coupled processes share one DOF table hold by _local_to_global_index_map. Other cases can be considered by overloading this member function in the derived class.

Reimplemented from ProcessLib::Process.

Definition at line 94 of file HMPhaseFieldProcess.cpp.

95{
96 // For displacement equation.
98 _process_data._mechanics_related_process_id);
99
100 // TODO move the two data members somewhere else.
101 // for extrapolation of secondary variables of stress or strain
107 // by location order is needed for output
109
111
112 // For the H (hydro) or PF (phasefield) process in the staggered scheme.
115}
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_process_id)
Definition Process.cpp:345
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:365
MeshLib::Mesh & _mesh
Definition Process.h:364
GlobalSparsityPattern computeSparsityPattern(LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
Computes a sparsity pattern for the given inputs.

References _local_to_global_index_map_single_component, ProcessLib::Process::_mesh, ProcessLib::Process::_mesh_subset_all_nodes, _process_data, _sparsity_pattern_with_single_component, NumLib::BY_LOCATION, NumLib::computeSparsityPattern(), and ProcessLib::Process::constructDofTableOfSpecifiedProcessStaggeredScheme().

◆ getDOFTable()

template<int DisplacementDim>
NumLib::LocalToGlobalIndexMap const & ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::getDOFTable ( const int process_id) const
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 81 of file HMPhaseFieldProcess.cpp.

82{
83 // For the M process (deformation) in the staggered scheme.
84 if (process_id == _process_data._mechanics_related_process_id)
85 {
87 }
88
89 // For the H (hydro) or PF (phasefield) process in the staggered scheme.
91}

References ProcessLib::Process::_local_to_global_index_map, _local_to_global_index_map_single_component, and _process_data.

Referenced by preTimestepConcreteProcess().

◆ getMatrixSpecifications()

template<int DisplacementDim>
MathLib::MatrixSpecifications ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::getMatrixSpecifications ( const int process_id) const
override

Definition at line 62 of file HMPhaseFieldProcess.cpp.

64{
65 // For the M process (deformation) in the staggered scheme.
66 if (process_id == _process_data._mechanics_related_process_id)
67 {
68 auto const& l = *_local_to_global_index_map;
69 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
70 &l.getGhostIndices(), &this->_sparsity_pattern};
71 }
72
73 // For the H (hydro) or PF (phasefield) process in the staggered scheme.
75 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
76 &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
77}
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:391

References ProcessLib::Process::_local_to_global_index_map, _local_to_global_index_map_single_component, _process_data, ProcessLib::Process::_sparsity_pattern, and _sparsity_pattern_with_single_component.

◆ initializeBoundaryConditions()

template<int DisplacementDim>
void ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::initializeBoundaryConditions ( std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media)
overrideprivatevirtual

Member function to initialize the boundary conditions for all coupled processes. It is called by initialize().

Reimplemented from ProcessLib::Process.

Definition at line 163 of file HMPhaseFieldProcess.cpp.

165{
166 // Staggered scheme:
167 // for the phase field
170 _process_data._phasefield_process_id, media);
171 // for the pressure
174 _process_data._hydro_process_id, media);
175 // for the deformation.
178 _process_data._mechanics_related_process_id, media);
179}
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:83

References ProcessLib::Process::_local_to_global_index_map, _local_to_global_index_map_single_component, _process_data, and ProcessLib::Process::initializeProcessBoundaryConditionsAndSourceTerms().

◆ initializeConcreteProcess()

template<int DisplacementDim>
void ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::initializeConcreteProcess ( NumLib::LocalToGlobalIndexMap const & dof_table,
MeshLib::Mesh const & mesh,
unsigned const integration_order )
overrideprivatevirtual

Process specific initialization called by initialize().

Implements ProcessLib::Process.

Definition at line 118 of file HMPhaseFieldProcess.cpp.

122{
125 mesh.getElements(), dof_table, _local_assemblers,
126 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
128
129 _secondary_variables.addSecondaryVariable(
130 "sigma",
135
136 _secondary_variables.addSecondaryVariable(
137 "epsilon",
142
143 _secondary_variables.addSecondaryVariable(
144 "width_node",
147
149 const_cast<MeshLib::Mesh&>(mesh), "damage", MeshLib::MeshItemType::Cell,
150 1);
151
153 const_cast<MeshLib::Mesh&>(mesh), "width", MeshLib::MeshItemType::Cell,
154 1);
155
156 // Initialize local assemblers after all variables have been set.
160}
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
SecondaryVariableCollection _secondary_variables
Definition Process.h:369
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:201
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
virtual std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtWidth(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0

References _local_assemblers, ProcessLib::Process::_local_to_global_index_map, _process_data, ProcessLib::Process::_secondary_variables, MeshLib::Cell, ProcessLib::SmallDeformation::createLocalAssemblers(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::HMPhaseField::HMPhaseFieldLocalAssemblerInterface::getIntPtEpsilon(), ProcessLib::HMPhaseField::HMPhaseFieldLocalAssemblerInterface::getIntPtSigma(), ProcessLib::HMPhaseField::HMPhaseFieldLocalAssemblerInterface::getIntPtWidth(), MeshLib::getOrCreateMeshProperty(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::Mesh::isAxiallySymmetric(), and ProcessLib::makeExtrapolator().

◆ isLinear()

template<int DisplacementDim>
bool ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::isLinear ( ) const
override

Definition at line 55 of file HMPhaseFieldProcess.cpp.

56{
57 return false;
58}

◆ postNonLinearSolverConcreteProcess()

template<int DisplacementDim>
void ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::postNonLinearSolverConcreteProcess ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
const double t,
double const dt,
int const process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 299 of file HMPhaseFieldProcess.cpp.

303{
305
308 dof_tables.emplace_back(_local_to_global_index_map.get());
309
310 if (process_id == _process_data._phasefield_process_id)
311 {
312 INFO("Update fracture width and porous properties");
316 }
317
318 if (process_id == _process_data._hydro_process_id)
319 {
320 INFO("PostNonLinearSolver for Hydro Process");
321
325 process_id);
326 }
327}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void postNonLinearSolver(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:376
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
virtual void approximateFractureWidth(std::size_t mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)=0

References _local_assemblers, ProcessLib::Process::_local_to_global_index_map, _local_to_global_index_map_single_component, _process_data, ProcessLib::HMPhaseField::HMPhaseFieldLocalAssemblerInterface::approximateFractureWidth(), NumLib::SerialExecutor::executeMemberOnDereferenced(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), ProcessLib::Process::getDOFTables(), INFO(), and ProcessLib::LocalAssemblerInterface::postNonLinearSolver().

◆ postTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::postTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
const double t,
const double dt,
int const process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 260 of file HMPhaseFieldProcess.cpp.

264{
265 if (process_id == _process_data._phasefield_process_id)
266 {
267 DBUG("PostTimestep HMPhaseFieldProcess.");
268
269 _process_data.elastic_energy = 0.0;
270 _process_data.surface_energy = 0.0;
271 _process_data.pressure_work = 0.0;
272
274
275 dof_tables.emplace_back(
277 dof_tables.emplace_back(
279 dof_tables.emplace_back(_local_to_global_index_map.get());
280
284 _process_data.elastic_energy, _process_data.surface_energy,
285 _process_data.pressure_work);
286
287 showEnergyAndWork(t, _process_data.elastic_energy,
288 _process_data.surface_energy,
289 _process_data.pressure_work);
290 }
291
295 process_id);
296}
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
void showEnergyAndWork(const double t, double &_elastic_energy, double &_surface_energy, double &_pressure_work)
virtual void computeEnergy(std::size_t mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &elastic_energy, double &surface_energy, double &pressure_work)=0

References _local_assemblers, ProcessLib::Process::_local_to_global_index_map, _local_to_global_index_map_single_component, _process_data, ProcessLib::HMPhaseField::HMPhaseFieldLocalAssemblerInterface::computeEnergy(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), ProcessLib::Process::getDOFTables(), ProcessLib::LocalAssemblerInterface::postTimestep(), and ProcessLib::HMPhaseField::showEnergyAndWork().

◆ preTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::preTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
double const t,
double const dt,
const int process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 240 of file HMPhaseFieldProcess.cpp.

243{
244 DBUG("PreTimestep HMPhaseFieldProcess {}.", process_id);
245
246 if (process_id == _process_data._phasefield_process_id)
247 {
248 DBUG("Store the value of phase field at previous time step.");
251 *x[process_id]);
252 }
253
257}
std::unique_ptr< GlobalVector > _x_previous_timestep
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
virtual void preTimestep(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)

References _local_assemblers, _process_data, _x_previous_timestep, DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), getDOFTable(), and ProcessLib::LocalAssemblerInterface::preTimestep().

◆ updateConstraints()

template<int DisplacementDim>
void ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::updateConstraints ( GlobalVector & lower,
GlobalVector & upper,
int const process_id )
overrideprivate

Definition at line 330 of file HMPhaseFieldProcess.cpp.

332{
333 lower.setZero();
336
337 GlobalIndexType const x_begin = _x_previous_timestep->getRangeBegin();
338 GlobalIndexType const x_end = _x_previous_timestep->getRangeEnd();
339
340 for (GlobalIndexType i = x_begin; i < x_end; i++)
341 {
342 if ((*_x_previous_timestep)[i] > _process_data.irreversible_threshold)
343 {
344 upper.set(i, 1.0);
345 }
346 }
347}
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:30
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:20

References _process_data, _x_previous_timestep, MathLib::LinAlg::copy(), MathLib::EigenVector::set(), MathLib::LinAlg::setLocalAccessibleVector(), and MathLib::EigenVector::setZero().

Member Data Documentation

◆ _local_assemblers

template<int DisplacementDim>
std::vector<std::unique_ptr<LocalAssemblerInterface> > ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::_local_assemblers
private

◆ _local_to_global_index_map_single_component

◆ _nodal_forces

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::_nodal_forces = nullptr
private

◆ _process_data

◆ _sparsity_pattern_with_single_component

template<int DisplacementDim>
GlobalSparsityPattern ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::_sparsity_pattern_with_single_component
private

Sparsity pattern for the phase field equation, and it is initialized only if the staggered scheme is used.

Definition at line 97 of file HMPhaseFieldProcess.h.

Referenced by constructDofTable(), and getMatrixSpecifications().

◆ _x_previous_timestep

template<int DisplacementDim>
std::unique_ptr<GlobalVector> ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::_x_previous_timestep
private

Definition at line 99 of file HMPhaseFieldProcess.h.

Referenced by preTimestepConcreteProcess(), and updateConstraints().


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