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.
flux_condition.h
Go to the documentation of this file.
1 // KRATOS ___ ___ _ ___ __ ___ ___ ___ ___
2 // / __/ _ \| \| \ \ / /__| \_ _| __| __|
3 // | (_| (_) | .` |\ V /___| |) | || _|| _|
4 // \___\___/|_|\_| \_/ |___/___|_| |_| APPLICATION
5 //
6 // License: BSD License
7 // Kratos default license: kratos/license.txt
8 //
9 // Main authors: Jordi Cotela
10 //
11 
12 #if !defined(KRATOS_FLUX_CONDITION_H_INCLUDED )
13 #define KRATOS_FLUX_CONDITION_H_INCLUDED
14 
15 // System includes
16 #include <string>
17 #include <iostream>
18 
19 
20 // External includes
21 
22 
23 // Project includes
24 #include "includes/define.h"
25 #include "includes/node.h"
26 #include "includes/process_info.h"
27 #include "includes/properties.h"
28 #include "includes/condition.h"
29 #include "geometries/geometry.h"
30 #include "includes/variables.h"
31 
32 namespace Kratos
33 {
34 
37 
38 namespace FluxConditionInternals
39 {
42 
43 template< unsigned int TNodeNumber >
45 {
46 public:
47 
49  Geometry< Node >& rGeometry,
50  const Variable<double>& rFluxVar
51  ):
52  mGaussPoint(0),
53  mNodalFluxes(TNodeNumber,0.0)
54  {
56  Vector DetJ = ZeroVector(NumGauss);
58 
59  mShapeFunctionValues.resize(NumGauss,TNodeNumber);
60  mIntegrationWeights.resize(NumGauss);
61 
63 
64  const auto& IntegrationPoints = rGeometry.IntegrationPoints(GeometryData::IntegrationMethod::GI_GAUSS_2);
65 
66  for (unsigned int g = 0; g < NumGauss; g++)
67  {
68  mIntegrationWeights[g] = DetJ[g]*IntegrationPoints[g].Weight();
69  }
70 
71  for (unsigned int i = 0; i < TNodeNumber; i++)
72  {
73  mNodalFluxes[i] = rGeometry[i].FastGetSolutionStepValue(rFluxVar);
74  }
75  }
76 
77  void SetCurrentGaussPoint(unsigned int g)
78  {
79  mGaussPoint = g;
80  }
81 
82  double N(unsigned int i) const
83  {
84  return mShapeFunctionValues(mGaussPoint,i);
85  }
86 
87  double NodalFlux(unsigned int i) const
88  {
89  return mNodalFluxes[i];
90  }
91 
92  double GaussPointFlux() const
93  {
94  double flux = mNodalFluxes[0]*mShapeFunctionValues(mGaussPoint,0);
95  for (unsigned int i = 1; i < TNodeNumber; i++)
96  {
97  flux += mNodalFluxes[i]*mShapeFunctionValues(mGaussPoint,i);
98  }
99  return flux;
100  }
101 
102  double IntegrationWeight() const
103  {
104  return mIntegrationWeights[mGaussPoint];
105  }
106 
107  unsigned int NumGauss;
108 
109 private:
110 
111  unsigned int mGaussPoint;
112 
113  array_1d<double,TNodeNumber> mNodalFluxes;
114 
115  Matrix mShapeFunctionValues;
116 
117  Vector mIntegrationWeights;
118 };
119 
121 
122 }
123 
127 
129 
133 template< unsigned int TNodeNumber >
135 {
136 public:
139 
142 
145 
149 
151  IndexType NewId,
152  Geometry< Node >::Pointer pGeometry);
153 
155  IndexType NewId,
156  Geometry< Node >::Pointer pGeometry,
157  Properties::Pointer pProperties);
158 
160  ~FluxCondition() override;
161 
165 
166  Condition::Pointer Create(
167  IndexType NewId,
168  NodesArrayType const& ThisNodes,
169  Properties::Pointer pProperties) const override;
170 
171  Condition::Pointer Create(
172  IndexType NewId,
173  GeometryType::Pointer pGeom,
174  Properties::Pointer pProperties) const override;
175 
177  MatrixType& rLeftHandSideMatrix,
178  VectorType& rRightHandSideVector,
179  const ProcessInfo& rCurrentProcessInfo) override;
180 
182  VectorType& rRightHandSideVector,
183  const ProcessInfo& rCurrentProcessInfo) override;
184 
185  void EquationIdVector(
186  EquationIdVectorType& rResult,
187  const ProcessInfo& rCurrentProcessInfo) const override;
188 
189  void GetDofList(
190  DofsVectorType& ConditionalDofList,
191  const ProcessInfo& CurrentProcessInfo) const override;
192 
194 
196  const Variable<array_1d<double, 3 > >& rVariable,
197  std::vector<array_1d<double, 3 > >& rValues,
198  const ProcessInfo& rCurrentProcessInfo) override;
199 
201  const Variable<double>& rVariable,
202  std::vector<double>& rValues,
203  const ProcessInfo& rCurrentProcessInfo) override;
204 
206  const Variable<array_1d<double, 6 > >& rVariable,
207  std::vector<array_1d<double, 6 > >& rValues,
208  const ProcessInfo& rCurrentProcessInfo) override;
209 
211  const Variable<Vector>& rVariable,
212  std::vector<Vector>& rValues,
213  const ProcessInfo& rCurrentProcessInfo) override;
214 
216  const Variable<Matrix>& rVariable,
217  std::vector<Matrix>& rValues,
218  const ProcessInfo& rCurrentProcessInfo) override;
219 
223 
225  std::string Info() const override;
226 
228  void PrintInfo(std::ostream& rOStream) const override;
229 
231  void PrintData(std::ostream& rOStream) const override;
232 
234 
235 protected:
236 
239 
240  // Default constructor necessary for serialization
241  FluxCondition();
242 
246 
248  VectorType& rRightHandSideVector,
250 
252 
254 
255 private:
256 
259  friend class Serializer;
260 
261  void save(Serializer& rSerializer) const override;
262 
263  void load(Serializer& rSerializer) override;
264 
268 
270  FluxCondition& operator=(FluxCondition const& rOther);
271 
273  FluxCondition(FluxCondition const& rOther);
274 
276 
277 }; // Class FluxCondition
278 
280 
282 
283 } // namespace Kratos.
284 
285 #endif // KRATOS_FLUX_CONDITION_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Base class for all Conditions.
Definition: condition.h:59
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
std::vector< DofType::Pointer > DofsVectorType
Definition: condition.h:100
std::size_t IndexType
Definition: flags.h:74
A basic Neumann condition for convection-diffusion problems.
Definition: flux_condition.h:135
Condition::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, Properties::Pointer pProperties) const override
Definition: flux_condition.cpp:45
Condition::VectorType VectorType
Definition: flux_condition.h:144
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(FluxCondition)
Pointer definition of FluxCondition.
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: flux_condition.cpp:82
Condition::MatrixType MatrixType
Definition: flux_condition.h:143
GeometryData::IntegrationMethod GetIntegrationMethod() const override
Definition: flux_condition.cpp:166
void CalculateOnIntegrationPoints(const Variable< array_1d< double, 3 > > &rVariable, std::vector< array_1d< double, 3 > > &rValues, const ProcessInfo &rCurrentProcessInfo) override
Definition: flux_condition.cpp:172
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: flux_condition.cpp:282
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: flux_condition.cpp:276
FluxCondition()
Definition: flux_condition.cpp:343
void CalculateNormal(array_1d< double, 3 > &An)
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: flux_condition.cpp:63
~FluxCondition() override
Destructor.
Definition: flux_condition.cpp:38
void AddIntegrationPointRHSContribution(VectorType &rRightHandSideVector, const FluxConditionInternals::IntegrationData< TNodeNumber > &rData)
Definition: flux_condition.cpp:292
std::string Info() const override
Turn back information as a string.
Definition: flux_condition.cpp:268
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
Definition: flux_condition.cpp:112
void GetDofList(DofsVectorType &ConditionalDofList, const ProcessInfo &CurrentProcessInfo) const override
Definition: flux_condition.cpp:139
double GaussPointFlux() const
Definition: flux_condition.h:92
void SetCurrentGaussPoint(unsigned int g)
Definition: flux_condition.h:77
IntegrationData(Geometry< Node > &rGeometry, const Variable< double > &rFluxVar)
Definition: flux_condition.h:48
double IntegrationWeight() const
Definition: flux_condition.h:102
double NodalFlux(unsigned int i) const
Definition: flux_condition.h:87
unsigned int NumGauss
Definition: flux_condition.h:107
double N(unsigned int i) const
Definition: flux_condition.h:82
IntegrationMethod
Definition: geometry_data.h:76
Geometry base class.
Definition: geometry.h:71
const Matrix & ShapeFunctionsValues() const
Definition: geometry.h:3393
Vector & DeterminantOfJacobian(Vector &rResult) const
Definition: geometry.h:3154
const IntegrationPointsArrayType & IntegrationPoints() const
Definition: geometry.h:2284
SizeType IntegrationPointsNumber() const
Definition: geometry.h:2257
void resize(std::size_t NewSize1, std::size_t NewSize2, bool preserve=0)
Definition: amatrix_interface.h:224
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
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
def load(f)
Definition: ode_solve.py:307
integer i
Definition: TensorModule.f:17