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 30 of file CreateHMPhaseFieldProcess.cpp.

32{
33 std::array const requiredMediumProperties = {
36 std::array const requiredFluidProperties = {MaterialPropertyLib::viscosity,
38 std::array const requiredSolidProperties = {
40
41 for (auto const& m : media)
42 {
43 checkRequiredProperties(*m.second, requiredMediumProperties);
44 checkRequiredProperties(m.second->phase("AqueousLiquid"),
45 requiredFluidProperties);
46 checkRequiredProperties(m.second->phase("Solid"),
47 requiredSolidProperties);
48 }
49}
void checkRequiredProperties(Component const &c, std::span< PropertyType const > const required_properties)
Definition Component.cpp:60

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 52 of file CreateHMPhaseFieldProcess.cpp.

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

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};
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35

References INFO().

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