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.
ulf_axisym.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-03-23 13:03:13 $
47 // Revision: $Revision: 1.3 $
48 //
49 //
50 
51 
52 #if !defined(KRATOS_ULF_AXISYM_H_INCLUDED )
53 #define KRATOS_ULF_AXISYM_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 
96 class UlfAxisym
97  : public Element
98 {
99 public:
102 
105 
109 
111  UlfAxisym(IndexType NewId, GeometryType::Pointer pGeometry);
112  UlfAxisym(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
113 
115  ~UlfAxisym() 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& rDampMatrix, 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 "UlfAxisym #" ;
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:
224 
229  double mA0; //stores area at the beginning of time-step
230  double r0;// stores the radius at the beginning of time-step
231 
235  friend class Serializer;
236 
237  // A private default constructor necessary for serialization
238  UlfAxisym() : Element()
239  {
240  }
241 
242  void save(Serializer& rSerializer) const override
243  {
245  }
246 
247  void load(Serializer& rSerializer) override
248  {
250  }
251 
252 
253 
254 
258 
259 
260 
261  //inline void CalculateGeometryData(Matrix& msDN_DX, Vector& N, double& Area)
262  inline void CalculateGeometryData(BoundedMatrix<double,3,2>& DN_DX, array_1d<double,3>& N, double& Area);
263 
267 
268 
272 
273 
277 
278 
282 
284  //UlfAxisym& operator=(const UlfAxisym& rOther);
285 
287  //UlfAxisym(const UlfAxisym& rOther);
288 
289 
291 
292 }; // Class UlfAxisym
293 
295 
298 
299 
303 
304 
306 /* inline std::istream& operator >> (std::istream& rIStream,
307  UlfAxisym& rThis);
308 */
310 /* inline std::ostream& operator << (std::ostream& rOStream,
311  const UlfAxisym& rThis)
312  {
313  rThis.PrintInfo(rOStream);
314  rOStream << std::endl;
315  rThis.PrintData(rOStream);
316 
317  return rOStream;
318  }*/
320 
321 } // namespace Kratos.
322 
323 #endif // KRATOS_ULF_AXISYM_H_INCLUDED defined
324 
325 
Base class for all Elements.
Definition: element.h:60
Element(IndexType NewId=0)
Definition: element.h:121
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
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: ulf_axisym.h:98
void EquationIdVector(EquationIdVectorType &rResult, ProcessInfo &rCurrentProcessInfo)
Definition: ulf_axisym.cpp:458
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: ulf_axisym.h:169
void Calculate(const Variable< double > &rVariable, double &Output, const ProcessInfo &rCurrentProcessInfo)
Definition: ulf_axisym.cpp:544
void GetFirstDerivativesVector(Vector &values, int Step=0)
Definition: ulf_axisym.cpp:511
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, ProcessInfo &rCurrentProcessInfo)
Definition: ulf_axisym.cpp:102
void CalculateRightHandSide(VectorType &rRightHandSideVector, ProcessInfo &rCurrentProcessInfo)
Definition: ulf_axisym.cpp:273
std::string Info() const override
Turn back information as a string.
Definition: ulf_axisym.h:163
void GetValuesVector(Vector &values, int Step=0)
Definition: ulf_axisym.cpp:493
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(UlfAxisym)
Counted pointer of UlfAxisym.
void GetSecondDerivativesVector(Vector &values, int Step=0)
Definition: ulf_axisym.cpp:527
void InitializeSolutionStep(ProcessInfo &CurrentProcessInfo)
Definition: ulf_axisym.cpp:444
~UlfAxisym() override
Destructor.
Definition: ulf_axisym.cpp:96
void CalculateMassMatrix(MatrixType &rMassMatrix, ProcessInfo &rCurrentProcessInfo)
Definition: ulf_axisym.cpp:283
void CalculateDampingMatrix(MatrixType &rDampMatrix, ProcessInfo &rCurrentProcessInfo)
Definition: ulf_axisym.cpp:338
void GetDofList(DofsVectorType &ElementalDofList, ProcessInfo &CurrentProcessInfo)
Definition: ulf_axisym.cpp:476
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const
It creates a new element pointer.
Definition: ulf_axisym.cpp:88
friend class Serializer
Definition: ulf_axisym.h:235
#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
def load(f)
Definition: ode_solve.py:307
N
Definition: sensitivityMatrix.py:29