26template <
int DisplacementDim>
30 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
31 std::vector<std::unique_ptr<ParameterLib::ParameterBase>>
const& parameters,
32 unsigned const integration_order,
33 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
37 bool const use_monolithic_scheme)
38 :
Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
39 integration_order, std::move(process_variables),
40 std::move(secondary_variables), use_monolithic_scheme),
41 _process_data(std::move(process_data))
43 if (use_monolithic_scheme)
46 "Monolithic scheme is not implemented for the PhaseField process.");
53 std::make_unique<MeshLib::IntegrationPointWriter>(
60template <
int DisplacementDim>
66template <
int DisplacementDim>
69 const int process_id)
const
74 auto const& l = *_local_to_global_index_map;
75 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
76 &l.getGhostIndices(), &this->_sparsity_pattern};
80 auto const& l = *_local_to_global_index_map_single_component;
81 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
82 &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
85template <
int DisplacementDim>
92 return *_local_to_global_index_map;
96 return *_local_to_global_index_map_single_component;
99template <
int DisplacementDim>
103 const int mechanics_process_id = 0;
104 constructDofTableOfSpecifiedProcessStaggeredScheme(mechanics_process_id);
108 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
109 *_mesh_subset_all_nodes};
110 _local_to_global_index_map_single_component =
111 std::make_unique<NumLib::LocalToGlobalIndexMap>(
112 std::move(all_mesh_subsets_single_component),
116 assert(_local_to_global_index_map_single_component);
120 *_local_to_global_index_map_single_component, _mesh);
123template <
int DisplacementDim>
127 unsigned const integration_order)
135 _secondary_variables.addSecondaryVariable(
138 DisplacementDim>::RowsAtCompileTime,
139 getExtrapolator(), _local_assemblers,
140 &LocalAssemblerInterface::getIntPtSigma));
142 _secondary_variables.addSecondaryVariable(
145 DisplacementDim>::RowsAtCompileTime,
146 getExtrapolator(), _local_assemblers,
147 &LocalAssemblerInterface::getIntPtEpsilon));
149 _secondary_variables.addSecondaryVariable(
152 DisplacementDim>::RowsAtCompileTime,
153 getExtrapolator(), _local_assemblers,
154 &LocalAssemblerInterface::getIntPtSigmaTensile));
156 _secondary_variables.addSecondaryVariable(
159 DisplacementDim>::RowsAtCompileTime,
160 getExtrapolator(), _local_assemblers,
161 &LocalAssemblerInterface::getIntPtSigmaCompressive));
163 _secondary_variables.addSecondaryVariable(
166 DisplacementDim>::RowsAtCompileTime,
167 getExtrapolator(), _local_assemblers,
168 &LocalAssemblerInterface::getIntPtEpsilonTensile));
176 *_local_to_global_index_map);
179template <
int DisplacementDim>
181 std::map<
int, std::shared_ptr<MaterialPropertyLib::Medium>>
const& media)
185 const int mechanical_process_id = 0;
186 initializeProcessBoundaryConditionsAndSourceTerms(
187 *_local_to_global_index_map, mechanical_process_id, media);
189 const int phasefield_process_id = 1;
190 initializeProcessBoundaryConditionsAndSourceTerms(
191 *_local_to_global_index_map_single_component, phasefield_process_id,
195template <
int DisplacementDim>
197 const double t,
double const dt, std::vector<GlobalVector*>
const& x,
198 std::vector<GlobalVector*>
const& x_prev,
int const process_id,
201 DBUG(
"Assemble PhaseFieldProcess.");
203 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
209 "Assemble the equations of phase field in "
210 "PhaseFieldProcess for the staggered scheme.");
215 "Assemble the equations of deformation in "
216 "PhaseFieldProcess for the staggered scheme.");
218 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
219 dof_tables.emplace_back(_local_to_global_index_map.get());
224 getActiveElementIDs(), dof_tables, t, dt, x, x_prev, process_id, &M, &K,
228template <
int DisplacementDim>
230 const double t,
double const dt, std::vector<GlobalVector*>
const& x,
231 std::vector<GlobalVector*>
const& x_prev,
int const process_id,
234 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
240 "Assemble the Jacobian equations of phase field in "
241 "PhaseFieldProcess for the staggered scheme.");
246 "Assemble the Jacobian equations of deformation in "
247 "PhaseFieldProcess for the staggered scheme.");
249 dof_tables.emplace_back(_local_to_global_index_map.get());
250 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
255 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
256 process_id, &b, &Jac);
261 std::transform(_nodal_forces->begin(), _nodal_forces->end(),
262 _nodal_forces->begin(), [](
double val) { return -val; });
266template <
int DisplacementDim>
268 std::vector<GlobalVector*>
const& x,
double const t,
double const dt,
269 const int process_id)
271 DBUG(
"PreTimestep PhaseFieldProcess {:d}.", process_id);
273 _process_data.injected_volume = t;
275 _x_previous_timestep =
280 getActiveElementIDs(), getDOFTable(process_id), *x[process_id], t, dt);
283template <
int DisplacementDim>
285 std::vector<GlobalVector*>
const& x,
286 std::vector<GlobalVector*>
const& ,
const double t,
287 const double ,
int const process_id)
289 if (isPhaseFieldProcess(process_id))
291 DBUG(
"PostTimestep PhaseFieldProcess.");
293 _process_data.elastic_energy = 0.0;
294 _process_data.surface_energy = 0.0;
295 _process_data.pressure_work = 0.0;
297 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
299 dof_tables.emplace_back(_local_to_global_index_map.get());
300 dof_tables.emplace_back(
301 _local_to_global_index_map_single_component.get());
304 &LocalAssemblerInterface::computeEnergy, _local_assemblers,
305 getActiveElementIDs(), dof_tables, x, t,
306 _process_data.elastic_energy, _process_data.surface_energy,
307 _process_data.pressure_work);
310 double const elastic_energy = _process_data.elastic_energy;
311 MPI_Allreduce(&elastic_energy, &_process_data.elastic_energy, 1,
312 MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
313 double const surface_energy = _process_data.surface_energy;
314 MPI_Allreduce(&surface_energy, &_process_data.surface_energy, 1,
315 MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
316 double const pressure_work = _process_data.pressure_work;
317 MPI_Allreduce(&pressure_work, &_process_data.pressure_work, 1,
318 MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
322 "Elastic energy: {:g} Surface energy: {:g} Pressure work: {:g} at "
324 _process_data.elastic_energy, _process_data.surface_energy,
325 _process_data.pressure_work, t);
326 if (_process_data.propagating_pressurized_crack)
328 INFO(
"Pressure: {:g} at time: {:g} ", _process_data.pressure, t);
333template <
int DisplacementDim>
335 std::vector<GlobalVector*>
const& x,
336 std::vector<GlobalVector*>
const& ,
const double t,
337 double const ,
const int process_id)
339 _process_data.crack_volume = 0.0;
341 if (isPhaseFieldProcess(process_id))
343 if (_process_data.propagating_pressurized_crack)
347 1 / _process_data.pressure);
352 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
354 dof_tables.emplace_back(_local_to_global_index_map.get());
355 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
357 DBUG(
"PostNonLinearSolver crack volume computation.");
360 &LocalAssemblerInterface::computeCrackIntegral, _local_assemblers,
361 getActiveElementIDs(), dof_tables, x, t, _process_data.crack_volume);
364 double const crack_volume = _process_data.crack_volume;
365 MPI_Allreduce(&crack_volume, &_process_data.crack_volume, 1, MPI_DOUBLE,
366 MPI_SUM, PETSC_COMM_WORLD);
369 INFO(
"Integral of crack: {:g}", _process_data.crack_volume);
371 if (_process_data.propagating_pressurized_crack)
373 _process_data.pressure_old = _process_data.pressure;
374 _process_data.pressure =
375 _process_data.injected_volume / _process_data.crack_volume;
376 _process_data.pressure_error =
377 std::abs(_process_data.pressure_old - _process_data.pressure) /
378 _process_data.pressure;
379 INFO(
"Internal pressure: {:g} and Pressure error: {:.4e}",
380 _process_data.pressure, _process_data.pressure_error);
384 _process_data.pressure);
388template <
int DisplacementDim>
401 if ((*_x_previous_timestep)[i] > _process_data.irreversible_threshold)
408template <
int DisplacementDim>
410 int const process_id)
const
412 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)
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