OGS
PhaseFieldProcess.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
4#include "PhaseFieldProcess.h"
5
6#include <cassert>
7
8#include "BaseLib/MPI.h"
11#include "PhaseFieldFEM.h"
12#include "ProcessLib/Process.h"
15
16namespace ProcessLib
17{
18namespace PhaseField
19{
20template <int DisplacementDim>
22 std::string name,
23 MeshLib::Mesh& mesh,
24 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
25 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
26 unsigned const integration_order,
27 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
28 process_variables,
30 SecondaryVariableCollection&& secondary_variables,
31 bool const use_monolithic_scheme)
32 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
33 integration_order, std::move(process_variables),
34 std::move(secondary_variables), use_monolithic_scheme),
35 _process_data(std::move(process_data))
36{
37 // For numerical Jacobian assembler
38 if (this->_jacobian_assembler->isPerturbationEnabled())
39 {
40 OGS_FATAL(
41 "The numericial Jacobian assembler is not supported for the "
42 "PhaseField process.");
43 }
44
45 if (use_monolithic_scheme)
46 {
48 "Monolithic scheme is not implemented for the PhaseField process.");
49 }
50
52 mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
53
54 _integration_point_writer.emplace_back(
55 std::make_unique<MeshLib::IntegrationPointWriter>(
56 "sigma_ip",
57 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
58 integration_order, _local_assemblers,
60}
61
62template <int DisplacementDim>
64{
65 return false;
66}
67
68template <int DisplacementDim>
71 const int process_id) const
72{
73 // For the M process (deformation) in the staggered scheme.
74 if (process_id == 0)
75 {
76 auto const& l = *_local_to_global_index_map;
77 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
78 &l.getGhostIndices(), &this->_sparsity_pattern};
79 }
80
81 // For staggered scheme and phase field process.
83 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
84 &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
85}
86
87template <int DisplacementDim>
90{
91 // For the M process (deformation) in the staggered scheme.
92 if (process_id == 0)
93 {
95 }
96
97 // For the equation of phasefield
99}
100
101template <int DisplacementDim>
103{
104 // For displacement equation.
105 const int mechanics_process_id = 0;
107
108 // TODO move the two data members somewhere else.
109 // for extrapolation of secondary variables of stress or strain
110 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
113 std::make_unique<NumLib::LocalToGlobalIndexMap>(
114 std::move(all_mesh_subsets_single_component),
115 // by location order is needed for output
117
119
120 // For phase field equation.
123}
124
125template <int DisplacementDim>
127 NumLib::LocalToGlobalIndexMap const& dof_table,
128 MeshLib::Mesh const& mesh,
129 unsigned const integration_order)
130{
132 DisplacementDim, PhaseFieldLocalAssembler>(
133 mesh.getElements(), dof_table, _local_assemblers,
134 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
136
137 _secondary_variables.addSecondaryVariable(
138 "sigma",
140 DisplacementDim>::RowsAtCompileTime,
143
144 _secondary_variables.addSecondaryVariable(
145 "epsilon",
147 DisplacementDim>::RowsAtCompileTime,
150
151 _secondary_variables.addSecondaryVariable(
152 "sigma_tensile",
154 DisplacementDim>::RowsAtCompileTime,
157
158 _secondary_variables.addSecondaryVariable(
159 "sigma_compressive",
161 DisplacementDim>::RowsAtCompileTime,
164
165 _secondary_variables.addSecondaryVariable(
166 "eps_tensile",
168 DisplacementDim>::RowsAtCompileTime,
171
174
175 // Initialize local assemblers after all variables have been set.
179}
180
181template <int DisplacementDim>
183 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
184{
185 // Staggered scheme:
186 // for the equations of deformation.
187 const int mechanical_process_id = 0;
189 *_local_to_global_index_map, mechanical_process_id, media);
190 // for the phase field
191 const int phasefield_process_id = 1;
193 *_local_to_global_index_map_single_component, phasefield_process_id,
194 media);
195}
196
197template <int DisplacementDim>
199 const double t, double const dt, std::vector<GlobalVector*> const& x,
200 std::vector<GlobalVector*> const& x_prev, int const process_id,
202{
203 DBUG("Assemble PhaseFieldProcess.");
204
205 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
206
207 // For the staggered scheme
208 if (process_id == 1)
209 {
210 DBUG(
211 "Assemble the equations of phase field in "
212 "PhaseFieldProcess for the staggered scheme.");
213 }
214 else
215 {
216 DBUG(
217 "Assemble the equations of deformation in "
218 "PhaseFieldProcess for the staggered scheme.");
219 }
220 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
221 dof_tables.emplace_back(_local_to_global_index_map.get());
222
223 // Call global assembler for each local assembly item.
226 getActiveElementIDs(), dof_tables, t, dt, x, x_prev, process_id, &M, &K,
227 &b);
228}
229
230template <int DisplacementDim>
232 const double t, double const dt, std::vector<GlobalVector*> const& x,
233 std::vector<GlobalVector*> const& x_prev, int const process_id,
234 GlobalVector& b, GlobalMatrix& Jac)
235{
236 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
237
238 // For the staggered scheme
239 if (process_id == 1)
240 {
241 DBUG(
242 "Assemble the Jacobian equations of phase field in "
243 "PhaseFieldProcess for the staggered scheme.");
244 }
245 else
246 {
247 DBUG(
248 "Assemble the Jacobian equations of deformation in "
249 "PhaseFieldProcess for the staggered scheme.");
250 }
251 dof_tables.emplace_back(_local_to_global_index_map.get());
252 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
253
254 // Call global assembler for each local assembly item.
257 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
258 process_id, &b, &Jac);
259
260 if (process_id == 0)
261 {
263 std::transform(_nodal_forces->begin(), _nodal_forces->end(),
264 _nodal_forces->begin(), [](double val) { return -val; });
265 }
266}
267
268template <int DisplacementDim>
270 std::vector<GlobalVector*> const& x, double const t, double const dt,
271 const int process_id)
272{
273 DBUG("PreTimestep PhaseFieldProcess {:d}.", process_id);
274
275 _process_data.injected_volume = t;
276
279
282 getActiveElementIDs(), getDOFTable(process_id), *x[process_id], t, dt);
283}
284
285template <int DisplacementDim>
287 std::vector<GlobalVector*> const& x,
288 std::vector<GlobalVector*> const& /*x_prev*/, const double t,
289 const double /*delta_t*/, int const process_id)
290{
291 if (isPhaseFieldProcess(process_id))
292 {
293 DBUG("PostTimestep PhaseFieldProcess.");
294
295 _process_data.elastic_energy = 0.0;
296 _process_data.surface_energy = 0.0;
297 _process_data.pressure_work = 0.0;
298
299 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
300
301 dof_tables.emplace_back(_local_to_global_index_map.get());
302 dof_tables.emplace_back(
304
307 getActiveElementIDs(), dof_tables, x, t,
308 _process_data.elastic_energy, _process_data.surface_energy,
309 _process_data.pressure_work);
310
311#ifdef USE_PETSC
312 BaseLib::MPI::Mpi mpi{};
313 _process_data.elastic_energy =
314 BaseLib::MPI::allreduce(_process_data.elastic_energy, MPI_SUM, mpi);
315 _process_data.surface_energy =
316 BaseLib::MPI::allreduce(_process_data.surface_energy, MPI_SUM, mpi);
317 _process_data.pressure_work =
318 BaseLib::MPI::allreduce(_process_data.pressure_work, MPI_SUM, mpi);
319#endif
320
321 INFO(
322 "Elastic energy: {:g} Surface energy: {:g} Pressure work: {:g} at "
323 "time: {:g} ",
324 _process_data.elastic_energy, _process_data.surface_energy,
325 _process_data.pressure_work, t);
326 if (_process_data.propagating_pressurized_crack)
327 {
328 INFO("Pressure: {:g} at time: {:g} ", _process_data.pressure, t);
329 }
330 }
331}
332
333template <int DisplacementDim>
335 std::vector<GlobalVector*> const& x,
336 std::vector<GlobalVector*> const& /*x_prev*/, const double t,
337 double const /*dt*/, const int process_id)
338{
339 _process_data.crack_volume = 0.0;
340
341 if (isPhaseFieldProcess(process_id))
342 {
343 if (_process_data.propagating_pressurized_crack)
344 {
345 auto& u = *x[0];
346 MathLib::LinAlg::scale(const_cast<GlobalVector&>(u),
347 1 / _process_data.pressure);
348 }
349 return;
350 }
351
352 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
353
354 dof_tables.emplace_back(_local_to_global_index_map.get());
355 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
356
357 DBUG("PostNonLinearSolver crack volume computation.");
358
361 getActiveElementIDs(), dof_tables, x, t, _process_data.crack_volume);
362
363#ifdef USE_PETSC
365 _process_data.crack_volume, MPI_SUM, BaseLib::MPI::Mpi{});
366#endif
367
368 INFO("Integral of crack: {:g}", _process_data.crack_volume);
369
370 if (_process_data.propagating_pressurized_crack)
371 {
372 _process_data.pressure_old = _process_data.pressure;
373 _process_data.pressure =
374 _process_data.injected_volume / _process_data.crack_volume;
375 _process_data.pressure_error =
376 std::abs(_process_data.pressure_old - _process_data.pressure) /
377 _process_data.pressure;
378 INFO("Internal pressure: {:g} and Pressure error: {:.4e}",
379 _process_data.pressure, _process_data.pressure_error);
380
381 auto& u = *x[0];
382 MathLib::LinAlg::scale(const_cast<GlobalVector&>(u),
383 _process_data.pressure);
384 }
385}
386
387template <int DisplacementDim>
389 GlobalVector& lower, GlobalVector& upper, int const /*process_id*/)
390{
391 lower.setZero();
394
395 GlobalIndexType const x_begin = _x_previous_timestep->getRangeBegin();
396 GlobalIndexType const x_end = _x_previous_timestep->getRangeEnd();
397
398 for (GlobalIndexType i = x_begin; i < x_end; i++)
399 {
400 if ((*_x_previous_timestep)[i] > _process_data.irreversible_threshold)
401 {
402 upper.set(i, 1.0);
403 }
404 }
405}
406
407template <int DisplacementDim>
409 int const process_id) const
410{
411 return process_id == 1;
412}
413
414template class PhaseFieldProcess<2>;
415template class PhaseFieldProcess<3>;
416
417} // namespace PhaseField
418} // 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
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:79
Properties & getProperties()
Definition Mesh.h:125
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)
bool isPhaseFieldProcess(int const process_id) const
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
void updateConstraints(GlobalVector &lower, GlobalVector &upper, int const process_id) override
std::unique_ptr< GlobalVector > _x_previous_timestep
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
MeshLib::PropertyVector< double > * _nodal_forces
void postNonLinearSolverConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id) override
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const 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
GlobalSparsityPattern _sparsity_pattern_with_single_component
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
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().
PhaseFieldProcessData< DisplacementDim > _process_data
PhaseFieldProcess(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, PhaseFieldProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
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
std::string const name
Definition Process.h:361
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:389
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 assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const 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)
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 T allreduce(T const &value, MPI_Op const &mpi_op, Mpi const &mpi)
Definition MPI.h:122
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
void scale(PETScVector &x, PetscScalar const a)
Definition LinAlg.cpp:37
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 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)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
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 & getIntPtSigmaTensile(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 & getIntPtEpsilonTensile(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
virtual std::vector< double > getSigma() const =0
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 computeCrackIntegral(std::size_t mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &crack_volume)=0
virtual std::vector< double > const & getIntPtSigmaCompressive(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =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