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.
updated_lagrangian_fluid_inc.h
Go to the documentation of this file.
1 /*
2 ==============================================================================
3 KratosULFApplication
4 A library based on:
5 Kratos
6 A General Purpose Software for Multi-Physics Finite Element Analysis
7 Version 1.0 (Released on march 05, 2007).
8 
9 Copyright 2007
10 Pooyan Dadvand, Riccardo Rossi, Pawel Ryzhakov
11 pooyan@cimne.upc.edu
12 rrossi@cimne.upc.edu
13 - CIMNE (International Center for Numerical Methods in Engineering),
14 Gran Capita' s/n, 08034 Barcelona, Spain
15 
16 
17 Permission is hereby granted, free of charge, to any person obtaining
18 a copy of this software and associated documentation files (the
19 "Software"), to deal in the Software without restriction, including
20 without limitation the rights to use, copy, modify, merge, publish,
21 distribute, sublicense and/or sell copies of the Software, and to
22 permit persons to whom the Software is furnished to do so, subject to
23 the following condition:
24 
25 Distribution of this code for any commercial purpose is permissible
26 ONLY BY DIRECT ARRANGEMENT WITH THE COPYRIGHT OWNERS.
27 
28 The above copyright notice and this permission notice shall be
29 included in all copies or substantial portions of the Software.
30 
31 THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
32 EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
33 MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.
34 IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY
35 CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT,
36 TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE
37 SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
38 
39 ==============================================================================
40 */
41 
42 
43 //
44 // Project Name: Kratos
45 // Last Modified by: $Author: anonymous $
46 // Date: $Date: 2007-11-06 12:28:50 $
47 // Revision: $Revision: 1.1 $
48 //
49 //
50 
51 
52 #if !defined(KRATOS_UPDATED_LAGRANGIAN__FLUID_INC_H_INCLUDED )
53 #define KRATOS_UPDATED_LAGRANGIAN__FLUID_INC_H_INCLUDED
54 
55 
56 
57 // System includes
58 
59 
60 // External includes
61 #include "boost/smart_ptr.hpp"
62 
63 
64 // Project includes
65 #include "includes/define.h"
66 #include "includes/element.h"
68 #include "includes/variables.h"
69 #include "includes/serializer.h"
70 
71 namespace Kratos
72 {
73 
76 
80 
84 
88 
92 
94 
97  : public Element
98 {
99 public:
102 
105 
109 
111  UpdatedLagrangianFluidInc(IndexType NewId, GeometryType::Pointer pGeometry);
112  UpdatedLagrangianFluidInc(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
113 
115  ~UpdatedLagrangianFluidInc() override;
116 
117 
121 
122 
126 
127  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const;
128 
129  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, ProcessInfo& rCurrentProcessInfo);
130 
131  void CalculateRightHandSide(VectorType& rRightHandSideVector, ProcessInfo& rCurrentProcessInfo);
132  //virtual void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix, ProcessInfo& rCurrentProcessInfo);
133 
134  void EquationIdVector(EquationIdVectorType& rResult, ProcessInfo& rCurrentProcessInfo);
135 
136  void GetDofList(DofsVectorType& ElementalDofList,ProcessInfo& CurrentProcessInfo);
137 
138  void InitializeSolutionStep(ProcessInfo& CurrentProcessInfo);
139 
140  void Calculate(const Variable<double >& rVariable, double& Output, const ProcessInfo& rCurrentProcessInfo);
141  void CalculateMassMatrix(MatrixType& rMassMatrix, ProcessInfo& rCurrentProcessInfo);
142  void CalculateDampingMatrix(MatrixType& rDampingMatrix, ProcessInfo& rCurrentProcessInfo);
143  void GetValuesVector(Vector& values, int Step = 0);
144  void GetFirstDerivativesVector(Vector& values, int Step = 0);
145  void GetSecondDerivativesVector(Vector& values, int Step = 0);
146 
147 
151 
152 
156 
157 
161 
163  std::string Info() const override
164  {
165  return "UpdatedLagrangianFluidInc #" ;
166  }
167 
169  void PrintInfo(std::ostream& rOStream) const override
170  {
171  rOStream << Info() << Id();
172  }
173 
175 // virtual void PrintData(std::ostream& rOStream) const;
176 
177 
181 
182 
184 
185 protected:
188 
189 
193 
194 
198 
199 
203 
204 
208 
209 
213 
214 
218 
219 
221 
222 private:
225  static BoundedMatrix<double,3,3> msMassFactors;
226  static BoundedMatrix<double,3,2> msDN_Dx;
227  static array_1d<double,3> msN; //dimension = number of nodes
228  //static Matrix msDN_DX;
229  //static Matrix msMassFactors;
230  static array_1d<double,2> ms_vel_gauss; //dimesion coincides with space dimension
231  static array_1d<double,3> ms_temp_vec_np; //dimension = number of nodes
232  static array_1d<double,3> ms_u_DN;
233 
234 
235  static BoundedMatrix<double,3,6> msB;
236  static BoundedMatrix<double,3,3> ms_constitutive_matrix;
237  static BoundedMatrix<double,3,6> ms_temp;
238  static array_1d<double,6> ms_temp_vec;
239 
240 
244  double mA0; //stores area at the beginning of time-step
248  friend class Serializer;
249 
250  // A private default constructor necessary for serialization
252  {
253  }
254 
255  void save(Serializer& rSerializer) const override
256  {
257 
259  }
260 
261  void load(Serializer& rSerializer) override
262  {
264  }
265 
266 
270  void Stage1(MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, ProcessInfo& rCurrentProcessInfo, unsigned int ComponentIndex);
271  void Stage2(MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, ProcessInfo& rCurrentProcessInfo);
272 
273  void MeshMovingStep(MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, ProcessInfo& rCurrentProcessInfo, unsigned int ComponentIndex);
274 
275  //inline void CalculateGeometryData(Matrix& msDN_DX, Vector& N, double& Area)
276  inline void CalculateGeometryData(BoundedMatrix<double,3,2>& DN_DX, array_1d<double,3>& N, double& Area);
277 
281 
282 
286 
287 
291 
292 
296 
298  //UpdatedLagrangianFluidInc& operator=(const UpdatedLagrangianFluidInc& rOther);
299 
301  //UpdatedLagrangianFluidInc(const UpdatedLagrangianFluidInc& rOther);
302 
303 
305 
306 }; // Class UpdatedLagrangianFluidInc
307 
309 
312 
313 
317 
318 
320 /* inline std::istream& operator >> (std::istream& rIStream,
321  UpdatedLagrangianFluidInc& rThis);
322 */
324 /* inline std::ostream& operator << (std::ostream& rOStream,
325  const UpdatedLagrangianFluidInc& rThis)
326  {
327  rThis.PrintInfo(rOStream);
328  rOStream << std::endl;
329  rThis.PrintData(rOStream);
330 
331  return rOStream;
332  }*/
334 
335 } // namespace Kratos.
336 
337 #endif // KRATOS_UPDATED_LAGRANGIAN__FLUID_INC_H_INCLUDED defined
338 
339 
Base class for all Elements.
Definition: element.h:60
Vector VectorType
Definition: element.h:88
Element(IndexType NewId=0)
Definition: element.h:121
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
Matrix MatrixType
Definition: element.h:90
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
std::size_t IndexType
Definition: flags.h:74
IndexType Id() const
Definition: indexed_object.h:107
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
Short class definition.
Definition: updated_lagrangian_fluid_inc.h:98
void GetValuesVector(Vector &values, int Step=0)
Definition: updated_lagrangian_fluid_inc.cpp:385
void CalculateDampingMatrix(MatrixType &rDampingMatrix, ProcessInfo &rCurrentProcessInfo)
Definition: updated_lagrangian_fluid_inc.cpp:279
void GetFirstDerivativesVector(Vector &values, int Step=0)
Definition: updated_lagrangian_fluid_inc.cpp:403
~UpdatedLagrangianFluidInc() override
Destructor.
Definition: updated_lagrangian_fluid_inc.cpp:100
void CalculateMassMatrix(MatrixType &rMassMatrix, ProcessInfo &rCurrentProcessInfo)
Definition: updated_lagrangian_fluid_inc.cpp:246
void GetDofList(DofsVectorType &ElementalDofList, ProcessInfo &CurrentProcessInfo)
Definition: updated_lagrangian_fluid_inc.cpp:368
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(UpdatedLagrangianFluidInc)
Counted pointer of UpdatedLagrangianFluidInc.
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: updated_lagrangian_fluid_inc.h:169
void GetSecondDerivativesVector(Vector &values, int Step=0)
Definition: updated_lagrangian_fluid_inc.cpp:419
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const
It creates a new element pointer.
Definition: updated_lagrangian_fluid_inc.cpp:92
void InitializeSolutionStep(ProcessInfo &CurrentProcessInfo)
Definition: updated_lagrangian_fluid_inc.cpp:339
void CalculateRightHandSide(VectorType &rRightHandSideVector, ProcessInfo &rCurrentProcessInfo)
Definition: updated_lagrangian_fluid_inc.cpp:189
void EquationIdVector(EquationIdVectorType &rResult, ProcessInfo &rCurrentProcessInfo)
Definition: updated_lagrangian_fluid_inc.cpp:350
friend class Serializer
Definition: updated_lagrangian_fluid_inc.h:248
std::string Info() const override
Turn back information as a string.
Definition: updated_lagrangian_fluid_inc.h:163
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, ProcessInfo &rCurrentProcessInfo)
Definition: updated_lagrangian_fluid_inc.cpp:106
void Calculate(const Variable< double > &rVariable, double &Output, const ProcessInfo &rCurrentProcessInfo)
Definition: updated_lagrangian_fluid_inc.cpp:436
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
void CalculateGeometryData(const GeometryType &rGeometry, const GeometryData::IntegrationMethod &rIntegrationMethod, Vector &rGaussWeights, Matrix &rNContainer, GeometryType::ShapeFunctionsGradientsType &rDN_DX)
Definition: rans_calculation_utilities.cpp:31
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
list values
Definition: bombardelli_test.py:42
ProcessInfo
Definition: edgebased_PureConvection.py:116
def load(f)
Definition: ode_solve.py:307
N
Definition: sensitivityMatrix.py:29