27template <
int DisplacementDim>
31 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
32 std::vector<std::unique_ptr<ParameterLib::ParameterBase>>
const& parameters,
33 unsigned const integration_order,
34 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
38 bool const use_monolithic_scheme)
39 :
Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
40 integration_order, std::move(process_variables),
41 std::move(secondary_variables), use_monolithic_scheme),
42 _process_data(std::move(process_data))
44 if (use_monolithic_scheme)
47 "Monolithic scheme is not implemented for the PhaseField process.");
54 std::make_unique<MeshLib::IntegrationPointWriter>(
61template <
int DisplacementDim>
67template <
int DisplacementDim>
70 const int process_id)
const
75 auto const& l = *_local_to_global_index_map;
76 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
77 &l.getGhostIndices(), &this->_sparsity_pattern};
81 auto const& l = *_local_to_global_index_map_single_component;
82 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
83 &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
86template <
int DisplacementDim>
93 return *_local_to_global_index_map;
97 return *_local_to_global_index_map_single_component;
100template <
int DisplacementDim>
104 const int mechanics_process_id = 0;
105 constructDofTableOfSpecifiedProcessStaggeredScheme(mechanics_process_id);
109 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
110 *_mesh_subset_all_nodes};
111 _local_to_global_index_map_single_component =
112 std::make_unique<NumLib::LocalToGlobalIndexMap>(
113 std::move(all_mesh_subsets_single_component),
117 assert(_local_to_global_index_map_single_component);
121 *_local_to_global_index_map_single_component, _mesh);
124template <
int DisplacementDim>
128 unsigned const integration_order)
136 _secondary_variables.addSecondaryVariable(
139 DisplacementDim>::RowsAtCompileTime,
140 getExtrapolator(), _local_assemblers,
141 &LocalAssemblerInterface::getIntPtSigma));
143 _secondary_variables.addSecondaryVariable(
146 DisplacementDim>::RowsAtCompileTime,
147 getExtrapolator(), _local_assemblers,
148 &LocalAssemblerInterface::getIntPtEpsilon));
150 _secondary_variables.addSecondaryVariable(
153 DisplacementDim>::RowsAtCompileTime,
154 getExtrapolator(), _local_assemblers,
155 &LocalAssemblerInterface::getIntPtSigmaTensile));
157 _secondary_variables.addSecondaryVariable(
160 DisplacementDim>::RowsAtCompileTime,
161 getExtrapolator(), _local_assemblers,
162 &LocalAssemblerInterface::getIntPtSigmaCompressive));
164 _secondary_variables.addSecondaryVariable(
167 DisplacementDim>::RowsAtCompileTime,
168 getExtrapolator(), _local_assemblers,
169 &LocalAssemblerInterface::getIntPtEpsilonTensile));
177 *_local_to_global_index_map);
180template <
int DisplacementDim>
182 std::map<
int, std::shared_ptr<MaterialPropertyLib::Medium>>
const& media)
186 const int mechanical_process_id = 0;
187 initializeProcessBoundaryConditionsAndSourceTerms(
188 *_local_to_global_index_map, mechanical_process_id, media);
190 const int phasefield_process_id = 1;
191 initializeProcessBoundaryConditionsAndSourceTerms(
192 *_local_to_global_index_map_single_component, phasefield_process_id,
196template <
int DisplacementDim>
198 const double t,
double const dt, std::vector<GlobalVector*>
const& x,
199 std::vector<GlobalVector*>
const& x_prev,
int const process_id,
202 DBUG(
"Assemble PhaseFieldProcess.");
204 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
210 "Assemble the equations of phase field in "
211 "PhaseFieldProcess for the staggered scheme.");
216 "Assemble the equations of deformation in "
217 "PhaseFieldProcess for the staggered scheme.");
219 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
220 dof_tables.emplace_back(_local_to_global_index_map.get());
225 getActiveElementIDs(), dof_tables, t, dt, x, x_prev, process_id, &M, &K,
229template <
int DisplacementDim>
231 const double t,
double const dt, std::vector<GlobalVector*>
const& x,
232 std::vector<GlobalVector*>
const& x_prev,
int const process_id,
235 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
241 "Assemble the Jacobian equations of phase field in "
242 "PhaseFieldProcess for the staggered scheme.");
247 "Assemble the Jacobian equations of deformation in "
248 "PhaseFieldProcess for the staggered scheme.");
250 dof_tables.emplace_back(_local_to_global_index_map.get());
251 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
256 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
257 process_id, &b, &Jac);
262 std::transform(_nodal_forces->begin(), _nodal_forces->end(),
263 _nodal_forces->begin(), [](
double val) { return -val; });
267template <
int DisplacementDim>
269 std::vector<GlobalVector*>
const& x,
double const t,
double const dt,
270 const int process_id)
272 DBUG(
"PreTimestep PhaseFieldProcess {:d}.", process_id);
274 _process_data.injected_volume = t;
276 _x_previous_timestep =
281 getActiveElementIDs(), getDOFTable(process_id), *x[process_id], t, dt);
284template <
int DisplacementDim>
286 std::vector<GlobalVector*>
const& x,
287 std::vector<GlobalVector*>
const& ,
const double t,
288 const double ,
int const process_id)
290 if (isPhaseFieldProcess(process_id))
292 DBUG(
"PostTimestep PhaseFieldProcess.");
294 _process_data.elastic_energy = 0.0;
295 _process_data.surface_energy = 0.0;
296 _process_data.pressure_work = 0.0;
298 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
300 dof_tables.emplace_back(_local_to_global_index_map.get());
301 dof_tables.emplace_back(
302 _local_to_global_index_map_single_component.get());
305 &LocalAssemblerInterface::computeEnergy, _local_assemblers,
306 getActiveElementIDs(), dof_tables, x, t,
307 _process_data.elastic_energy, _process_data.surface_energy,
308 _process_data.pressure_work);
312 _process_data.elastic_energy =
314 _process_data.surface_energy =
316 _process_data.pressure_work =
321 "Elastic energy: {:g} Surface energy: {:g} Pressure work: {:g} at "
323 _process_data.elastic_energy, _process_data.surface_energy,
324 _process_data.pressure_work, t);
325 if (_process_data.propagating_pressurized_crack)
327 INFO(
"Pressure: {:g} at time: {:g} ", _process_data.pressure, t);
332template <
int DisplacementDim>
334 std::vector<GlobalVector*>
const& x,
335 std::vector<GlobalVector*>
const& ,
const double t,
336 double const ,
const int process_id)
338 _process_data.crack_volume = 0.0;
340 if (isPhaseFieldProcess(process_id))
342 if (_process_data.propagating_pressurized_crack)
346 1 / _process_data.pressure);
351 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
353 dof_tables.emplace_back(_local_to_global_index_map.get());
354 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
356 DBUG(
"PostNonLinearSolver crack volume computation.");
359 &LocalAssemblerInterface::computeCrackIntegral, _local_assemblers,
360 getActiveElementIDs(), dof_tables, x, t, _process_data.crack_volume);
367 INFO(
"Integral of crack: {:g}", _process_data.crack_volume);
369 if (_process_data.propagating_pressurized_crack)
371 _process_data.pressure_old = _process_data.pressure;
372 _process_data.pressure =
373 _process_data.injected_volume / _process_data.crack_volume;
374 _process_data.pressure_error =
375 std::abs(_process_data.pressure_old - _process_data.pressure) /
376 _process_data.pressure;
377 INFO(
"Internal pressure: {:g} and Pressure error: {:.4e}",
378 _process_data.pressure, _process_data.pressure_error);
382 _process_data.pressure);
386template <
int DisplacementDim>
399 if ((*_x_previous_timestep)[i] > _process_data.irreversible_threshold)
406template <
int DisplacementDim>
408 int const process_id)
const
410 return process_id == 1;
GlobalMatrix::IndexType GlobalIndexType
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Global vector based on Eigen vector.
void copyValues(std::vector< double > &u) const
void set(IndexType rowId, double v)
set entry
bool isAxiallySymmetric() const
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Properties & getProperties()
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
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
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 constructDofTable() 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
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
bool isLinear() const 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().
PhaseFieldProcess(std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const ¶meters, 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::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
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)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
void copy(PETScVector const &x, PETScVector &y)
void setLocalAccessibleVector(PETScVector const &x)
void scale(PETScVector &x, PetscScalar const a)
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.
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 > getSigma() const =0