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.
monolithic_wall_condition.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ \.
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Jordi Cotela
11 //
12 
13 #ifndef KRATOS_MONOLITHIC_WALL_CONDITION_H
14 #define KRATOS_MONOLITHIC_WALL_CONDITION_H
15 
16 // System includes
17 #include <string>
18 #include <iostream>
19 
20 
21 // External includes
22 
23 
24 // Project includes
25 #include "includes/define.h"
26 #include "includes/serializer.h"
27 #include "includes/condition.h"
28 #include "includes/process_info.h"
29 
30 // Application includes
33 #include "includes/cfd_variables.h"
34 
35 namespace Kratos
36 {
39 
42 
46 
50 
54 
58 
60 
68 template< unsigned int TDim, unsigned int TNumNodes = TDim >
69 class KRATOS_API(FLUID_DYNAMICS_APPLICATION) MonolithicWallCondition : public Condition
70 {
71 public:
74 
77 
78  typedef Node NodeType;
79 
81 
83 
85 
86  typedef Vector VectorType;
87 
88  typedef Matrix MatrixType;
89 
90  typedef std::size_t IndexType;
91 
92  typedef std::size_t SizeType;
93 
94  typedef std::vector<std::size_t> EquationIdVectorType;
95 
96  typedef std::vector< Dof<double>::Pointer > DofsVectorType;
97 
99 
103 
105 
109  Condition(NewId)
110  {
111  }
112 
114 
119  const NodesArrayType& ThisNodes):
120  Condition(NewId,ThisNodes)
121  {
122  }
123 
125 
130  GeometryType::Pointer pGeometry):
131  Condition(NewId,pGeometry)
132  {
133  }
134 
136 
142  GeometryType::Pointer pGeometry,
143  PropertiesType::Pointer pProperties):
144  Condition(NewId,pGeometry,pProperties)
145  {
146  }
147 
150  Condition(rOther)
151  {
152  }
153 
156 
157 
161 
164  {
165  Condition::operator=(rOther);
166 
167  return *this;
168  }
169 
173 
175 
180  Condition::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override
181  {
182  return Kratos::make_intrusive<MonolithicWallCondition>(NewId, GetGeometry().Create(ThisNodes), pProperties);
183  }
184 
185  Condition::Pointer Create(IndexType NewId,
186  GeometryType::Pointer pGeom,
187  PropertiesType::Pointer pProperties) const override
188  {
189  return Kratos::make_intrusive< MonolithicWallCondition >(NewId, pGeom, pProperties);
190  }
191 
200  Condition::Pointer Clone(IndexType NewId, NodesArrayType const& rThisNodes) const override
201  {
202  Condition::Pointer pNewCondition = Create(NewId, GetGeometry().Create( rThisNodes ), pGetProperties() );
203 
204  pNewCondition->SetData(this->GetData());
205  pNewCondition->SetFlags(this->GetFlags());
206 
207  return pNewCondition;
208  }
209 
211 
214  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
215  VectorType& rRightHandSideVector,
216  const ProcessInfo& rCurrentProcessInfo) override
217  {
218  const SizeType BlockSize = TDim + 1;
219  const SizeType LocalSize = BlockSize * TNumNodes;
220 
221  if (rLeftHandSideMatrix.size1() != LocalSize)
222  rLeftHandSideMatrix.resize(LocalSize,LocalSize);
223 
224  if (rRightHandSideVector.size() != LocalSize)
225  rRightHandSideVector.resize(LocalSize);
226 
227  noalias(rLeftHandSideMatrix) = ZeroMatrix(LocalSize,LocalSize);
228  noalias(rRightHandSideVector) = ZeroVector(LocalSize);
229  }
230 
232 
236  MatrixType& rLeftHandSideMatrix,
237  const ProcessInfo& rCurrentProcessInfo) override
238  {
239  const SizeType BlockSize = TDim + 1;
240  const SizeType LocalSize = BlockSize * TNumNodes;
241 
242  if (rLeftHandSideMatrix.size1() != LocalSize)
243  rLeftHandSideMatrix.resize(LocalSize,LocalSize);
244 
245  noalias(rLeftHandSideMatrix) = ZeroMatrix(LocalSize,LocalSize);
246  }
247 
249 
252  void CalculateRightHandSide(VectorType& rRightHandSideVector,
253  const ProcessInfo& rCurrentProcessInfo) override
254  {
255  const SizeType BlockSize = TDim + 1;
256  const SizeType LocalSize = BlockSize * TNumNodes;
257 
258  if (rRightHandSideVector.size() != LocalSize)
259  rRightHandSideVector.resize(LocalSize);
260 
261  noalias(rRightHandSideVector) = ZeroVector(LocalSize);
262  }
263 
264 
265 
266  void CalculateDampingMatrix(MatrixType& rDampingMatrix,
267  const ProcessInfo& rCurrentProcessInfo) override
268  {
269  VectorType RHS;
270  this->CalculateLocalVelocityContribution(rDampingMatrix,RHS,rCurrentProcessInfo);
271  }
272 
273 
274 
276 
281  void CalculateLocalVelocityContribution(MatrixType& rDampingMatrix,
282  VectorType& rRightHandSideVector,
283  const ProcessInfo& rCurrentProcessInfo) override;
284 
285 
287  int Check(const ProcessInfo& rCurrentProcessInfo) const override
288  {
289  KRATOS_TRY;
290 
291  int Check = Condition::Check(rCurrentProcessInfo); // Checks id > 0 and area > 0
292 
293  if (Check != 0)
294  {
295  return Check;
296  }
297  else
298  {
299  // Checks on nodes
300  // Check that the element's nodes contain all required SolutionStepData and Degrees of freedom
301  for(unsigned int i=0; i<this->GetGeometry().size(); ++i)
302  {
303  if(this->GetGeometry()[i].SolutionStepsDataHas(VELOCITY) == false)
304  KRATOS_THROW_ERROR(std::invalid_argument,"missing VELOCITY variable on solution step data for node ",this->GetGeometry()[i].Id());
305  if(this->GetGeometry()[i].SolutionStepsDataHas(PRESSURE) == false)
306  KRATOS_THROW_ERROR(std::invalid_argument,"missing PRESSURE variable on solution step data for node ",this->GetGeometry()[i].Id());
307  if(this->GetGeometry()[i].SolutionStepsDataHas(MESH_VELOCITY) == false)
308  KRATOS_THROW_ERROR(std::invalid_argument,"missing MESH_VELOCITY variable on solution step data for node ",this->GetGeometry()[i].Id());
309  if(this->GetGeometry()[i].SolutionStepsDataHas(ACCELERATION) == false)
310  KRATOS_THROW_ERROR(std::invalid_argument,"missing ACCELERATION variable on solution step data for node ",this->GetGeometry()[i].Id());
311  if(this->GetGeometry()[i].SolutionStepsDataHas(EXTERNAL_PRESSURE) == false)
312  KRATOS_THROW_ERROR(std::invalid_argument,"missing EXTERNAL_PRESSURE variable on solution step data for node ",this->GetGeometry()[i].Id());
313  if(this->GetGeometry()[i].HasDofFor(VELOCITY_X) == false ||
314  this->GetGeometry()[i].HasDofFor(VELOCITY_Y) == false ||
315  this->GetGeometry()[i].HasDofFor(VELOCITY_Z) == false)
316  KRATOS_THROW_ERROR(std::invalid_argument,"missing VELOCITY component degree of freedom on node ",this->GetGeometry()[i].Id());
317  if(this->GetGeometry()[i].HasDofFor(PRESSURE) == false)
318  KRATOS_THROW_ERROR(std::invalid_argument,"missing PRESSURE component degree of freedom on node ",this->GetGeometry()[i].Id());
319  }
320 
321  return Check;
322  }
323 
324  KRATOS_CATCH("");
325  }
326 
327 
329 
334  const ProcessInfo& rCurrentProcessInfo) const override;
335 
336 
338 
342  void GetDofList(DofsVectorType& ConditionDofList,
343  const ProcessInfo& CurrentProcessInfo) const override;
344 
345 
347 
352  int Step = 0) const override
353  {
354  const SizeType LocalSize = (TDim + 1) * TNumNodes;
355  unsigned int LocalIndex = 0;
356 
357  if (Values.size() != LocalSize)
358  Values.resize(LocalSize, false);
359 
360  for (unsigned int iNode = 0; iNode < TNumNodes; ++iNode)
361  {
362  const array_1d<double,3>& rVelocity = this->GetGeometry()[iNode].FastGetSolutionStepValue(VELOCITY, Step);
363  for (unsigned int d = 0; d < TDim; ++d)
364  Values[LocalIndex++] = rVelocity[d];
365  Values[LocalIndex++] = this->GetGeometry()[iNode].FastGetSolutionStepValue(PRESSURE, Step);
366  }
367  }
368 
369 
371 
376  int Step = 0) const override
377  {
378  const SizeType LocalSize = (TDim + 1) * TNumNodes;
379  unsigned int LocalIndex = 0;
380 
381  if (Values.size() != LocalSize)
382  Values.resize(LocalSize, false);
383 
384  for (unsigned int iNode = 0; iNode < TNumNodes; ++iNode)
385  {
386  const array_1d<double,3>& rVelocity = this->GetGeometry()[iNode].FastGetSolutionStepValue(ACCELERATION, Step);
387  for (unsigned int d = 0; d < TDim; ++d)
388  Values[LocalIndex++] = rVelocity[d];
389  Values[LocalIndex++] = 0.0; // No value on pressure positions
390  }
391  }
392 
393 
395  const Variable<array_1d<double, 3 > >& rVariable,
396  std::vector<array_1d<double, 3 > >& rValues,
397  const ProcessInfo& rCurrentProcessInfo) override;
398 
400  const Variable<double>& rVariable,
401  std::vector<double>& rValues,
402  const ProcessInfo& rCurrentProcessInfo) override;
403 
405  const Variable<array_1d<double, 6 > >& rVariable,
406  std::vector<array_1d<double, 6 > >& rValues,
407  const ProcessInfo& rCurrentProcessInfo) override;
408 
410  const Variable<Vector>& rVariable,
411  std::vector<Vector>& rValues,
412  const ProcessInfo& rCurrentProcessInfo) override;
413 
415  const Variable<Matrix>& rVariable,
416  std::vector<Matrix>& rValues,
417  const ProcessInfo& rCurrentProcessInfo) override;
418 
422 
423 
427 
428 
432 
434  std::string Info() const override
435  {
436  std::stringstream buffer;
437  buffer << "MonolithicWallCondition" << TDim << "D";
438  return buffer.str();
439  }
440 
442  void PrintInfo(std::ostream& rOStream) const override
443  {
444  rOStream << "MonolithicWallCondition";
445  }
446 
448  void PrintData(std::ostream& rOStream) const override {}
449 
450 
454 
455 
457 
458 protected:
461 
462 
466 
467 
471 
472 
476 
478 
482  virtual void ApplyWallLaw(MatrixType& rLocalMatrix,
483  VectorType& rLocalVector,
484  const ProcessInfo& rCurrentProcessInfo)
485  {
486  GeometryType& rGeometry = this->GetGeometry();
487  const size_t BlockSize = TDim + 1;
488  const double NodalFactor = 1.0 / double(TDim);
489 
490  double area = NodalFactor * rGeometry.DomainSize();
491  // DomainSize() is the way to ask the geometry's length/area/volume (whatever is relevant for its dimension) without asking for the number of spatial dimensions first
492 
493  for(size_t itNode = 0; itNode < rGeometry.PointsNumber(); ++itNode)
494  {
495  const NodeType& rConstNode = rGeometry[itNode];
496  const double y = rConstNode.GetValue(Y_WALL); // wall distance to use in stress calculation
497  if( y > 0.0 && rConstNode.Is(SLIP) )
498  {
499  array_1d<double,3> Vel = rGeometry[itNode].FastGetSolutionStepValue(VELOCITY);
500  const array_1d<double,3>& VelMesh = rGeometry[itNode].FastGetSolutionStepValue(MESH_VELOCITY);
501  Vel -= VelMesh;
502  const double Ikappa = 1.0/0.41; // inverse of Von Karman's kappa
503  const double B = 5.2;
504  const double limit_yplus = 10.9931899; // limit between linear and log regions
505 
506  const double rho = rGeometry[itNode].FastGetSolutionStepValue(DENSITY);
507  const double nu = rGeometry[itNode].FastGetSolutionStepValue(VISCOSITY);
508 
509  double wall_vel = 0.0;
510  for (size_t d = 0; d < TDim; d++)
511  {
512  wall_vel += Vel[d]*Vel[d];
513  }
514  wall_vel = sqrt(wall_vel);
515 
516  if (wall_vel > 1e-12) // do not bother if velocity is zero
517  {
518 
519  // linear region
520  double utau = sqrt(wall_vel * nu / y);
521  double yplus = y * utau / nu;
522 
523  // log region
524  if (yplus > limit_yplus)
525  {
526 
527  // wall_vel / utau = 1/kappa * log(yplus) + B
528  // this requires solving a nonlinear problem:
529  // f(utau) = utau*(1/kappa * log(y*utau/nu) + B) - wall_vel = 0
530  // note that f'(utau) = 1/kappa * log(y*utau/nu) + B + 1/kappa
531 
532  unsigned int iter = 0;
533  double dx = 1e10;
534  const double tol = 1e-6;
535  double uplus = Ikappa * log(yplus) + B;
536 
537  while(iter < 100 && fabs(dx) > tol * utau)
538  {
539  // Newton-Raphson iteration
540  double f = utau * uplus - wall_vel;
541  double df = uplus + Ikappa;
542  dx = f/df;
543 
544  // Update variables
545  utau -= dx;
546  yplus = y * utau / nu;
547  uplus = Ikappa * log(yplus) + B;
548  ++iter;
549  }
550  if (iter == 100)
551  {
552  std::cout << "Warning: wall condition Newton-Raphson did not converge. Residual is " << dx << std::endl;
553  }
554  }
555  const double Tmp = area * utau * utau * rho / wall_vel;
556  for (size_t d = 0; d < TDim; d++)
557  {
558  size_t k = itNode*BlockSize+d;
559  rLocalVector[k] -= Vel[d] * Tmp;
560  rLocalMatrix(k,k) += Tmp;
561  }
562  }
563  }
564  }
565  }
566 
567 
569 
570 
571  void ApplyNeumannCondition(MatrixType &rLocalMatrix, VectorType &rLocalVector);
572 
573 
577 
578 
582 
583 
587 
588 
590 
591 private:
594 
595 
599 
600 
604 
605  friend class Serializer;
606 
607  void save(Serializer& rSerializer) const override
608  {
610  }
611 
612  void load(Serializer& rSerializer) override
613  {
615  }
616 
620 
621 
625 
626 
630 
631 
635 
636 
640 
641 
643 
644 }; // Class MonolithicWallCondition
645 
646 
648 
651 
652 
656 
657 
659 template< unsigned int TDim, unsigned int TNumNodes >
660 inline std::istream& operator >> (std::istream& rIStream,
662 {
663  return rIStream;
664 }
665 
667 template< unsigned int TDim, unsigned int TNumNodes >
668 inline std::ostream& operator << (std::ostream& rOStream,
670 {
671  rThis.PrintInfo(rOStream);
672  rOStream << std::endl;
673  rThis.PrintData(rOStream);
674 
675  return rOStream;
676 }
677 
679 
681 
682 
683 } // namespace Kratos.
684 
685 #endif // KRATOS_MONOLITHIC_WALL_CONDITION_H
Base class for all Conditions.
Definition: condition.h:59
std::size_t SizeType
Definition: condition.h:94
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
std::vector< DofType::Pointer > DofsVectorType
Definition: condition.h:100
virtual int Check(const ProcessInfo &rCurrentProcessInfo) const
Definition: condition.h:854
Condition & operator=(Condition const &rOther)
Assignment operator.
Definition: condition.h:181
std::size_t IndexType
Definition: flags.h:74
bool Is(Flags const &rOther) const
Definition: flags.h:274
Geometry base class.
Definition: geometry.h:71
SizeType PointsNumber() const
Definition: geometry.h:528
virtual double DomainSize() const
This method calculate and return length, area or volume of this geometry depending to it's dimension.
Definition: geometry.h:1371
This object defines an indexed object.
Definition: indexed_object.h:54
void resize(std::size_t NewSize1, std::size_t NewSize2, bool preserve=0)
Definition: amatrix_interface.h:224
Implements a wall condition for the monolithic formulation.
Definition: monolithic_wall_condition.h:70
MonolithicWallCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Constructor using Properties.
Definition: monolithic_wall_condition.h:141
Properties PropertiesType
Definition: monolithic_wall_condition.h:80
Condition::Pointer Clone(IndexType NewId, NodesArrayType const &rThisNodes) const override
Definition: monolithic_wall_condition.h:200
~MonolithicWallCondition() override
Destructor.
Definition: monolithic_wall_condition.h:155
Condition::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
It creates a new condition pointer.
Definition: monolithic_wall_condition.h:185
virtual void ApplyWallLaw(MatrixType &rLocalMatrix, VectorType &rLocalVector, const ProcessInfo &rCurrentProcessInfo)
Commpute the wall stress and add corresponding terms to the system contributions.
Definition: monolithic_wall_condition.h:482
void CalculateNormal(array_1d< double, 3 > &An)
PointerVectorSet< Dof< double >, IndexedObject > DofsArrayType
Definition: monolithic_wall_condition.h:98
void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Return a matrix of the correct size, filled with zeros (for compatibility with time schemes).
Definition: monolithic_wall_condition.h:235
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
Provides the global indices for each one of this element's local rows.
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: monolithic_wall_condition.h:448
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Check that all data required by this condition is available and reasonable.
Definition: monolithic_wall_condition.h:287
MonolithicWallCondition(IndexType NewId=0)
Default constructor.
Definition: monolithic_wall_condition.h:108
MonolithicWallCondition & operator=(MonolithicWallCondition const &rOther)
Assignment operator.
Definition: monolithic_wall_condition.h:163
Geometry< NodeType > GeometryType
Definition: monolithic_wall_condition.h:82
void CalculateDampingMatrix(MatrixType &rDampingMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: monolithic_wall_condition.h:266
void GetSecondDerivativesVector(Vector &Values, int Step=0) const override
Returns ACCELERATION_X, ACCELERATION_Y, (ACCELERATION_Z,) 0 for each node.
Definition: monolithic_wall_condition.h:375
std::size_t SizeType
Definition: monolithic_wall_condition.h:92
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(MonolithicWallCondition)
Pointer definition of MonolithicWallCondition.
Condition::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Create a new MonolithicWallCondition object.
Definition: monolithic_wall_condition.h:180
Geometry< NodeType >::PointsArrayType NodesArrayType
Definition: monolithic_wall_condition.h:84
Node NodeType
Definition: monolithic_wall_condition.h:78
MonolithicWallCondition(MonolithicWallCondition const &rOther)
Copy constructor.
Definition: monolithic_wall_condition.h:149
void GetDofList(DofsVectorType &ConditionDofList, const ProcessInfo &CurrentProcessInfo) const override
Returns a list of the element's Dofs.
std::vector< std::size_t > EquationIdVectorType
Definition: monolithic_wall_condition.h:94
void GetFirstDerivativesVector(Vector &Values, int Step=0) const override
Returns VELOCITY_X, VELOCITY_Y, (VELOCITY_Z,) PRESSURE for each node.
Definition: monolithic_wall_condition.h:351
std::size_t IndexType
Definition: monolithic_wall_condition.h:90
Matrix MatrixType
Definition: monolithic_wall_condition.h:88
MonolithicWallCondition(IndexType NewId, const NodesArrayType &ThisNodes)
Constructor using an array of nodes.
Definition: monolithic_wall_condition.h:118
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Return local contributions of the correct size, filled with zeros (for compatibility with time scheme...
Definition: monolithic_wall_condition.h:214
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: monolithic_wall_condition.h:442
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Return local right hand side of the correct size, filled with zeros (for compatibility with time sche...
Definition: monolithic_wall_condition.h:252
std::string Info() const override
Turn back information as a string.
Definition: monolithic_wall_condition.h:434
std::vector< Dof< double >::Pointer > DofsVectorType
Definition: monolithic_wall_condition.h:96
Vector VectorType
Definition: monolithic_wall_condition.h:86
MonolithicWallCondition(IndexType NewId, GeometryType::Pointer pGeometry)
Constructor using Geometry.
Definition: monolithic_wall_condition.h:129
This class defines the node.
Definition: node.h:65
TVariableType::Type & GetValue(const TVariableType &rThisVariable)
Definition: node.h:466
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
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
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#define KRATOS_THROW_ERROR(ExceptionType, ErrorMessage, MoreInfo)
Definition: define.h:77
#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
Matrix MatrixType
Definition: geometrical_transformation_utilities.h:55
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
pybind11::list CalculateOnIntegrationPoints(TObject &dummy, const Variable< TDataType > &rVariable, const ProcessInfo &rProcessInfo)
Definition: add_mesh_to_python.cpp:142
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
KratosZeroMatrix< double > ZeroMatrix
Definition: amatrix_interface.h:559
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
TABLE_NUMBER_ANGULAR_VELOCITY TABLE_NUMBER_MOMENT I33 BEAM_INERTIA_ROT_UNIT_LENGHT_Y KRATOS_DEFINE_APPLICATION_VARIABLE(DEM_APPLICATION, double, BEAM_INERTIA_ROT_UNIT_LENGHT_Z) typedef std double
Definition: DEM_application_variables.h:182
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
y
Other simbols definition.
Definition: generate_axisymmetric_navier_stokes_element.py:54
f
Definition: generate_convection_diffusion_explicit_element.py:112
rho
Definition: generate_droplet_dynamics.py:86
int tol
Definition: hinsberg_optimization.py:138
nu
Definition: isotropic_damage_automatic_differentiation.py:135
def load(f)
Definition: ode_solve.py:307
int d
Definition: ode_solve.py:397
int k
Definition: quadrature.py:595
B
Definition: sensitivityMatrix.py:76
integer i
Definition: TensorModule.f:17
e
Definition: run_cpp_mpi_tests.py:31