Loading [MathJax]/extensions/MathZoom.js
OGS
HMPhaseFieldProcessData.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <Eigen/Core>
14#include <memory>
15#include <utility>
16#ifdef USE_PETSC
17#include <petscsys.h>
18#endif
19
24
25namespace MaterialLib
26{
27namespace Solids
28{
29template <int DisplacementDim>
30struct MechanicsBase;
31}
32} // namespace MaterialLib
33namespace ProcessLib
34{
35template <typename T>
36struct Parameter;
37
38namespace HMPhaseField
39{
40inline void showEnergyAndWork(const double t, double& _elastic_energy,
41 double& _surface_energy, double& _pressure_work)
42{
43#ifdef USE_PETSC
44 double const elastic_energy = _elastic_energy;
45 MPI_Allreduce(&elastic_energy, &_elastic_energy, 1, MPI_DOUBLE, MPI_SUM,
46 PETSC_COMM_WORLD);
47 double const surface_energy = _surface_energy;
48 MPI_Allreduce(&surface_energy, &_surface_energy, 1, MPI_DOUBLE, MPI_SUM,
49 PETSC_COMM_WORLD);
50 double const pressure_work = _pressure_work;
51 MPI_Allreduce(&pressure_work, &_pressure_work, 1, MPI_DOUBLE, MPI_SUM,
52 PETSC_COMM_WORLD);
53#endif
54
55 INFO("Elastic energy: {} Surface energy: {} Pressure work: {} at time: {} ",
56 _elastic_energy, _surface_energy, _pressure_work, t);
57};
58template <int DisplacementDim>
60{
62
64
65 std::map<int, std::shared_ptr<
71 Eigen::Matrix<double, DisplacementDim, 1> const specific_body_force;
72 Eigen::Matrix<double, DisplacementDim, 1> const specific_fracture_direction;
78 std::unique_ptr<MaterialLib::Solids::Phasefield::DegradationDerivative>
82 double const fracture_threshold;
87
90
91 double elastic_energy = 0.0;
92 double surface_energy = 0.0;
93 double pressure_work = 0.0;
95 int const _hydro_process_id = 1;
97};
98
99} // namespace HMPhaseField
100} // namespace ProcessLib
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void showEnergyAndWork(const double t, double &_elastic_energy, double &_surface_energy, double &_pressure_work)
ParameterLib::Parameter< double > const & crack_length_scale
ParameterLib::Parameter< double > const & residual_stiffness
MaterialLib::Solids::Phasefield::PhaseFieldModel phasefield_model
MaterialLib::Solids::Phasefield::EnergySplitModel energy_split_model
ParameterLib::Parameter< double > const & width_init
ParameterLib::Parameter< double > const & crack_resistance
Eigen::Matrix< double, DisplacementDim, 1 > const specific_fracture_direction
Eigen::Matrix< double, DisplacementDim, 1 > const specific_body_force
MaterialLib::Solids::Phasefield::SofteningCurve softening_curve
std::map< int, std::shared_ptr< MaterialLib::Solids::MechanicsBase< DisplacementDim > > > solid_materials
std::unique_ptr< MaterialLib::Solids::Phasefield::DegradationDerivative > degradation_derivative
MeshLib::PropertyVector< int > const *const material_ids
MaterialPropertyLib::MaterialSpatialDistributionMap media_map