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.
dof_updater.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 #if !defined(KRATOS_DOF_UPDATER_H_INCLUDED )
14 #define KRATOS_DOF_UPDATER_H_INCLUDED
15 
16 // Project includes
17 #include "includes/define.h"
18 #include "includes/model_part.h"
20 
21 namespace Kratos
22 {
25 
28 
30 
38 template< class TSparseSpace >
40 {
41 public:
44 
47 
50 
52 
56 
59 
61  DofUpdater(DofUpdater const& rOther) = delete;
62 
64  virtual ~DofUpdater(){}
65 
67  DofUpdater& operator=(DofUpdater const& rOther) = delete;
68 
72 
74 
79  virtual typename DofUpdater::UniquePointer Create() const
80  {
81  return Kratos::make_unique<DofUpdater>();
82  }
83 
85 
89  virtual void Initialize(
90  const DofsArrayType& rDofSet,
91  const SystemVectorType& rDx)
92  {}
93 
95 
97  virtual void Clear() {}
98 
100 
106  virtual void UpdateDofs(
107  DofsArrayType& rDofSet,
108  const SystemVectorType& rDx)
109  {
111  rDofSet,
112  [&rDx](DofType& rDof)
113  {
114  if (rDof.IsFree()) {
115  rDof.GetSolutionStepValue() += TSparseSpace::GetValue(rDx, rDof.EquationId());
116  }
117  }
118  );
119  }
120 
122 
128  virtual void AssignDofs(DofsArrayType& rDofSet, const SystemVectorType& rX)
129  {
131  rDofSet,
132  [&rX](DofType& rDof)
133  {
134  if (rDof.IsFree()) {
135  rDof.GetSolutionStepValue() = TSparseSpace::GetValue(rX, rDof.EquationId());
136  }
137  }
138  );
139  }
140 
144 
146  virtual std::string Info() const
147  {
148  std::stringstream buffer;
149  buffer << "DofUpdater" ;
150  return buffer.str();
151  }
152 
154  virtual void PrintInfo(std::ostream& rOStream) const
155  {
156  rOStream << this->Info() << std::endl;
157  }
158 
160  virtual void PrintData(std::ostream& rOStream) const
161  {
162  rOStream << this->Info() << std::endl;
163  }
164 
166 
167 }; // Class DofUpdater
168 
169 
173 
175 template< class TSparseSpace >
176 inline std::istream& operator >> (
177  std::istream& rIStream,
179 {
180  return rIStream;
181 }
182 
184 template< class TSparseSpace >
185 inline std::ostream& operator << (
186  std::ostream& rOStream,
187  const DofUpdater<TSparseSpace>& rThis)
188 {
189  rThis.PrintInfo(rOStream);
190  rOStream << std::endl;
191  rThis.PrintData(rOStream);
192 
193  return rOStream;
194 }
195 
197 
199 
200 } // namespace Kratos.
201 
202 #endif // KRATOS_DOF_UPDATER_H_INCLUDED defined
Dof represents a degree of freedom (DoF).
Definition: dof.h:86
bool IsFree() const
Definition: dof.h:382
Utility class to update the values of degree of freedom (Dof) variables after solving the system.
Definition: dof_updater.h:40
virtual void Clear()
Free internal storage to reset the instance and/or optimize memory consumption.
Definition: dof_updater.h:97
virtual void Initialize(const DofsArrayType &rDofSet, const SystemVectorType &rDx)
Initialize the DofUpdater in preparation for a subsequent UpdateDofs call.
Definition: dof_updater.h:89
DofUpdater(DofUpdater const &rOther)=delete
Deleted copy constructor.
DofUpdater()
Default constructor.
Definition: dof_updater.h:58
virtual void UpdateDofs(DofsArrayType &rDofSet, const SystemVectorType &rDx)
Calculate new values for the problem's degrees of freedom using the update vector rDx.
Definition: dof_updater.h:106
typename TSparseSpace::VectorType SystemVectorType
Definition: dof_updater.h:51
virtual DofUpdater::UniquePointer Create() const
Create a new instance of this class.
Definition: dof_updater.h:79
virtual void PrintData(std::ostream &rOStream) const
Print object's data.
Definition: dof_updater.h:160
virtual void AssignDofs(DofsArrayType &rDofSet, const SystemVectorType &rX)
Assign new values for the problem's degrees of freedom using the vector rX.
Definition: dof_updater.h:128
virtual void PrintInfo(std::ostream &rOStream) const
Print information about this object.
Definition: dof_updater.h:154
KRATOS_CLASS_POINTER_DEFINITION(DofUpdater)
Pointer definition of DofUpdater.
virtual ~DofUpdater()
Destructor.
Definition: dof_updater.h:64
DofUpdater & operator=(DofUpdater const &rOther)=delete
Deleted assignment operator.
virtual std::string Info() const
Turn back information as a string.
Definition: dof_updater.h:146
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
Vector VectorType
Definition: geometrical_transformation_utilities.h:56
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
void block_for_each(TIterator itBegin, TIterator itEnd, TFunction &&rFunction)
Execute a functor on all items of a range in parallel.
Definition: parallel_utilities.h:299
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