Loading [MathJax]/extensions/tex2jax.js
OGS
HMPhaseFieldProcess.cpp
Go to the documentation of this file.
1
11#include "HMPhaseFieldProcess.h"
12
13#include <cassert>
14
15#include "HMPhaseFieldFEM.h"
18#include "ProcessLib/Process.h"
20
21namespace ProcessLib
22{
23namespace HMPhaseField
24{
25template <int DisplacementDim>
27 std::string name,
28 MeshLib::Mesh& mesh,
29 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
30 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
31 unsigned const integration_order,
32 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
33 process_variables,
35 SecondaryVariableCollection&& secondary_variables,
36 bool const use_monolithic_scheme)
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}
52
53template <int DisplacementDim>
55{
56 return false;
57}
58
59template <int DisplacementDim>
62 const int process_id) const
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.
73 auto const& l = *_local_to_global_index_map_single_component;
74 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
75 &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
76}
77
78template <int DisplacementDim>
81{
82 // For the M process (deformation) in the staggered scheme.
83 if (process_id == _process_data._mechanics_related_process_id)
84 {
85 return *_local_to_global_index_map;
86 }
87
88 // For the H (hydro) or PF (phasefield) process in the staggered scheme.
89 return *_local_to_global_index_map_single_component;
90}
91
92template <int DisplacementDim>
94{
95 // For displacement equation.
96 constructDofTableOfSpecifiedProcessStaggeredScheme(
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{
102 *_mesh_subset_all_nodes};
103 _local_to_global_index_map_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
109 assert(_local_to_global_index_map_single_component);
110
111 // For the H (hydro) or PF (phasefield) process in the staggered scheme.
112 _sparsity_pattern_with_single_component = NumLib::computeSparsityPattern(
113 *_local_to_global_index_map_single_component, _mesh);
114}
115
116template <int DisplacementDim>
118 NumLib::LocalToGlobalIndexMap const& dof_table,
119 MeshLib::Mesh const& mesh,
120 unsigned const integration_order)
121{
123 DisplacementDim, HMPhaseFieldLocalAssembler>(
124 mesh.getElements(), dof_table, _local_assemblers,
125 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
126 _process_data);
127
128 _secondary_variables.addSecondaryVariable(
129 "sigma",
131 DisplacementDim>::RowsAtCompileTime,
132 getExtrapolator(), _local_assemblers,
133 &LocalAssemblerInterface::getIntPtSigma));
134
135 _secondary_variables.addSecondaryVariable(
136 "epsilon",
138 DisplacementDim>::RowsAtCompileTime,
139 getExtrapolator(), _local_assemblers,
140 &LocalAssemblerInterface::getIntPtEpsilon));
141
142 _secondary_variables.addSecondaryVariable(
143 "width_node",
144 makeExtrapolator(1, getExtrapolator(), _local_assemblers,
145 &LocalAssemblerInterface::getIntPtWidth));
146
147 _process_data.ele_d = MeshLib::getOrCreateMeshProperty<double>(
148 const_cast<MeshLib::Mesh&>(mesh), "damage", MeshLib::MeshItemType::Cell,
149 1);
150
151 _process_data.width = MeshLib::getOrCreateMeshProperty<double>(
152 const_cast<MeshLib::Mesh&>(mesh), "width", MeshLib::MeshItemType::Cell,
153 1);
154
155 // Initialize local assemblers after all variables have been set.
157 &LocalAssemblerInterface::initialize, _local_assemblers,
158 *_local_to_global_index_map);
159}
160
161template <int DisplacementDim>
163 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
164{
165 // Staggered scheme:
166 // for the phase field
167 initializeProcessBoundaryConditionsAndSourceTerms(
168 *_local_to_global_index_map_single_component,
169 _process_data._phasefield_process_id, media);
170 // for the pressure
171 initializeProcessBoundaryConditionsAndSourceTerms(
172 *_local_to_global_index_map_single_component,
173 _process_data._hydro_process_id, media);
174 // for the deformation.
175 initializeProcessBoundaryConditionsAndSourceTerms(
176 *_local_to_global_index_map,
177 _process_data._mechanics_related_process_id, media);
178}
179
180template <int DisplacementDim>
182 const double /*t*/, double const /*dt*/,
183 std::vector<GlobalVector*> const& /*x*/,
184 std::vector<GlobalVector*> const& /*x_prev*/, int const /*process_id*/,
185 GlobalMatrix& /*M*/, GlobalMatrix& /*K*/, GlobalVector& /*b*/)
186{
187 OGS_FATAL(
188 "HMPhaseFieldLocalAssembler: assembly for the Picard non-linear solver "
189 "is not implemented.");
190}
191
192template <int DisplacementDim>
194 const double t, double const dt, std::vector<GlobalVector*> const& x,
195 std::vector<GlobalVector*> const& x_prev, int const process_id,
196 GlobalVector& b, GlobalMatrix& Jac)
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 {
232 b.copyValues(*_nodal_forces);
233 std::transform(_nodal_forces->begin(), _nodal_forces->end(),
234 _nodal_forces->begin(), [](double val) { return -val; });
235 }
236}
237
238template <int DisplacementDim>
240 std::vector<GlobalVector*> const& x, double const t, double const dt,
241 const int process_id)
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.");
248 _x_previous_timestep =
250 *x[process_id]);
251 }
252
254 &LocalAssemblerInterface::preTimestep, _local_assemblers,
255 getActiveElementIDs(), getDOFTable(process_id), *x[process_id], t, dt);
256}
257
258template <int DisplacementDim>
260 std::vector<GlobalVector*> const& x,
261 std::vector<GlobalVector*> const& x_prev, const double t, const double dt,
262 int const process_id)
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(
275 _local_to_global_index_map_single_component.get());
276 dof_tables.emplace_back(
277 _local_to_global_index_map_single_component.get());
278 dof_tables.emplace_back(_local_to_global_index_map.get());
279
281 &LocalAssemblerInterface::computeEnergy, _local_assemblers,
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
292 &LocalAssemblerInterface::postTimestep, _local_assemblers,
293 getActiveElementIDs(), getDOFTables(x.size()), x, x_prev, t, dt,
294 process_id);
295}
296
297template <int DisplacementDim>
299 std::vector<GlobalVector*> const& x,
300 std::vector<GlobalVector*> const& x_prev, const double t, double const dt,
301 const int process_id)
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");
313 &LocalAssemblerInterface::approximateFractureWidth,
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}
327
328template <int DisplacementDim>
330 GlobalVector& lower, GlobalVector& upper, int const /*process_id*/)
331{
332 lower.setZero();
333 MathLib::LinAlg::setLocalAccessibleVector(*_x_previous_timestep);
334 MathLib::LinAlg::copy(*_x_previous_timestep, upper);
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}
347
348template class HMPhaseFieldProcess<2>;
349template class HMPhaseFieldProcess<3>;
350
351} // namespace HMPhaseField
352} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
GlobalMatrix::IndexType GlobalIndexType
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Global vector based on Eigen vector.
Definition EigenVector.h:25
void copyValues(std::vector< double > &u) const
void set(IndexType rowId, double v)
set entry
Definition EigenVector.h:73
bool isAxiallySymmetric() const
Definition Mesh.h:139
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:111
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)
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
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)
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:37
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:27
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)