Loading [MathJax]/extensions/tex2jax.js
OGS
ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

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

Definition at line 22 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
 
- 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 48 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 26 of file HMPhaseFieldProcess.cpp.

37 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
38 integration_order, std::move(process_variables),
39 std::move(secondary_variables), use_monolithic_scheme),
40 _process_data(std::move(process_data))
41{
42 if (use_monolithic_scheme)
43 {
45 "Monolithic scheme is not implemented for the HMPhaseField "
46 "process.");
47 }
48
50 mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
51}
#define OGS_FATAL(...)
Definition Error.h:26
HMPhaseFieldProcessData< DisplacementDim > _process_data
MeshLib::PropertyVector< double > * _nodal_forces
std::string const name
Definition Process.h:362
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:44
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)

References ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::_nodal_forces, MeshLib::getOrCreateMeshProperty(), 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 181 of file HMPhaseFieldProcess.cpp.

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

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 193 of file HMPhaseFieldProcess.cpp.

197{
198 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
199
200 // For the staggered scheme
201 if (process_id == _process_data._phasefield_process_id)
202 {
203 DBUG(
204 "Assemble the Jacobian equations of phase field in "
205 "HMPhaseFieldProcess for the staggered scheme.");
206 }
207 else if (process_id == _process_data._hydro_process_id)
208 {
209 DBUG(
210 "Assemble the Jacobian equations of pressure in "
211 "HMPhaseFieldProcess for the staggered scheme.");
212 }
213 else
214 {
215 DBUG(
216 "Assemble the Jacobian equations of deformation in "
217 "HMPhaseFieldProcess for the staggered scheme.");
218 }
219
220 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
221 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
222 dof_tables.emplace_back(_local_to_global_index_map.get());
223
224 // Call global assembler for each local assembly item.
227 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
228 process_id, &b, &Jac);
229
230 if (process_id == _process_data._mechanics_related_process_id)
231 {
233 std::transform(_nodal_forces->begin(), _nodal_forces->end(),
234 _nodal_forces->begin(), [](double val) { return -val; });
235 }
236}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
void copyValues(std::vector< double > &u) const
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:167
VectorMatrixAssembler _global_assembler
Definition Process.h:377
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:368
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, 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)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), MathLib::EigenVector::copyValues(), DBUG(), and NumLib::SerialExecutor::executeSelectedMemberDereferenced().

◆ 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 93 of file HMPhaseFieldProcess.cpp.

94{
95 // For displacement equation.
97 _process_data._mechanics_related_process_id);
98
99 // TODO move the two data members somewhere else.
100 // for extrapolation of secondary variables of stress or strain
101 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
104 std::make_unique<NumLib::LocalToGlobalIndexMap>(
105 std::move(all_mesh_subsets_single_component),
106 // by location order is needed for output
108
110
111 // For the H (hydro) or PF (phasefield) process in the staggered scheme.
114}
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_process_id)
Definition Process.cpp:353
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:366
MeshLib::Mesh & _mesh
Definition Process.h:365
@ BY_LOCATION
Ordering data by spatial location.
GlobalSparsityPattern computeSparsityPattern(LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
Computes a sparsity pattern for the given inputs.

References NumLib::BY_LOCATION, and NumLib::computeSparsityPattern().

◆ 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 80 of file HMPhaseFieldProcess.cpp.

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

◆ getMatrixSpecifications()

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

Definition at line 61 of file HMPhaseFieldProcess.cpp.

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

◆ 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 162 of file HMPhaseFieldProcess.cpp.

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

◆ 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 117 of file HMPhaseFieldProcess.cpp.

121{
123 DisplacementDim, HMPhaseFieldLocalAssembler>(
124 mesh.getElements(), dof_table, _local_assemblers,
125 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
127
129 "sigma",
131 DisplacementDim>::RowsAtCompileTime,
134
136 "epsilon",
138 DisplacementDim>::RowsAtCompileTime,
141
143 "width_node",
146
148 const_cast<MeshLib::Mesh&>(mesh), "damage", MeshLib::MeshItemType::Cell,
149 1);
150
152 const_cast<MeshLib::Mesh&>(mesh), "width", MeshLib::MeshItemType::Cell,
153 1);
154
155 // Initialize local assemblers after all variables have been set.
159}
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
SecondaryVariableCollection _secondary_variables
Definition Process.h:370
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
void createLocalAssemblers(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, ProviderOrOrder const &provider_or_order, ExtraCtorArgs &&... extra_ctor_args)
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 MeshLib::Cell, ProcessLib::SmallDeformation::createLocalAssemblers(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getElements(), 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 54 of file HMPhaseFieldProcess.cpp.

55{
56 return false;
57}

◆ 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 298 of file HMPhaseFieldProcess.cpp.

302{
303 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
304
305 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
306 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
307 dof_tables.emplace_back(_local_to_global_index_map.get());
308
309 if (process_id == _process_data._phasefield_process_id)
310 {
311 INFO("Update fracture width and porous properties");
314 _local_assemblers, dof_tables, x, t, dt);
315 }
316
317 if (process_id == _process_data._hydro_process_id)
318 {
319 INFO("PostNonLinearSolver for Hydro Process");
320
323 getActiveElementIDs(), getDOFTables(x.size()), x, x_prev, t, dt,
324 process_id);
325 }
326}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
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:384
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 NumLib::SerialExecutor::executeMemberOnDereferenced(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), 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 259 of file HMPhaseFieldProcess.cpp.

