OGS
HMPhaseFieldProcess.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
5
6#include <cassert>
7
8#include "HMPhaseFieldFEM.h"
11#include "ProcessLib/Process.h"
13
14namespace ProcessLib
15{
16namespace HMPhaseField
17{
18template <int DisplacementDim>
20 std::string name,
21 MeshLib::Mesh& mesh,
22 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
23 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
24 unsigned const integration_order,
25 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
26 process_variables,
28 SecondaryVariableCollection&& secondary_variables,
29 bool const use_monolithic_scheme)
30 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
31 integration_order, std::move(process_variables),
32 std::move(secondary_variables), use_monolithic_scheme),
33 _process_data(std::move(process_data))
34{
35 // For numerical Jacobian assembler
36 if (this->_jacobian_assembler->isPerturbationEnabled())
37 {
38 OGS_FATAL(
39 "The numericial Jacobian assembler is not supported for the "
40 "HMPhaseField process.");
41 }
42
43 if (use_monolithic_scheme)
44 {
46 "Monolithic scheme is not implemented for the HMPhaseField "
47 "process.");
48 }
49
51 mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
52}
53
54template <int DisplacementDim>
56{
57 return false;
58}
59
60template <int DisplacementDim>
63 const int process_id) const
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}
78
79template <int DisplacementDim>
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}
92
93template <int DisplacementDim>
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
102 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
105 std::make_unique<NumLib::LocalToGlobalIndexMap>(
106 std::move(all_mesh_subsets_single_component),
107 // by location order is needed for output
109
111
112 // For the H (hydro) or PF (phasefield) process in the staggered scheme.
115}
116
117template <int DisplacementDim>
119 NumLib::LocalToGlobalIndexMap const& dof_table,
120 MeshLib::Mesh const& mesh,
121 unsigned const integration_order)
122{
124 DisplacementDim, HMPhaseFieldLocalAssembler>(
125 mesh.getElements(), dof_table, _local_assemblers,
126 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
128
129 _secondary_variables.addSecondaryVariable(
130 "sigma",
132 DisplacementDim>::RowsAtCompileTime,
135
136 _secondary_variables.addSecondaryVariable(
137 "epsilon",
139 DisplacementDim>::RowsAtCompileTime,
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}
161
162template <int DisplacementDim>
164 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
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}
180
181template <int DisplacementDim>
183 const double /*t*/, double const /*dt*/,
184 std::vector<GlobalVector*> const& /*x*/,
185 std::vector<GlobalVector*> const& /*x_prev*/, int const /*process_id*/,
186 GlobalMatrix& /*M*/, GlobalMatrix& /*K*/, GlobalVector& /*b*/)
187{
188 OGS_FATAL(
189 "HMPhaseFieldLocalAssembler: assembly for the Picard non-linear solver "
190 "is not implemented.");
191}
192
193template <int DisplacementDim>
195 const double t, double const dt, std::vector<GlobalVector*> const& x,
196 std::vector<GlobalVector*> const& x_prev, int const process_id,
197 GlobalVector& b, GlobalMatrix& Jac)
198{
199 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
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
221 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
222 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
223 dof_tables.emplace_back(_local_to_global_index_map.get());
224
225 // Call global assembler for each local assembly item.
228 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
229 process_id, &b, &Jac);
230
231 if (process_id == _process_data._mechanics_related_process_id)
232 {
234 std::transform(_nodal_forces->begin(), _nodal_forces->end(),
235 _nodal_forces->begin(), [](double val) { return -val; });
236 }
237}
238
239template <int DisplacementDim>
241 std::vector<GlobalVector*> const& x, double const t, double const dt,
242 const int process_id)
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
256 getActiveElementIDs(), getDOFTable(process_id), *x[process_id], t, dt);
257}
258
259template <int DisplacementDim>
261 std::vector<GlobalVector*> const& x,
262 std::vector<GlobalVector*> const& x_prev, const double t, const double dt,
263 int const process_id)
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
273 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
274
275 dof_tables.emplace_back(
277 dof_tables.emplace_back(
279 dof_tables.emplace_back(_local_to_global_index_map.get());
280
283 getActiveElementIDs(), dof_tables, x, t,
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
294 getActiveElementIDs(), getDOFTables(x.size()), x, x_prev, t, dt,
295 process_id);
296}
297
298template <int DisplacementDim>
300 std::vector<GlobalVector*> const& x,
301 std::vector<GlobalVector*> const& x_prev, const double t, double const dt,
302 const int process_id)
303{
304 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
305
306 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
307 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
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");
315 _local_assemblers, dof_tables, x, t, dt);
316 }
317
318 if (process_id == _process_data._hydro_process_id)
319 {
320 INFO("PostNonLinearSolver for Hydro Process");
321
324 getActiveElementIDs(), getDOFTables(x.size()), x, x_prev, t, dt,
325 process_id);
326 }
327}
328
329template <int DisplacementDim>
331 GlobalVector& lower, GlobalVector& upper, int const /*process_id*/)
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}
348
349template class HMPhaseFieldProcess<2>;
350template class HMPhaseFieldProcess<3>;
351
352} // namespace HMPhaseField
353} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:19
MathLib::EigenMatrix GlobalMatrix
MathLib::EigenVector GlobalVector
GlobalMatrix::IndexType GlobalIndexType
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
void copyValues(std::vector< double > &u) const
void set(IndexType rowId, double v)
set entry
Definition EigenVector.h:67
bool isAxiallySymmetric() const
Definition Mesh.h:128
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:100
std::unique_ptr< GlobalVector > _x_previous_timestep
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &, const double t, const double dt, int const process_id) override
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)
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) 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 initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const 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
HMPhaseFieldProcessData< DisplacementDim > _process_data
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
void updateConstraints(GlobalVector &lower, GlobalVector &upper, int const process_id) override
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
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
MeshLib::PropertyVector< double > * _nodal_forces
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)
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)
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)
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
std::string const name
Definition Process.h:361
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:376
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_process_id)
Definition Process.cpp:345
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
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< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:365
MeshLib::Mesh & _mesh
Definition Process.h:364
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
SecondaryVariableCollection _secondary_variables
Definition Process.h:369
VectorMatrixAssembler _global_assembler
Definition Process.h:376
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:367
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:375
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:201
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:391
Handles configuration of several secondary variables from the project file.
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)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:30
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:20
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
@ 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.
void showEnergyAndWork(const double t, double &_elastic_energy, double &_surface_energy, double &_pressure_work)
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 executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
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 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
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
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