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.
newtonian_plane_strain_2D_law.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosConstitutiveModelsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: April 2017 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_NEWTONIAN_FLUID_PLANE_STRAIN_2D_LAW_H_INCLUDED )
11 #define KRATOS_NEWTONIAN_FLUID_PLANE_STRAIN_2D_LAW_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
19 
20 namespace Kratos
21 {
24 
27 
31 
35 
39 
43 
45 
48  {
49  public:
52 
55 
59 
62 
63 
66 
69  {
71  return *this;
72  }
73 
75  ConstitutiveLaw::Pointer Clone() const override
76  {
77  return Kratos::make_shared<NewtonianFluidPlaneStrain2DLaw>(*this);
78  }
79 
82 
83 
87 
89  SizeType WorkingSpaceDimension() override { return 2; }
90 
92  SizeType GetStrainSize() const override { return 3; }
93 
95  void GetLawFeatures(Features& rFeatures) override
96  {
98 
99  //Set the type of law
100  rFeatures.mOptions.Set( PLANE_STRAIN_LAW );
101  rFeatures.mOptions.Set( INFINITESIMAL_STRAINS );
102  rFeatures.mOptions.Set( ISOTROPIC );
103 
104  //Set strain measure required by the consitutive law
106 
107  //Set the strain size
108  rFeatures.mStrainSize = GetStrainSize();
109 
110  //Set the spacedimension
112 
113  KRATOS_CATCH(" ")
114  }
115 
119 
120 
124 
125 
129 
130 
134 
136  std::string Info() const override
137  {
138  std::stringstream buffer;
139  buffer << "NewtonianFluidPlaneStrain2DLaw" ;
140  return buffer.str();
141  }
142 
144  void PrintInfo(std::ostream& rOStream) const override {rOStream << "NewtonianFluidPlaneStrain2DLaw";}
145 
147  void PrintData(std::ostream& rOStream) const override {}
148 
149 
153 
154 
156 
157  protected:
160 
161 
165 
166 
170 
171 
178  void CalculateStress(Vector& rStressVector,
179  const Vector & rStrainVector,
180  const Properties& rProperties) override
181  {
182  KRATOS_TRY
183 
184  const double& rViscosity = rProperties[DYNAMIC_VISCOSITY];
185 
186  const double pressure = (rStrainVector[0]+rStrainVector[1])/3.0;
187 
188  // Cauchy StressVector
189  rStressVector[0] = 2.0*rViscosity*(rStrainVector[0] - pressure);
190  rStressVector[1] = 2.0*rViscosity*(rStrainVector[1] - pressure);
191  rStressVector[2] = rViscosity*rStrainVector[2];
192 
193  KRATOS_CATCH(" ")
194 
195  }
196 
197 
203  void CalculateConstitutiveMatrix(Matrix& rConstitutiveMatrix,
204  const Properties& rProperties) override
205  {
206  KRATOS_TRY
207 
208  // Viscosity
209  const double& rViscosity = rProperties[DYNAMIC_VISCOSITY];
210 
211  const double diagonal_component = 4.0 * rViscosity / 3.0;
212  const double side_component = -0.5 * diagonal_component;
213 
214  // 3D linear elastic constitutive matrix
215  rConstitutiveMatrix ( 0 , 0 ) = diagonal_component;
216  rConstitutiveMatrix ( 1 , 1 ) = diagonal_component;
217 
218  rConstitutiveMatrix ( 2 , 2 ) = rViscosity;
219 
220  rConstitutiveMatrix ( 0 , 1 ) = side_component;
221  rConstitutiveMatrix ( 1 , 0 ) = side_component;
222 
223  //initialize to zero other values
224  rConstitutiveMatrix ( 0 , 2 ) = 0.0;
225  rConstitutiveMatrix ( 2 , 0 ) = 0.0;
226  rConstitutiveMatrix ( 1 , 2 ) = 0.0;
227  rConstitutiveMatrix ( 2 , 1 ) = 0.0;
228 
229  KRATOS_CATCH(" ")
230  }
231 
232 
236 
237 
241 
242 
246 
247 
251 
252 
254 
255  private:
258 
259 
263 
264 
268 
269 
273 
274 
278 
279 
283 
284 
288  friend class Serializer;
289 
290  void save(Serializer& rSerializer) const override
291  {
293  }
294 
295  void load(Serializer& rSerializer) override
296  {
298  }
299 
300 
304 
306 
307  }; // Class NewtonianFluidPlaneStrain2DLaw
308 
310 
313 
314 
318 
319 
321 
323 
324 } // namespace Kratos.
325 
326 #endif // KRATOS_NEWTONIAN_FLUID_PLANE_STRAIN_2D_LAW_H_INCLUDED defined
@ StrainMeasure_Velocity_Gradient
Definition: constitutive_law.h:65
std::size_t SizeType
Definition: constitutive_law.h:82
void Set(const Flags ThisFlag)
Definition: flags.cpp:33
Definition: newtonian_3D_law.hpp:30
NewtonianFluid3DLaw & operator=(const NewtonianFluid3DLaw &rOther)
Assignment operator.
Definition: newtonian_3D_law.cpp:43
Short class definition.
Definition: newtonian_plane_strain_2D_law.hpp:48
SizeType GetStrainSize() const override
Law Voigt Strain Size.
Definition: newtonian_plane_strain_2D_law.hpp:92
void GetLawFeatures(Features &rFeatures) override
Law Features.
Definition: newtonian_plane_strain_2D_law.hpp:95
void CalculateConstitutiveMatrix(Matrix &rConstitutiveMatrix, const Properties &rProperties) override
Definition: newtonian_plane_strain_2D_law.hpp:203
NewtonianFluidPlaneStrain2DLaw(const NewtonianFluidPlaneStrain2DLaw &rOther)
Copy constructor.
Definition: newtonian_plane_strain_2D_law.hpp:65
void CalculateStress(Vector &rStressVector, const Vector &rStrainVector, const Properties &rProperties) override
Definition: newtonian_plane_strain_2D_law.hpp:178
NewtonianFluidPlaneStrain2DLaw()
Default constructor.
Definition: newtonian_plane_strain_2D_law.hpp:61
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: newtonian_plane_strain_2D_law.hpp:144
std::string Info() const override
Turn back information as a string.
Definition: newtonian_plane_strain_2D_law.hpp:136
ConstitutiveLaw::Pointer Clone() const override
Clone.
Definition: newtonian_plane_strain_2D_law.hpp:75
NewtonianFluidPlaneStrain2DLaw & operator=(NewtonianFluidPlaneStrain2DLaw const &rOther)
Assignment operator.
Definition: newtonian_plane_strain_2D_law.hpp:68
~NewtonianFluidPlaneStrain2DLaw() override
Destructor.
Definition: newtonian_plane_strain_2D_law.hpp:81
KRATOS_CLASS_POINTER_DEFINITION(NewtonianFluidPlaneStrain2DLaw)
Pointer definition of NewtonianFluidPlaneStrain2DLawg.
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: newtonian_plane_strain_2D_law.hpp:147
SizeType WorkingSpaceDimension() override
Law Dimension.
Definition: newtonian_plane_strain_2D_law.hpp:89
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
def load(f)
Definition: ode_solve.py:307
float pressure
Definition: test_pureconvectionsolver_benchmarking.py:101
Definition: constitutive_law.h:137
SizeType mStrainSize
Definition: constitutive_law.h:152
std::vector< StrainMeasure > mStrainMeasures
Definition: constitutive_law.h:154
SizeType mSpaceDimension
Definition: constitutive_law.h:153
Flags mOptions
Definition: constitutive_law.h:151