263{
264 if (process_id == _process_data._phasefield_process_id)
265 {
266 DBUG("PostTimestep HMPhaseFieldProcess.");
267
268 _process_data.elastic_energy = 0.0;
269 _process_data.surface_energy = 0.0;
270 _process_data.pressure_work = 0.0;
271
272 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
273
274 dof_tables.emplace_back(
276 dof_tables.emplace_back(
278 dof_tables.emplace_back(_local_to_global_index_map.get());
279
282 getActiveElementIDs(), dof_tables, x, t,
283 _process_data.elastic_energy, _process_data.surface_energy,
284 _process_data.pressure_work);
285
286 showEnergyAndWork(t, _process_data.elastic_energy,
287 _process_data.surface_energy,
288 _process_data.pressure_work);
289 }
290
293 getActiveElementIDs(), getDOFTables(x.size()), x, x_prev, t, dt,
294 process_id);
295}
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 DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), 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 239 of file HMPhaseFieldProcess.cpp.

242{
243 DBUG("PreTimestep HMPhaseFieldProcess {}.", process_id);
244
245 if (process_id == _process_data._phasefield_process_id)
246 {
247 DBUG("Store the value of phase field at previous time step.");
250 *x[process_id]);
251 }
252
255 getActiveElementIDs(), getDOFTable(process_id), *x[process_id], t, dt);
256}
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 DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), 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 329 of file HMPhaseFieldProcess.cpp.

331{
332 lower.setZero();
335
336 GlobalIndexType const x_begin = _x_previous_timestep->getRangeBegin();
337 GlobalIndexType const x_end = _x_previous_timestep->getRangeEnd();
338
339 for (GlobalIndexType i = x_begin; i < x_end; i++)
340 {
341 if ((*_x_previous_timestep)[i] > _process_data.irreversible_threshold)
342 {
343 upper.set(i, 1.0);
344 }
345 }
346}
GlobalMatrix::IndexType GlobalIndexType
void set(IndexType rowId, double v)
set entry
Definition EigenVector.h:73
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:37
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:27

References 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

Definition at line 95 of file HMPhaseFieldProcess.h.

◆ _local_to_global_index_map_single_component

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::_local_to_global_index_map_single_component
private

Definition at line 98 of file HMPhaseFieldProcess.h.

◆ _nodal_forces

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

◆ _process_data

template<int DisplacementDim>
HMPhaseFieldProcessData<DisplacementDim> ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::_process_data
private

Definition at line 93 of file HMPhaseFieldProcess.h.

◆ _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 104 of file HMPhaseFieldProcess.h.

◆ _x_previous_timestep

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

Definition at line 106 of file HMPhaseFieldProcess.h.


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