OGS
ProcessLib::HMPhaseField Namespace Reference

Classes

class  HMPhaseFieldLocalAssembler
struct  HMPhaseFieldLocalAssemblerInterface
class  HMPhaseFieldProcess
struct  HMPhaseFieldProcessData
struct  IntegrationPointData
struct  SecondaryData

Functions

void checkMPLProperties (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template<int DisplacementDim>
std::unique_ptr< ProcesscreateHMPhaseFieldProcess (std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< ProcessVariable > const &variables, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< ProcesscreateHMPhaseFieldProcess< 2 > (std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< ProcessVariable > const &variables, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< ProcesscreateHMPhaseFieldProcess< 3 > (std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< ProcessVariable > const &variables, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
void showEnergyAndWork (const double t, double &_elastic_energy, double &_surface_energy, double &_pressure_work)

Function Documentation

◆ checkMPLProperties()

void ProcessLib::HMPhaseField::checkMPLProperties ( std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media)

Definition at line 23 of file CreateHMPhaseFieldProcess.cpp.

25{
26 std::array const requiredMediumProperties = {
29 std::array const requiredFluidProperties = {MaterialPropertyLib::viscosity,
31 std::array const requiredSolidProperties = {
33
34 for (auto const& m : media)
35 {
36 checkRequiredProperties(*m.second, requiredMediumProperties);
37 checkRequiredProperties(m.second->phase("AqueousLiquid"),
38 requiredFluidProperties);
39 checkRequiredProperties(m.second->phase("Solid"),
40 requiredSolidProperties);
41 }
42}
void checkRequiredProperties(Component const &c, std::span< PropertyType const > const required_properties)
Definition Component.cpp:51

References MaterialPropertyLib::biot_coefficient, MaterialPropertyLib::density, MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, MaterialPropertyLib::reference_temperature, and MaterialPropertyLib::viscosity.

Referenced by createHMPhaseFieldProcess().

◆ createHMPhaseFieldProcess()

template<int DisplacementDim>
std::unique_ptr< Process > ProcessLib::HMPhaseField::createHMPhaseFieldProcess ( std::string name,
MeshLib::Mesh & mesh,
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && jacobian_assembler,
std::vector< ProcessVariable > const & variables,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
std::optional< ParameterLib::CoordinateSystem > const & local_coordinate_system,
unsigned const integration_order,
BaseLib::ConfigTree const & config,
std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media )
Input File Parameter
prj__processes__process__type
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__coupling_scheme
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__process_variables

Primary process variables as they appear in the global component vector:

Input File Parameter
prj__processes__process__HM_PHASE_FIELD__process_variables__phasefield
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__process_variables__pressure
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__process_variables__displacement

Process Parameters

Input File Parameter
prj__processes__process__HM_PHASE_FIELD__phasefield_parameters
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__phasefield_parameters__residual_stiffness
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__phasefield_parameters__crack_resistance
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__phasefield_parameters__crack_length_scale
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__characteristic_length
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__fluid_compressibility
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__specific_body_force
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__specific_fracture_direction
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__irreversible_threshold
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__phasefield_model
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__softening_curve
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__energy_split_model
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__diffused_range_parameter
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__fracture_threshold
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__fracture_permeability_parameter
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__fixed_stress_stabilization_parameter
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__spatial_stabilization_parameter
Input File Parameter
prj__processes__process__HM_PHASE_FIELD__width_init

Definition at line 45 of file CreateHMPhaseFieldProcess.cpp.

54{
56 config.checkConfigParameter("type", "HM_PHASE_FIELD");
57 DBUG("Create HMPhaseFieldProcess.");
58
59 auto const coupling_scheme =
61 config.getConfigParameter<std::string>("coupling_scheme", "staggered");
62 const bool use_monolithic_scheme = (coupling_scheme != "staggered");
63
65 auto const pv_config = config.getConfigSubtree("process_variables");
66
67 ProcessVariable* variable_ph;
68 ProcessVariable* variable_p;
69 ProcessVariable* variable_u;
70 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>
71 process_variables;
72 if (use_monolithic_scheme) // monolithic scheme.
73 {
74 OGS_FATAL("Monolithic implementation is not available.");
75 }
76 else // staggered scheme.
77 {
78 using namespace std::string_literals;
79 for (
82 auto const& variable_name :
83 {
84 "phasefield"s,
86 "pressure"s,
88 "displacement"s})
89 {
90 auto per_process_variables =
91 findProcessVariables(variables, pv_config, {variable_name});
92 process_variables.push_back(std::move(per_process_variables));
93 }
94 variable_ph = &process_variables[0][0].get();
95 variable_p = &process_variables[1][0].get();
96 variable_u = &process_variables[2][0].get();
97 }
98
99 DBUG("Associate phase field with process variable '{}'.",
100 variable_ph->getName());
101 if (variable_ph->getNumberOfGlobalComponents() != 1)
102 {
103 OGS_FATAL(
104 "HMPhasefield process variable '{}' is not a scalar variable but "
105 "has {:d} components.",
106 variable_ph->getName(),
107 variable_ph->getNumberOfGlobalComponents());
108 }
109
110 DBUG("Associate pressure with process variable '{}'.",
111 variable_p->getName());
112 if (variable_p->getNumberOfGlobalComponents() != 1)
113 {
114 OGS_FATAL(
115 "HMPhasefield process variable '{}' is not a scalar variable but "
116 "has {:d} components.",
117 variable_p->getName(),
118 variable_p->getNumberOfGlobalComponents());
119 }
120
121 DBUG("Associate displacement with process variable '{}'.",
122 variable_u->getName());
123
124 if (variable_u->getNumberOfGlobalComponents() != DisplacementDim)
125 {
126 OGS_FATAL(
127 "The component number of the process variable '{}' is different "
128 "from the displacement dimension: got {:d}, expected {:d}",
129 variable_u->getName(),
130 variable_u->getNumberOfGlobalComponents(),
131 DisplacementDim);
132 }
133
135 auto solid_constitutive_relations =
137 parameters, local_coordinate_system, materialIDs(mesh), config);
138
139 auto const phasefield_parameters_config =
141 config.getConfigSubtree("phasefield_parameters");
142
143 // Residual stiffness
144 auto const& residual_stiffness = ParameterLib::findParameter<double>(
145 phasefield_parameters_config,
147 "residual_stiffness", parameters, 1);
148 DBUG("Use '{}' as residual stiffness.", residual_stiffness.name);
149
150 // Crack resistance
151 auto const& crack_resistance = ParameterLib::findParameter<double>(
152 phasefield_parameters_config,
154 "crack_resistance", parameters, 1);
155 DBUG("Use '{}' as crack resistance.", crack_resistance.name);
156
157 // Crack length scale
158 auto const& crack_length_scale = ParameterLib::findParameter<double>(
159 phasefield_parameters_config,
161 "crack_length_scale", parameters, 1);
162 DBUG("Use '{}' as crack length scale.", crack_length_scale.name);
163
164 // Characteristic_length
165 auto const characteristic_length =
167 config.getConfigParameter<double>("characteristic_length", 1.0);
168
169 // Fluid compression modulus. Default value is 0.0.
170 auto const fluid_compressibility =
172 config.getConfigParameter<double>("fluid_compressibility", 0.0);
173
174 // Specific body force
175 Eigen::Matrix<double, DisplacementDim, 1> specific_body_force;
176 {
177 std::vector<double> const b =
179 config.getConfigParameter<std::vector<double>>(
180 "specific_body_force");
181 if (b.size() != DisplacementDim)
182 {
183 OGS_FATAL(
184 "The size of the specific body force vector does not match the "
185 "displacement dimension. Vector size is {:d}, displacement "
186 "dimension is {:d}",
187 b.size(), DisplacementDim);
188 }
189
190 std::copy_n(b.data(), b.size(), specific_body_force.data());
191 }
192
193 // Specific fracture direction
194 Eigen::Vector<double, DisplacementDim> specific_fracture_direction;
195 {
196 std::vector<double> const fracture_normal =
198 config.getConfigParameter<std::vector<double>>(
199 "specific_fracture_direction");
200 if (fracture_normal.size() != DisplacementDim)
201 {
202 OGS_FATAL(
203 "The size of the specific fracture direction vector does not "
204 "match the displacement dimension. Vector size is {:d}, "
205 "displacement dimension is {:d}",
206 fracture_normal.size(), DisplacementDim);
207 }
208
209 std::copy_n(fracture_normal.data(), fracture_normal.size(),
210 specific_fracture_direction.data());
211 }
212
213 auto const irreversible_threshold =
215 config.getConfigParameter<double>("irreversible_threshold", 0.05);
216
217 auto const phasefield_model_string =
219 config.getConfigParameter<std::string>("phasefield_model");
220 auto const phasefield_model =
222 DisplacementDim>(phasefield_model_string);
223
224 auto const softening_curve_string =
226 config.getConfigParameterOptional<std::string>("softening_curve");
227 auto const softening_curve =
229 DisplacementDim>(softening_curve_string);
230
231 auto const energy_split_model_string =
233 config.getConfigParameter<std::string>("energy_split_model");
234 auto const energy_split_model =
236 DisplacementDim>(energy_split_model_string);
237
238 auto degradation_derivative = creatDegradationDerivative<DisplacementDim>(
239 phasefield_model, characteristic_length, softening_curve);
240
241 auto const diffused_range_parameter =
243 config.getConfigParameter<double>("diffused_range_parameter", 2.0);
244
245 auto const fracture_threshold =
247 config.getConfigParameter<double>("fracture_threshold", 1.0);
248
249 auto const fracture_permeability_parameter =
251 config.getConfigParameter<double>("fracture_permeability_parameter",
252 0.0);
253
254 // Default value is 0.5, which is recommended by [Mikelic & Wheeler].
255 double const fixed_stress_stabilization_parameter =
257 config.getConfigParameter<double>(
258 "fixed_stress_stabilization_parameter", 0.5);
259
260 DBUG("Using value {:g} for coupling parameter of staggered scheme.",
261 fixed_stress_stabilization_parameter);
262
263 auto spatial_stabilization_parameter_optional =
265 config.getConfigParameterOptional<double>(
266 "spatial_stabilization_parameter");
267 double const spatial_stabilization_parameter =
268 spatial_stabilization_parameter_optional
269 ? spatial_stabilization_parameter_optional.value()
270 : 0.0;
271 DBUG("Using value {} for spatial stabilization coupling parameter.",
272 spatial_stabilization_parameter);
273
274 // Initial width
275 auto const& width_init = ParameterLib::findParameter<double>(
276 config,
278 "width_init", parameters, 0.0);
279 DBUG("Use '{}' as initial width.", width_init.name);
280
281 auto media_map =
283 DBUG("Check the media properties of HMPhaseField process ...");
284 checkMPLProperties(media);
285 DBUG("Media properties verified.");
286
288 materialIDs(mesh),
289 std::move(media_map),
290 std::move(solid_constitutive_relations),
291 residual_stiffness,
292 crack_resistance,
293 crack_length_scale,
294 specific_body_force,
295 specific_fracture_direction,
296 irreversible_threshold,
297 phasefield_model,
298 energy_split_model,
299 softening_curve,
300 characteristic_length,
301 std::move(degradation_derivative),
302 diffused_range_parameter,
303 fluid_compressibility,
304 fracture_threshold,
305 fracture_permeability_parameter,
306 fixed_stress_stabilization_parameter,
307 spatial_stabilization_parameter,
308 width_init};
309
310 SecondaryVariableCollection secondary_variables;
311
312 ProcessLib::createSecondaryVariables(config, secondary_variables);
313
314 return std::make_unique<HMPhaseFieldProcess<DisplacementDim>>(
315 std::move(name), mesh, std::move(jacobian_assembler), parameters,
316 integration_order, std::move(process_variables),
317 std::move(process_data), std::move(secondary_variables),
318 use_monolithic_scheme);
319}
#define OGS_FATAL(...)
Definition Error.h:19
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
std::string const & getName() const
int getNumberOfGlobalComponents() const
Returns the number of components of the process variable.
PhaseFieldModel convertStringToPhaseFieldModel(std::string const &phasefield_model)
SofteningCurve convertStringToSofteningCurve(std::optional< std::string > const &softening_curve)
EnergySplitModel convertStringToEnergySplitModel(std::string const &energy_split_model)
std::map< int, std::shared_ptr< MaterialLib::Solids::MechanicsBase< DisplacementDim > > > createConstitutiveRelations(std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, MeshLib::PropertyVector< int > const *const material_ids, BaseLib::ConfigTree const &config)
MaterialSpatialDistributionMap createMaterialSpatialDistributionMap(std::map< int, std::shared_ptr< Medium > > const &media, MeshLib::Mesh const &mesh)
PropertyVector< int > const * materialIDs(Mesh const &mesh)
Definition Mesh.cpp:258
OGS_NO_DANGLING Parameter< ParameterDataType > & findParameter(std::string const &parameter_name, std::vector< std::unique_ptr< ParameterBase > > const &parameters, int const num_components, MeshLib::Mesh const *const mesh=nullptr)
void checkMPLProperties(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
std::vector< std::reference_wrapper< ProcessVariable > > findProcessVariables(std::vector< ProcessVariable > const &variables, BaseLib::ConfigTree const &pv_config, std::initializer_list< std::string > tags)
void createSecondaryVariables(BaseLib::ConfigTree const &config, SecondaryVariableCollection &secondary_variables)

References BaseLib::ConfigTree::checkConfigParameter(), checkMPLProperties(), MaterialLib::Solids::Phasefield::convertStringToEnergySplitModel(), MaterialLib::Solids::Phasefield::convertStringToPhaseFieldModel(), MaterialLib::Solids::Phasefield::convertStringToSofteningCurve(), MaterialLib::Solids::createConstitutiveRelations(), MaterialPropertyLib::createMaterialSpatialDistributionMap(), ProcessLib::createSecondaryVariables(), DBUG(), ParameterLib::findParameter(), ProcessLib::findProcessVariables(), BaseLib::ConfigTree::getConfigParameter(), BaseLib::ConfigTree::getConfigParameterOptional(), BaseLib::ConfigTree::getConfigSubtree(), ProcessLib::ProcessVariable::getName(), ProcessLib::ProcessVariable::getNumberOfGlobalComponents(), and OGS_FATAL.

◆ createHMPhaseFieldProcess< 2 >()

template std::unique_ptr< Process > ProcessLib::HMPhaseField::createHMPhaseFieldProcess< 2 > ( std::string name,
MeshLib::Mesh & mesh,
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && jacobian_assembler,
std::vector< ProcessVariable > const & variables,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
std::optional< ParameterLib::CoordinateSystem > const & local_coordinate_system,
unsigned const integration_order,
BaseLib::ConfigTree const & config,
std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media )

◆ createHMPhaseFieldProcess< 3 >()

template std::unique_ptr< Process > ProcessLib::HMPhaseField::createHMPhaseFieldProcess< 3 > ( std::string name,
MeshLib::Mesh & mesh,
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && jacobian_assembler,
std::vector< ProcessVariable > const & variables,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
std::optional< ParameterLib::CoordinateSystem > const & local_coordinate_system,
unsigned const integration_order,
BaseLib::ConfigTree const & config,
std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media )

◆ showEnergyAndWork()

void ProcessLib::HMPhaseField::showEnergyAndWork ( const double t,
double & _elastic_energy,
double & _surface_energy,
double & _pressure_work )
inline

Definition at line 33 of file HMPhaseFieldProcessData.h.

35{
36#ifdef USE_PETSC
37 double const elastic_energy = _elastic_energy;
38 MPI_Allreduce(&elastic_energy, &_elastic_energy, 1, MPI_DOUBLE, MPI_SUM,
39 PETSC_COMM_WORLD);
40 double const surface_energy = _surface_energy;
41 MPI_Allreduce(&surface_energy, &_surface_energy, 1, MPI_DOUBLE, MPI_SUM,
42 PETSC_COMM_WORLD);
43 double const pressure_work = _pressure_work;
44 MPI_Allreduce(&pressure_work, &_pressure_work, 1, MPI_DOUBLE, MPI_SUM,
45 PETSC_COMM_WORLD);
46#endif
47
48 INFO("Elastic energy: {} Surface energy: {} Pressure work: {} at time: {} ",
49 _elastic_energy, _surface_energy, _pressure_work, t);
50};
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28

References INFO().

Referenced by ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::postTimestepConcreteProcess().