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.
two_step_updated_lagrangian_V_P_implicit_element.h
Go to the documentation of this file.
1 //
2 // Project Name: KratosFluidDynamicsApplication $
3 // Last modified by: $Author: AFranci $
4 // Date: $Date: January 2016 $
5 // Revision: $Revision: 0.0 $
6 //
7 //
8 
9 #if !defined(KRATOS_TWO_STEP_UPDATED_LAGRANGIAN_V_P_IMPLICIT_ELEMENT_H_INCLUDED)
10 #define KRATOS_TWO_STEP_UPDATED_LAGRANGIAN_V_P_IMPLICIT_ELEMENT_H_INCLUDED
11 
12 // System includes
13 #include <string>
14 #include <iostream>
15 
16 // External includes
17 
18 // Project includes
19 #include "containers/array_1d.h"
20 #include "includes/define.h"
21 #include "includes/serializer.h"
22 #include "geometries/geometry.h"
23 #include "utilities/math_utils.h"
24 #include "utilities/geometry_utilities.h"
25 
27 
29 
30 #include "includes/model_part.h"
31 
32 namespace Kratos
33 {
34 
37 
40 
44 
48 
52 
56 
58 
61  template <unsigned int TDim>
63  {
64 
65  public:
68 
71 
74 
76  typedef Node NodeType;
77 
80 
83 
86 
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 
102 
105 
108 
109  /* typedef Element::PropertiesType::Pointer PropertiesType::Pointer; */
110 
112 
113  typedef typename BaseType::PropertiesType::Pointer pPropertiesType;
114 
116 
119 
121  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
122 
126 
127  // Constructors.
128 
130 
134  {
135  }
136 
138 
142  TwoStepUpdatedLagrangianVPImplicitElement(IndexType NewId, const NodesArrayType &ThisNodes) : BaseType(NewId, ThisNodes)
143  {
144  }
145 
147 
151  TwoStepUpdatedLagrangianVPImplicitElement(IndexType NewId, GeometryType::Pointer pGeometry) : BaseType(NewId, pGeometry)
152  {
153  }
154 
156 
161  TwoStepUpdatedLagrangianVPImplicitElement(IndexType NewId, GeometryType::Pointer pGeometry, pPropertiesType pProperties) : BaseType(NewId, pGeometry, pProperties)
162  {
163  }
164 
166 
168  {
169  }
170 
173  {
174  }
175 
179 
183 
185 
192  Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes,
193  pPropertiesType pProperties) const override
194  {
195  return Element::Pointer(new TwoStepUpdatedLagrangianVPImplicitElement(NewId, this->GetGeometry().Create(ThisNodes), pProperties));
196  }
197 
198  Element::Pointer Clone(IndexType NewId, NodesArrayType const &ThisNodes) const override;
199 
200  void Initialize(const ProcessInfo &rCurrentProcessInfo) override{};
201 
203  void InitializeSolutionStep(const ProcessInfo &rCurrentProcessInfo) override{};
204 
205  void InitializeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo) override{};
206 
208  void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix,
209  VectorType &rRightHandSideVector,
210  const ProcessInfo &rCurrentProcessInfo) override;
211 
212  void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix,
213  const ProcessInfo &rCurrentProcessInfo) override
214  {
215  KRATOS_TRY;
216  KRATOS_THROW_ERROR(std::logic_error, "TwoStepUpdatedLagrangianVPImplicitElement::CalculateLeftHandSide not implemented", "");
217  KRATOS_CATCH("");
218  }
219 
220  void CalculateRightHandSide(VectorType &rRightHandSideVector,
221  const ProcessInfo &rCurrentProcessInfo) override
222  {
223  KRATOS_TRY;
224  KRATOS_THROW_ERROR(std::logic_error, "TwoStepUpdatedLagrangianVPImplicitElement::CalculateRightHandSide not implemented", "");
225  KRATOS_CATCH("");
226  }
227 
228  void CalculateOnIntegrationPoints(const Variable<bool> &rVariable,
229  std::vector<bool> &rOutput,
230  const ProcessInfo &rCurrentProcessInfo) override;
231 
232  void CalculateOnIntegrationPoints(const Variable<double> &rVariable,
233  std::vector<double> &rOutput,
234  const ProcessInfo &rCurrentProcessInfo) override;
235 
236  void CalculateOnIntegrationPoints(const Variable<Vector> &rVariable,
237  std::vector<Vector> &rOutput,
238  const ProcessInfo &rCurrentProcessInfo) override;
239 
240  /* // The following methods have different implementations depending on TDim */
241  /* /// Provides the global indices for each one of this element's local rows */
242  /* /\** */
243  /* * this determines the elemental equation ID vector for all elemental */
244  /* * DOFs */
245  /* * @param rResult A vector containing the global Id of each row */
246  /* * @param rCurrentProcessInfo the current process info object (unused) */
247  /* *\/ */
248 
249  /* virtual GeometryData::IntegrationMethod GetIntegrationMethod() const; */
250 
251  void UpdateCauchyStress(unsigned int g, const ProcessInfo &rCurrentProcessInfo) override{};
252 
253  void InitializeElementalVariables(ElementalVariables &rElementalVariables) override
254  {
255  KRATOS_TRY;
256  KRATOS_THROW_ERROR(std::logic_error, "InitializeElementalVariables", "");
257  KRATOS_CATCH("");
258  };
259 
260  /* void CalculateDeltaPosition (Matrix & rDeltaPosition); */
261 
265 
269 
271 
279  int Check(const ProcessInfo &rCurrentProcessInfo) const override;
280 
284 
288 
290  std::string Info() const override
291  {
292  std::stringstream buffer;
293  buffer << "TwoStepUpdatedLagrangianVPImplicitElement #" << BaseType::Id();
294  return buffer.str();
295  }
296 
298  void PrintInfo(std::ostream &rOStream) const override
299  {
300  rOStream << "TwoStepUpdatedLagrangianVPImplicitElement" << TDim << "D";
301  }
302 
303  // /// Print object's data.
304  // virtual void PrintData(std::ostream& rOStream) const;
305 
309 
311  protected:
314 
318 
319  ConstitutiveLaw::Pointer mpConstitutiveLaw = nullptr;
320 
324 
328 
329  void CalculateLocalMomentumEquations(MatrixType &rLeftHandSideMatrix,
330  VectorType &rRightHandSideVector,
331  const ProcessInfo &rCurrentProcessInfo) override;
332 
334  VectorType &rRightHandSideVector,
335  const ProcessInfo &rCurrentProcessInfo) override{};
336 
337  double GetThetaMomentum() override
338  {
339  std::cout << "I SHOULD NOT ENTER HERE!" << std::endl;
340  return 1.0;
341  };
342 
343  double GetThetaContinuity() override
344  {
345  std::cout << "I SHOULD NOT ENTER HERE!" << std::endl;
346  return 1.0;
347  };
348 
350  double &MeanValue,
351  const ShapeFunctionDerivativesType &rShapeDeriv,
352  const double secondLame,
353  double &bulkModulus,
354  const double Weight,
355  double &MeanValueMass,
356  const double TimeStep){};
357 
359  MatrixType StiffnessMatrix,
360  double &meanValueStiff,
361  double &bulkCoefficient,
362  double timeStep){};
363 
365  MatrixType &rDampingMatrix,
366  const ShapeFunctionDerivativesType &rShapeDeriv,
367  const double secondLame,
368  const double bulkModulus,
369  const double theta,
370  const double Weight);
371 
372  virtual void ComputeBulkMatrixLump(MatrixType &BulkMatrix,
373  const double Weight) override{};
374 
375  virtual void ComputeBoundLHSMatrix(MatrixType &BoundLHSMatrix,
376  const ShapeFunctionsType &rN,
377  const double Weight){};
378 
379  virtual void ComputeBoundRHSVector(VectorType &BoundRHSVector,
380  const ShapeFunctionsType &rN,
381  const double TimeStep,
382  const double BoundRHSCoeffAcc,
383  const double BoundRHSCoeffDev){};
384 
385  virtual void ComputeStabLaplacianMatrix(MatrixType &StabLaplacianMatrix,
386  const ShapeFunctionDerivativesType &rShapeDeriv,
387  const double Weight){};
388 
390  ElementalVariables &rElementalVariables,
391  const unsigned int g,
392  const Vector& rN,
393  const ProcessInfo &rCurrentProcessInfo,
394  double &Density,
395  double &DeviatoricCoeff,
396  double &VolumetricCoeff) override
397  {};
398 
399  virtual void CalculateTauFIC(double &TauOne,
400  double ElemSize,
401  const double Density,
402  const double Viscosity,
403  const ProcessInfo &rCurrentProcessInfo){};
404 
405  virtual void AddStabilizationMatrixLHS(MatrixType &rLeftHandSideMatrix,
406  Matrix &BulkAccMatrix,
407  const ShapeFunctionsType &rN,
408  const double Weight){};
409 
410  virtual void AddStabilizationNodalTermsLHS(MatrixType &rLeftHandSideMatrix,
411  const double Tau,
412  const double Weight,
413  const ShapeFunctionDerivativesType &rDN_DX,
414  const SizeType i){};
415 
416  virtual void AddStabilizationNodalTermsRHS(VectorType &rRightHandSideVector,
417  const double Tau,
418  const double Density,
419  const double Weight,
420  const ShapeFunctionDerivativesType &rDN_DX,
421  const SizeType i){};
422 
426 
430 
434 
436 
437  private:
440 
444 
448 
449  friend class Serializer;
450 
451  void save(Serializer &rSerializer) const override
452  {
454  }
455 
456  void load(Serializer &rSerializer) override
457  {
459  }
460 
464 
468 
472 
476 
480 
483 
484  /* /// Copy constructor. */
485  /* TwoStepUpdatedLagrangianVPImplicitElement(TwoStepUpdatedLagrangianVPImplicitElement const& rOther); */
486 
488 
489  }; // Class TwoStepUpdatedLagrangianVPImplicitElement
490 
492 
495 
499 
501  template <unsigned int TDim>
502  inline std::istream &operator>>(std::istream &rIStream,
504  {
505  return rIStream;
506  }
507 
509  template <unsigned int TDim>
510  inline std::ostream &operator<<(std::ostream &rOStream,
512  {
513  rThis.PrintInfo(rOStream);
514  rOStream << std::endl;
515  rThis.PrintData(rOStream);
516 
517  return rOStream;
518  }
519 
520 } // namespace Kratos.
521 
522 #endif // KRATOS_TWO_STEP_UPDATED_LAGRANGIAN_V_P_ELEMENT defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Definition: constitutive_law.h:47
Base class for all Elements.
Definition: element.h:60
std::size_t SizeType
Definition: element.h:94
virtual void MassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:926
std::size_t IndexType
Definition: flags.h:74
This defines the geometrical object, base definition of the element and condition entities.
Definition: geometrical_object.h:58
GeometryType & GetGeometry()
Returns the reference of the geometry.
Definition: geometrical_object.h:158
Geometry base class.
Definition: geometry.h:71
This object defines an indexed object.
Definition: indexed_object.h:54
IndexType Id() const
Definition: indexed_object.h:107
This class defines the node.
Definition: node.h:65
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
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
A stabilized element for the incompressible Navier-Stokes equations.
Definition: two_step_updated_lagrangian_element.h:63
A stabilized element for the incompressible Navier-Stokes equations.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:63
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
void CalcElasticPlasticCauchySplitted(ElementalVariables &rElementalVariables, const unsigned int g, const Vector &rN, const ProcessInfo &rCurrentProcessInfo, double &Density, double &DeviatoricCoeff, double &VolumetricCoeff) override
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:389
std::size_t SizeType
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:92
std::vector< Dof< double >::Pointer > DofsVectorType
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:96
BaseType::PropertiesType::Pointer pPropertiesType
Definition: two_step_updated_lagrangian_element.h:116
Kratos::Matrix MatrixType
Matrix type for local contributions to the linear system.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:88
std::size_t IndexType
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:90
virtual void CalculateTauFIC(double &TauOne, double ElemSize, const double Density, const double Viscosity, const ProcessInfo &rCurrentProcessInfo)
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:399
std::string Info() const override
Turn back information as a string.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:290
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, pPropertiesType pProperties) const override
Create a new element of this type.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:192
virtual void ComputeBoundLHSMatrix(MatrixType &BoundLHSMatrix, const ShapeFunctionsType &rN, const double Weight)
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:375
void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:212
TwoStepUpdatedLagrangianVPImplicitElement(IndexType NewId, const NodesArrayType &ThisNodes)
Constructor using an array of nodes.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:142
Geometry< NodeType > GeometryType
Geometry type (using with given NodeType)
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:79
ConstitutiveLaw::Pointer mpConstitutiveLaw
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:319
void InitializeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo) override
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:205
virtual void ComputeBoundRHSVector(VectorType &BoundRHSVector, const ShapeFunctionsType &rN, const double TimeStep, const double BoundRHSCoeffAcc, const double BoundRHSCoeffDev)
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:379
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:121
Kratos::Vector VectorType
Vector type for local contributions to the linear system.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:85
Kratos::Matrix ShapeFunctionDerivativesType
Type for a matrix containing the shape function gradients.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:104
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:298
BaseType::PropertiesType PropertiesType
Definition: two_step_updated_lagrangian_element.h:114
PointerVectorSet< Dof< double >, IndexedObject > DofsArrayType
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:98
void CalculateOnIntegrationPoints(const Variable< bool > &rVariable, std::vector< bool > &rOutput, const ProcessInfo &rCurrentProcessInfo) override
Definition: two_step_updated_lagrangian_V_P_implicit_element.cpp:215
GeometryType::ShapeFunctionsGradientsType ShapeFunctionDerivativesArrayType
Type for an array of shape function gradient matrices.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:107
void CalculateLocalContinuityEqForPressure(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:333
virtual void AddStabilizationNodalTermsLHS(MatrixType &rLeftHandSideMatrix, const double Tau, const double Weight, const ShapeFunctionDerivativesType &rDN_DX, const SizeType i)
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:410
BaseType::ElementalVariables ElementalVariables
Definition: two_step_updated_lagrangian_element.h:118
BaseType::PropertiesType::Pointer pPropertiesType
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:113
void InitializeElementalVariables(ElementalVariables &rElementalVariables) override
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:253
Node NodeType
Node type (default is: Node)
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:76
TwoStepUpdatedLagrangianElement< TDim > BaseType
base type:
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:73
Kratos::Vector ShapeFunctionsType
Type for shape function values container.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:101
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:220
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(TwoStepUpdatedLagrangianVPImplicitElement)
Pointer definition of TwoStepUpdatedLagrangianVPImplicitElement.
void ComputeCompleteTangentTerm(ElementalVariables &rElementalVariables, MatrixType &rDampingMatrix, const ShapeFunctionDerivativesType &rShapeDeriv, const double secondLame, const double bulkModulus, const double theta, const double Weight)
ConstitutiveLaw ConstitutiveLawType
Reference type definition for constitutive laws.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:118
virtual void ComputeBulkMatrixLump(MatrixType &BulkMatrix, const double Weight) override
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:372
std::vector< std::size_t > EquationIdVectorType
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:94
BaseType::PropertiesType PropertiesType
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:111
virtual ~TwoStepUpdatedLagrangianVPImplicitElement()
Destructor.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:172
TwoStepUpdatedLagrangianVPImplicitElement(IndexType NewId=0)
Default constuctor.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:133
TwoStepUpdatedLagrangianVPImplicitElement(IndexType NewId, GeometryType::Pointer pGeometry, pPropertiesType pProperties)
Constuctor using geometry and properties.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:161
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:200
Element::Pointer Clone(IndexType NewId, NodesArrayType const &ThisNodes) const override
It creates a new element pointer and clones the previous element data.
Definition: two_step_updated_lagrangian_V_P_implicit_element.cpp:21
Geometry< NodeType >::PointsArrayType NodesArrayType
Definition of nodes container type, redefined from GeometryType.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:82
virtual void AddStabilizationMatrixLHS(MatrixType &rLeftHandSideMatrix, Matrix &BulkAccMatrix, const ShapeFunctionsType &rN, const double Weight)
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:405
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Checks the input and that all required Kratos variables have been registered.
Definition: two_step_updated_lagrangian_V_P_implicit_element.cpp:354
void ComputeMeanValueMaterialTangentMatrix(ElementalVariables &rElementalVariables, double &MeanValue, const ShapeFunctionDerivativesType &rShapeDeriv, const double secondLame, double &bulkModulus, const double Weight, double &MeanValueMass, const double TimeStep)
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:349
virtual void ComputeBulkReductionCoefficient(MatrixType MassMatrix, MatrixType StiffnessMatrix, double &meanValueStiff, double &bulkCoefficient, double timeStep)
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:358
double GetThetaContinuity() override
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:343
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Calculate the element's local contribution to the system for the current step.
Definition: two_step_updated_lagrangian_V_P_implicit_element.cpp:32
TwoStepUpdatedLagrangianVPImplicitElement(IndexType NewId, GeometryType::Pointer pGeometry)
Constructor using a geometry object.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:151
void UpdateCauchyStress(unsigned int g, const ProcessInfo &rCurrentProcessInfo) override
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:251
BaseType::ElementalVariables ElementalVariables
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:115
double GetThetaMomentum() override
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:337
virtual void ComputeStabLaplacianMatrix(MatrixType &StabLaplacianMatrix, const ShapeFunctionDerivativesType &rShapeDeriv, const double Weight)
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:385
void InitializeSolutionStep(const ProcessInfo &rCurrentProcessInfo) override
Initializes the element and all geometric information required for the problem.
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:203
TwoStepUpdatedLagrangianVPImplicitElement(TwoStepUpdatedLagrangianVPImplicitElement const &rOther)
copy constructor
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:167
virtual void AddStabilizationNodalTermsRHS(VectorType &rRightHandSideVector, const double Tau, const double Density, const double Weight, const ShapeFunctionDerivativesType &rDN_DX, const SizeType i)
Definition: two_step_updated_lagrangian_V_P_implicit_element.h:416
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: updated_lagrangian_element.hpp:177
void CalculateLocalMomentumEquations(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: two_step_updated_lagrangian_V_P_implicit_element.cpp:61
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
def load(f)
Definition: ode_solve.py:307
integer i
Definition: TensorModule.f:17
Definition: updated_lagrangian_element.h:74