OGS
ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim > Class Template Referencefinal

Detailed Description

template<int GlobalDim>
class ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >

Definition at line 26 of file HydroMechanicsProcess.h.

#include <HydroMechanicsProcess.h>

Inheritance diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >:
[legend]

Public Member Functions

 HydroMechanicsProcess (std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, HydroMechanicsProcessData< GlobalDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
 
void postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id) override
 
ODESystem interface
bool isLinear () const override
 
- Public Member Functions inherited from ProcessLib::Process
 Process (std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
 
void preTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id)
 Preprocessing before starting assembly for new timestep.
 
void postTimestep (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep.
 
void postNonLinearSolver (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id)
 
void preIteration (const unsigned iter, GlobalVector const &x) final
 
void computeSecondaryVariable (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
 compute secondary variables for the coupled equations or for output.
 
NumLib::IterationResult postIteration (GlobalVector const &x) final
 
void initialize (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
void setInitialConditions (std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
 
void updateDeactivatedSubdomains (double const time, const int process_id)
 
virtual bool isMonolithicSchemeUsed () const
 
virtual void extrapolateIntegrationPointValuesToNodes (const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &)
 
void preAssemble (const double t, double const dt, GlobalVector const &x) final
 
void assemble (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) final
 
void assembleWithJacobian (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) final
 
void preOutput (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
 
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * getKnownSolutions (double const t, GlobalVector const &x, int const process_id) const final
 
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable (const int) const
 
MeshLib::MeshgetMesh () const
 
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables () const
 
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables (const int process_id) const
 
std::vector< std::size_t > const & getActiveElementIDs () const
 
SecondaryVariableCollection const & getSecondaryVariables () const
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const & getIntegrationPointWriters () const
 
virtual Eigen::Vector3d getFlux (std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const
 
virtual void solveReactionEquation (std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
 
bool requiresNormalization () const override
 
- Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
virtual std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes)
 
virtual ~SubmeshAssemblySupport ()=default
 

Private Types

using LocalAssemblerInterface = HydroMechanicsLocalAssemblerInterface
 

Private Member Functions

void constructDofTable () override
 
void initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
 Process specific initialization called by initialize().
 
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 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
 
void preTimestepConcreteProcess (std::vector< GlobalVector * > const &x, double const t, double const dt, int const process_id) override
 

Private Attributes

HydroMechanicsProcessData< GlobalDim > _process_data
 
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
 
std::vector< MeshLib::Element * > _vec_matrix_elements
 
std::vector< MeshLib::Element * > _vec_fracture_elements
 
std::vector< MeshLib::Element * > _vec_fracture_matrix_elements
 
std::vector< MeshLib::Node * > _vec_fracture_nodes
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_fracture_nodes
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_matrix_nodes
 
std::vector< MeshLib::Node * > _mesh_nodes_p
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_nodes_p
 

Additional Inherited Members

- Public Attributes inherited from ProcessLib::Process
std::string const name
 
- Static Public Attributes inherited from ProcessLib::Process
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name = "constant_one"
 
- Protected Member Functions inherited from ProcessLib::Process
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables (int const number_of_processes) const
 
NumLib::ExtrapolatorgetExtrapolator () const
 
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
 
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
void constructMonolithicProcessDofTable ()
 
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_process_id)
 
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
 
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
 
- Protected Attributes inherited from ProcessLib::Process
MeshLib::Mesh_mesh
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map
 
SecondaryVariableCollection _secondary_variables
 
CellAverageData cell_average_data_
 
std::unique_ptr< ProcessLib::AbstractJacobianAssembler_jacobian_assembler
 
VectorMatrixAssembler _global_assembler
 
const bool _use_monolithic_scheme
 
unsigned const _integration_order
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
 
GlobalSparsityPattern _sparsity_pattern
 
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
 
std::vector< BoundaryConditionCollection_boundary_conditions
 

Member Typedef Documentation

◆ LocalAssemblerInterface

template<int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::LocalAssemblerInterface = HydroMechanicsLocalAssemblerInterface
private

Definition at line 58 of file HydroMechanicsProcess.h.

Constructor & Destructor Documentation

◆ HydroMechanicsProcess()

template<int GlobalDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::HydroMechanicsProcess ( std::string name,
MeshLib::Mesh & mesh,
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && jacobian_assembler,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
unsigned const integration_order,
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > && process_variables,
HydroMechanicsProcessData< GlobalDim > && process_data,
SecondaryVariableCollection && secondary_variables,
bool const use_monolithic_scheme )

Definition at line 38 of file HydroMechanicsProcess.cpp.

49 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
50 integration_order, std::move(process_variables),
51 std::move(secondary_variables), use_monolithic_scheme),
52 _process_data(std::move(process_data))
53{
54 INFO("[LIE/HM] looking for fracture elements in the given mesh");
55 std::vector<int> vec_fracture_mat_IDs;
56 std::vector<std::vector<MeshLib::Element*>> vec_vec_fracture_elements;
57 std::vector<std::vector<MeshLib::Element*>>
58 vec_vec_fracture_matrix_elements;
59 std::vector<std::vector<MeshLib::Node*>> vec_vec_fracture_nodes;
60 std::vector<std::pair<std::size_t, std::vector<int>>>
61 vec_branch_nodeID_matIDs;
62 std::vector<std::pair<std::size_t, std::vector<int>>>
63 vec_junction_nodeID_matIDs;
65 mesh, _vec_matrix_elements, vec_fracture_mat_IDs,
66 vec_vec_fracture_elements, vec_vec_fracture_matrix_elements,
67 vec_vec_fracture_nodes, vec_branch_nodeID_matIDs,
68 vec_junction_nodeID_matIDs);
70 vec_vec_fracture_elements[0].begin(),
71 vec_vec_fracture_elements[0].end());
74 vec_vec_fracture_matrix_elements[0].begin(),
75 vec_vec_fracture_matrix_elements[0].end());
77 vec_vec_fracture_nodes[0].begin(),
78 vec_vec_fracture_nodes[0].end());
79
80 if (!_vec_fracture_elements.empty())
81 {
82 // set fracture property assuming a fracture forms a straight line
83 setFractureProperty(GlobalDim,
85 *_process_data.fracture_property);
86 }
87
88 //
89 // If Neumann BCs for the displacement_jump variable are required they need
90 // special treatment because of the levelset function. The implementation
91 // exists in the version 6.1.0 (e54815cc07ee89c81f953a4955b1c788595dd725)
92 // and was removed due to lack of applications.
93 //
94
95 if (!_process_data.deactivate_matrix_in_flow)
96 {
97 _process_data.p_element_status =
98 std::make_unique<MeshLib::ElementStatus>(&mesh);
99 }
100 else
101 {
102 auto const range =
104 if (!range)
105 {
106 OGS_FATAL(
107 "Could not get minimum/maximum ranges values for the "
108 "MaterialIDs property in the mesh '{:s}'.",
109 mesh.getName());
110 }
111
112 std::vector<int> vec_p_inactive_matIDs;
113 for (int matID = range->first; matID <= range->second; matID++)
114 {
115 if (std::find(vec_fracture_mat_IDs.begin(),
116 vec_fracture_mat_IDs.end(),
117 matID) == vec_fracture_mat_IDs.end())
118 {
119 vec_p_inactive_matIDs.push_back(matID);
120 }
121 }
122 _process_data.p_element_status =
123 std::make_unique<MeshLib::ElementStatus>(&mesh,
124 vec_p_inactive_matIDs);
125
126 const int monolithic_process_id = 0;
127 ProcessVariable const& pv_p =
128 getProcessVariables(monolithic_process_id)[0];
129 _process_data.p0 = &pv_p.getInitialCondition();
130 }
131}
#define OGS_FATAL(...)
Definition Error.h:26
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
const std::string getName() const
Get name of the mesh.
Definition Mesh.h:103
static std::optional< std::pair< T, T > > const getValueBounds(MeshLib::PropertyVector< T > const &property)
std::vector< MeshLib::Element * > _vec_fracture_matrix_elements
std::string const name
Definition Process.h:362
Process(std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
Definition Process.cpp:44
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables() const
Definition Process.h:156
PropertyVector< int > const * materialIDs(Mesh const &mesh)
Definition Mesh.cpp:268
void setFractureProperty(int const dim, MeshLib::Element const &e, FractureProperty &frac_prop)
void getFractureMatrixDataInMesh(MeshLib::Mesh const &mesh, std::vector< MeshLib::Element * > &vec_matrix_elements, std::vector< int > &vec_fracture_mat_IDs, std::vector< std::vector< MeshLib::Element * > > &vec_fracture_elements, std::vector< std::vector< MeshLib::Element * > > &vec_fracture_matrix_elements, std::vector< std::vector< MeshLib::Node * > > &vec_fracture_nodes, std::vector< std::pair< std::size_t, std::vector< int > > > &vec_branch_nodeID_matIDs, std::vector< std::pair< std::size_t, std::vector< int > > > &vec_junction_nodeID_matIDs)

References ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_process_data, ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_vec_fracture_elements, ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_vec_fracture_matrix_elements, ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_vec_fracture_nodes, ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_vec_matrix_elements, ProcessLib::LIE::getFractureMatrixDataInMesh(), ProcessLib::ProcessVariable::getInitialCondition(), MeshLib::Mesh::getName(), ProcessLib::Process::getProcessVariables(), MeshToolsLib::MeshInformation::getValueBounds(), INFO(), OGS_FATAL, and ProcessLib::LIE::setFractureProperty().

Member Function Documentation

◆ assembleConcreteProcess()

template<int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::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 )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 445 of file HydroMechanicsProcess.cpp.

