KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
apply_constant_phreatic_multi_line_pressure_process.h
Go to the documentation of this file.
1 // KRATOS___
2 // // ) )
3 // // ___ ___
4 // // ____ //___) ) // ) )
5 // // / / // // / /
6 // ((____/ / ((____ ((___/ / MECHANICS
7 //
8 // License: geo_mechanics_application/license.txt
9 //
10 // Main authors: Vahid Galavi
11 // Jonathan Nuttall
12 
13 #pragma once
14 
15 #include <algorithm>
16 #include "includes/kratos_flags.h"
18 #include "processes/process.h"
19 
21 
22 namespace Kratos
23 {
24 
25 class KRATOS_API(GEO_MECHANICS_APPLICATION) ApplyConstantPhreaticMultiLinePressureProcess : public Process
26 {
27 
28 public:
30 
35 
38  void ExecuteInitialize() override;
39 
40  std::string Info() const override;
41  void PrintInfo(std::ostream& rOStream) const override;
42  const std::string& VariableName() const;
43  bool IsFixed() const;
44  bool IsSeepage() const;
45  unsigned int GravityDirection() const;
46  double SpecificWeight() const;
47  unsigned int OutOfPlaneDirection() const;
48  unsigned int HorizontalDirection() const;
49  const Vector& HorizontalDirectionCoordinates() const;
50  const Vector& GravityDirectionCoordinates() const;
51  const Vector& XCoordinates() const;
52  const Vector& YCoordinates() const;
53  const Vector& ZCoordinates() const;
54  double PressureTensionCutOff() const;
55  bool IsFixedProvided() const;
56 
57 protected:
59  int findIndex(const Node& rNode) const;
60  double CalculatePressure(const Node& rNode, std::vector<double> deltaH = {}) const;
61 
62 private:
63  std::string mVariableName;
64  bool mIsFixed;
65  bool mIsFixedProvided;
66  bool mIsSeepage;
67  unsigned int mGravityDirection;
68  double mSpecificWeight;
69  unsigned int mOutOfPlaneDirection;
70  unsigned int mHorizontalDirection;
71  Vector mHorizontalDirectionCoordinates;
72  Vector mGravityDirectionCoordinates;
73  Vector mXCoordinates;
74  Vector mYCoordinates;
75  Vector mZCoordinates;
76  double mPressureTensionCutOff;
77 
78  void InitializeHorizontalDirection();
79  void InitializeGravityDirection(const Parameters &rParameters);
80  void InitializeCoordinates(const Parameters &rParameters);
81  void ValidateCoordinates(const Parameters &rParameters) const;
82 
83  void InitializeParameters(Parameters &rParameters) const;
84 }; // Class ApplyConstantPhreaticMultiLinePressureProcess
85 
87 inline std::istream& operator >> (std::istream& rIStream,
89 
91 inline std::ostream& operator << (std::ostream& rOStream,
93 {
94  rThis.PrintInfo(rOStream);
95  rOStream << std::endl;
96  rThis.PrintData(rOStream);
97 
98  return rOStream;
99 }
100 
101 } // namespace Kratos.
void ExecuteInitialize() override
Definition: periodic_interface_process.hpp:37
std::string Info() const override
Turn back information as a string.
Definition: periodic_interface_process.hpp:93
Definition: apply_constant_phreatic_multi_line_pressure_process.h:26
ModelPart & mrModelPart
Definition: apply_constant_phreatic_multi_line_pressure_process.h:58
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: apply_constant_phreatic_multi_line_pressure_process.cpp:168
ApplyConstantPhreaticMultiLinePressureProcess & operator=(ApplyConstantPhreaticMultiLinePressureProcess const &)=delete
ApplyConstantPhreaticMultiLinePressureProcess(ApplyConstantPhreaticMultiLinePressureProcess const &)=delete
KRATOS_CLASS_POINTER_DEFINITION(ApplyConstantPhreaticMultiLinePressureProcess)
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
This class defines the node.
Definition: node.h:65
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
The base class for all processes in Kratos.
Definition: process.h:49
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: process.h:210
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
model_part
Definition: face_heat.py:14