449{
450 DBUG("Assemble HydroMechanicsProcess.");
451
452 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
454 // Call global assembler for each local assembly item.
457 dof_table, t, dt, x, x_prev, process_id, &M, &K, &b);
458}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
VectorMatrixAssembler _global_assembler
Definition Process.h:377
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:368
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)
static void executeMemberDereferenced(Object &object, Method method, Container const &container, Args &&... args)

References ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), and NumLib::SerialExecutor::executeMemberDereferenced().

◆ assembleWithJacobianConcreteProcess()

template<int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::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 )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 461 of file HydroMechanicsProcess.cpp.

465{
466 DBUG("AssembleWithJacobian HydroMechanicsProcess.");
467
468 // Call global assembler for each local assembly item.
469 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
473 _local_assemblers, getActiveElementIDs(), dof_table, t, dt, x, x_prev,
474 process_id, &b, &Jac);
475
476 auto copyRhs = [&](int const variable_id, auto& output_vector)
477 {
480 output_vector, std::negate<double>());
481 };
482 copyRhs(0, *_process_data.mesh_prop_hydraulic_flow);
483 copyRhs(1, *_process_data.mesh_prop_nodal_forces);
484 copyRhs(2, *_process_data.mesh_prop_nodal_forces_jump);
485}
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:167
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)
void transformVariableFromGlobalVector(GlobalVector const &input_vector, int const variable_id, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, MeshLib::PropertyVector< double > &output_vector, Functor map_function)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), and NumLib::SerialExecutor::executeSelectedMemberDereferenced().

◆ constructDofTable()

template<int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::constructDofTable ( )
overrideprivatevirtual

This function is for general cases, in which all equations of the coupled processes have the same number of unknowns. For the general cases with the staggered scheme, all equations of the coupled processes share one DOF table hold by _local_to_global_index_map. Other cases can be considered by overloading this member function in the derived class.

Reimplemented from ProcessLib::Process.

Definition at line 134 of file HydroMechanicsProcess.cpp.

135{
136 //------------------------------------------------------------
137 // prepare mesh subsets to define DoFs
138 //------------------------------------------------------------
139 // for extrapolation
141 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
142 // pressure
144 _process_data.p_element_status->getActiveElements());
146 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh_nodes_p);
147 // regular u
149 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
150 if (!_vec_fracture_nodes.empty())
151 {
152 // u jump
154 std::make_unique<MeshLib::MeshSubset>(_mesh, _vec_fracture_nodes);
155 }
156
157 // Collect the mesh subsets in a vector.
158 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
159 std::vector<int> vec_n_components;
160 std::vector<std::vector<MeshLib::Element*> const*> vec_var_elements;
161 // pressure
162 vec_n_components.push_back(1);
163 all_mesh_subsets.emplace_back(*_mesh_subset_nodes_p);
164 if (!_process_data.deactivate_matrix_in_flow)
165 {
166 vec_var_elements.push_back(&_mesh.getElements());
167 }
168 else
169 {
170 // TODO set elements including active nodes for pressure.
171 // cannot use ElementStatus
172 vec_var_elements.push_back(&_vec_fracture_matrix_elements);
173 }
174 // regular displacement
175 vec_n_components.push_back(GlobalDim);
176 std::generate_n(std::back_inserter(all_mesh_subsets), GlobalDim,
177 [&]() { return *_mesh_subset_matrix_nodes; });
178 vec_var_elements.push_back(&_vec_matrix_elements);
179 if (!_vec_fracture_nodes.empty())
180 {
181 // displacement jump
182 vec_n_components.push_back(GlobalDim);
183 std::generate_n(std::back_inserter(all_mesh_subsets), GlobalDim,
184 [&]() { return *_mesh_subset_fracture_nodes; });
185 vec_var_elements.push_back(&_vec_fracture_matrix_elements);
186 }
187
188 INFO("[LIE/HM] creating a DoF table");
190 std::make_unique<NumLib::LocalToGlobalIndexMap>(
191 std::move(all_mesh_subsets),
192 vec_n_components,
193 vec_var_elements,
195
196 DBUG("created {:d} DoF", _local_to_global_index_map->size());
197}
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_matrix_nodes
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_fracture_nodes
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_nodes_p
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:366
MeshLib::Mesh & _mesh
Definition Process.h:365
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
Definition Utils.h:26
@ BY_COMPONENT
Ordering data by component type.

References NumLib::BY_COMPONENT, DBUG(), MeshLib::getBaseNodes(), and INFO().

◆ initializeConcreteProcess()

template<int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::initializeConcreteProcess ( NumLib::LocalToGlobalIndexMap const & dof_table,
MeshLib::Mesh const & mesh,
unsigned const integration_order )
overrideprivatevirtual

Process specific initialization called by initialize().

Implements ProcessLib::Process.

Definition at line 200 of file HydroMechanicsProcess.cpp.

204{
205 assert(mesh.getDimension() == GlobalDim);
206 INFO("[LIE/HM] creating local assemblers");
208 GlobalDim, HydroMechanicsLocalAssemblerMatrix,
209 HydroMechanicsLocalAssemblerMatrixNearFracture,
210 HydroMechanicsLocalAssemblerFracture>(
211 mesh.getElements(), dof_table, _local_assemblers,
212 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
214
215 auto add_secondary_variable = [&](std::string const& name,
216 int const num_components,
217 auto get_ip_values_function)
218 {
220 name,
221 makeExtrapolator(num_components, getExtrapolator(),
223 std::move(get_ip_values_function)));
224 };
225
226 add_secondary_variable(
227 "sigma",
230
231 add_secondary_variable(
232 "epsilon",
235
236 add_secondary_variable("velocity", GlobalDim,
238
239 add_secondary_variable("fracture_velocity", GlobalDim,
241
242 add_secondary_variable("fracture_stress", GlobalDim,
244
245 add_secondary_variable("fracture_aperture", 1,
247
248 add_secondary_variable(
249 "fracture_permeability", 1,
251
253 const_cast<MeshLib::Mesh&>(mesh), "sigma_avg",
256
258 const_cast<MeshLib::Mesh&>(mesh), "velocity_avg",
259 MeshLib::MeshItemType::Cell, GlobalDim);
260
261 if (!_vec_fracture_elements.empty())
262 {
263 auto mesh_prop_levelset = MeshLib::getOrCreateMeshProperty<double>(
264 const_cast<MeshLib::Mesh&>(mesh), "levelset1",
266 mesh_prop_levelset->resize(mesh.getNumberOfElements());
267 for (MeshLib::Element const* e : _mesh.getElements())
268 {
269 if (e->getDimension() < GlobalDim)
270 {
271 continue;
272 }
273
274 std::vector<FractureProperty*> fracture_props(
275 {_process_data.fracture_property.get()});
276 std::vector<JunctionProperty*> junction_props;
277 std::unordered_map<int, int> fracID_to_local({{0, 0}});
278 std::vector<double> levelsets = uGlobalEnrichments(
279 fracture_props, junction_props, fracID_to_local,
280 Eigen::Vector3d(MeshLib::getCenterOfGravity(*e).data()));
281 (*mesh_prop_levelset)[e->getID()] = levelsets[0];
282 }
283
284 _process_data.element_local_jumps =
286 const_cast<MeshLib::Mesh&>(mesh), "local_jump_w_avg",
287 MeshLib::MeshItemType::Cell, GlobalDim);
288
289 _process_data.element_fracture_stresses =
291 const_cast<MeshLib::Mesh&>(mesh), "fracture_stress_avg",
292 MeshLib::MeshItemType::Cell, GlobalDim);
293
294 _process_data.element_fracture_velocities =
296 const_cast<MeshLib::Mesh&>(mesh), "fracture_velocity_avg",
297 MeshLib::MeshItemType::Cell, GlobalDim);
298
300 const_cast<MeshLib::Mesh&>(mesh), "fracture_aperture_avg",
302 mesh_prop_b->resize(mesh.getNumberOfElements());
303
304 auto const mesh_prop_matid = materialIDs(mesh);
305 if (!mesh_prop_matid)
306 {
307 OGS_FATAL("Could not access MaterialIDs property from mesh.");
308 }
309 auto const& frac = _process_data.fracture_property;
310 for (MeshLib::Element const* e : _mesh.getElements())
311 {
312 if (e->getDimension() == GlobalDim)
313 {
314 continue;
315 }
316 if ((*mesh_prop_matid)[e->getID()] != frac->mat_id)
317 {
318 continue;
319 }
320 // Mean value for the element. This allows usage of node based
321 // properties for aperture.
322 (*mesh_prop_b)[e->getID()] =
323 frac->aperture0
324 .getNodalValuesOnElement(*e, /*time independent*/ 0)
325 .mean();
326 }
327 _process_data.mesh_prop_b = mesh_prop_b;
328
329 auto mesh_prop_k_f = MeshLib::getOrCreateMeshProperty<double>(
330 const_cast<MeshLib::Mesh&>(mesh), "fracture_permeability_avg",
332 mesh_prop_k_f->resize(mesh.getNumberOfElements());
333 _process_data.mesh_prop_k_f = mesh_prop_k_f;
334
335 auto mesh_prop_fracture_shear_failure =
337 const_cast<MeshLib::Mesh&>(mesh), "f_shear_failure",
339 mesh_prop_fracture_shear_failure->resize(mesh.getNumberOfElements());
340 _process_data.mesh_prop_fracture_shear_failure =
341 mesh_prop_fracture_shear_failure;
342
343 auto mesh_prop_nodal_p = MeshLib::getOrCreateMeshProperty<double>(
344 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
346 mesh_prop_nodal_p->resize(mesh.getNumberOfNodes());
347 _process_data.mesh_prop_nodal_p = mesh_prop_nodal_p;
348
349 _process_data.mesh_prop_nodal_forces =
351 const_cast<MeshLib::Mesh&>(mesh), "NodalForces",
352 MeshLib::MeshItemType::Node, GlobalDim);
353 assert(_process_data.mesh_prop_nodal_forces->size() ==
354 GlobalDim * mesh.getNumberOfNodes());
355
356 _process_data.mesh_prop_nodal_forces_jump =
358 const_cast<MeshLib::Mesh&>(mesh), "NodalForcesJump",
359 MeshLib::MeshItemType::Node, GlobalDim);
360 assert(_process_data.mesh_prop_nodal_forces_jump->size() ==
361 GlobalDim * mesh.getNumberOfNodes());
362
363 _process_data.mesh_prop_hydraulic_flow =
365 const_cast<MeshLib::Mesh&>(mesh), "MassFlowRate",
367 assert(_process_data.mesh_prop_hydraulic_flow->size() ==
368 mesh.getNumberOfNodes());
369 }
370}
const double * data() const
Definition Point3d.h:60
virtual std::vector< double > const & getIntPtFracturePermeability(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtFractureAperture(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtFractureStress(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtFractureVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
SecondaryVariableCollection _secondary_variables
Definition Process.h:370
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
MathLib::Point3d getCenterOfGravity(Element const &element)
Calculates the center of gravity for the mesh element.
Definition Element.cpp:124
void createLocalAssemblers(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, NumLib::IntegrationOrder const integration_order, ExtraCtorArgs &&... extra_ctor_args)
std::vector< double > uGlobalEnrichments(std::vector< FractureProperty * > const &frac_props, std::vector< JunctionProperty * > const &junction_props, std::unordered_map< int, int > const &fracID_to_local, Eigen::Vector3d const &x)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)

References MeshLib::Cell, ProcessLib::LIE::HydroMechanics::createLocalAssemblers(), MathLib::Point3d::data(), MeshLib::getCenterOfGravity(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getNumberOfElements(), MeshLib::Mesh::getNumberOfNodes(), MeshLib::getOrCreateMeshProperty(), INFO(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), MeshLib::Node, OGS_FATAL, and ProcessLib::LIE::uGlobalEnrichments().

◆ isLinear()

template<int GlobalDim>
bool ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::isLinear ( ) const
override

Definition at line 439 of file HydroMechanicsProcess.cpp.

440{
441 return false;
442}

◆ postTimestepConcreteProcess()

template<int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::postTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
double const t,
double const dt,
int const process_id )
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 373 of file HydroMechanicsProcess.cpp.

377{
378 if (process_id == 0)
379 {
380 DBUG("PostTimestep HydroMechanicsProcess.");
381
385 x_prev, t, dt, process_id);
386 }
387
388 DBUG("Compute the secondary variables for HydroMechanicsProcess.");
389
390 const auto& dof_table = getDOFTable(process_id);
391
392 // Copy displacement jumps in a solution vector to mesh property
393 // Remark: the copy is required because mesh properties for primary
394 // variables are set during output and are not ready yet when this function
395 // is called.
396 int g_variable_id = 0;
397 {
398 const int monolithic_process_id = 0;
399 auto const& pvs = getProcessVariables(monolithic_process_id);
400 auto const it =
401 std::find_if(pvs.begin(), pvs.end(),
402 [](ProcessVariable const& pv)
403 { return pv.getName() == "displacement_jump1"; });
404 if (it == pvs.end())
405 {
406 OGS_FATAL(
407 "Didn't find expected 'displacement_jump1' process variable.");
408 }
409 g_variable_id = static_cast<int>(std::distance(pvs.begin(), it));
410 }
411
413
414 const int monolithic_process_id = 0;
415 ProcessVariable& pv_g =
416 this->getProcessVariables(monolithic_process_id)[g_variable_id];
417 auto const num_comp = pv_g.getNumberOfGlobalComponents();
418 auto& mesh_prop_g = *MeshLib::getOrCreateMeshProperty<double>(
419 _mesh, pv_g.getName(), MeshLib::MeshItemType::Node, num_comp);
420 for (int component_id = 0; component_id < num_comp; ++component_id)
421 {
422 auto const& mesh_subset =
423 dof_table.getMeshSubset(g_variable_id, component_id);
424 auto const mesh_id = mesh_subset.getMeshID();
425 for (auto const* node : mesh_subset.getNodes())
426 {
428 node->getID());
429
430 auto const global_index =
431 dof_table.getGlobalIndex(l, g_variable_id, component_id);
432 mesh_prop_g[node->getID() * num_comp + component_id] =
433 (*x[process_id])[global_index];
434 }
435 }
436}
std::vector< std::size_t > getNodes(GeoLib::Point const &pnt, std::vector< MeshLib::Node * > const &nodes, MeshLib::PropertyVector< int > const &mat_ids, std::pair< int, int > const &mat_limits, std::pair< double, double > const &elevation_limits, MeshLib::Mesh const &mesh)
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:385
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable(const int) const
Definition Process.h:147
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:27
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getName(), ProcessLib::ProcessVariable::getNumberOfGlobalComponents(), MeshLib::getOrCreateMeshProperty(), MeshLib::Node, OGS_FATAL, ProcessLib::LocalAssemblerInterface::postTimestep(), and MathLib::LinAlg::setLocalAccessibleVector().

◆ preTimestepConcreteProcess()

template<int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::preTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
double const t,
double const dt,
int const process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 488 of file HydroMechanicsProcess.cpp.

491{
492 DBUG("PreTimestep HydroMechanicsProcess.");
493
496 getActiveElementIDs(), *_local_to_global_index_map, *x[process_id], t,
497 dt);
498}
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)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), and ProcessLib::LocalAssemblerInterface::preTimestep().

Member Data Documentation

◆ _local_assemblers

template<int GlobalDim>
std::vector<std::unique_ptr<LocalAssemblerInterface> > ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_local_assemblers
private

Definition at line 84 of file HydroMechanicsProcess.h.

◆ _mesh_nodes_p

template<int GlobalDim>
std::vector<MeshLib::Node*> ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_mesh_nodes_p
private

Definition at line 94 of file HydroMechanicsProcess.h.

◆ _mesh_subset_fracture_nodes

template<int GlobalDim>
std::unique_ptr<MeshLib::MeshSubset const> ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_mesh_subset_fracture_nodes
private

Definition at line 91 of file HydroMechanicsProcess.h.

◆ _mesh_subset_matrix_nodes

template<int GlobalDim>
std::unique_ptr<MeshLib::MeshSubset const> ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_mesh_subset_matrix_nodes
private

Definition at line 92 of file HydroMechanicsProcess.h.

◆ _mesh_subset_nodes_p

template<int GlobalDim>
std::unique_ptr<MeshLib::MeshSubset const> ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_mesh_subset_nodes_p
private

Definition at line 95 of file HydroMechanicsProcess.h.

◆ _process_data

◆ _vec_fracture_elements

template<int GlobalDim>
std::vector<MeshLib::Element*> ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_vec_fracture_elements
private

◆ _vec_fracture_matrix_elements

template<int GlobalDim>
std::vector<MeshLib::Element*> ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_vec_fracture_matrix_elements
private

◆ _vec_fracture_nodes

template<int GlobalDim>
std::vector<MeshLib::Node*> ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_vec_fracture_nodes
private

◆ _vec_matrix_elements

template<int GlobalDim>
std::vector<MeshLib::Element*> ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::_vec_matrix_elements
private

The documentation for this class was generated from the following files: