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.
relaxed_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: Suneth Warnakulasuriya
11 //
12 
13 #if !defined(KRATOS_RELAXED_DOF_UPDATER_H_INCLUDED )
14 #define KRATOS_RELAXED_DOF_UPDATER_H_INCLUDED
15 
16 // Project includes
17 #include "includes/define.h"
18 #include "includes/model_part.h"
19 #include "utilities/dof_updater.h"
20 
21 #ifdef KRATOS_USING_MPI // mpi-parallel compilation
22 class Epetra_Import;
23 #endif
24 
25 
26 
27 namespace Kratos
28 {
31 
34 
36 
43 template< class TSparseSpace >
44 class KRATOS_API(RANS_APPLICATION) RelaxedDofUpdater
45  :public DofUpdater<TSparseSpace>
46 {
47 public:
50 
53 
55 
56  using DofType = typename BaseType::DofType;
57 
59 
61 
62 
63 
67 
69  RelaxedDofUpdater(const double RelaxationFactor)
70  : mRelaxationFactor(RelaxationFactor)
71  {}
72 
74  RelaxedDofUpdater(RelaxedDofUpdater const& rOther) = delete;
75 
77  virtual ~RelaxedDofUpdater() = default;
78 
80  RelaxedDofUpdater& operator=(RelaxedDofUpdater const& rOther) = delete;
81 
85 
87 
92  typename BaseType::UniquePointer Create() const override
93  {
94  return Kratos::make_unique<RelaxedDofUpdater>(mRelaxationFactor);
95  }
96 
98 
103  const DofsArrayType& rDofSet,
104  const SystemVectorType& rDx) override;
105 
107 
109  void Clear() override;
110 
112 
120  DofsArrayType& rDofSet,
121  const SystemVectorType& rDx) override;
122 
126 
128  std::string Info() const override;
129 
131  void PrintInfo(std::ostream& rOStream) const override
132  {
133  rOStream << this->Info() << std::endl;
134  }
135 
137  void PrintData(std::ostream& rOStream) const override
138  {
139  rOStream << this->Info() << std::endl;
140  }
141 
143 private:
144 
145  //@name Member Variables
147 
149  bool mImportIsInitialized = false;
150  const double mRelaxationFactor;
151 
152  #ifdef KRATOS_USING_MPI // mpi-parallel compilation
154  std::shared_ptr<Epetra_Import> mpDofImport = nullptr;
155  #endif
156 
158 
159 
160 }; // Class RelaxedDofUpdater
161 
162 
166 
168 template< class TSparseSpace >
169 inline std::istream& operator >> (
170  std::istream& rIStream,
172 {
173  return rIStream;
174 }
175 
177 template< class TSparseSpace >
178 inline std::ostream& operator << (
179  std::ostream& rOStream,
180  const RelaxedDofUpdater<TSparseSpace>& rThis)
181 {
182  rThis.PrintInfo(rOStream);
183  rOStream << std::endl;
184  rThis.PrintData(rOStream);
185 
186  return rOStream;
187 }
188 
190 
192 
193 } // namespace Kratos.
194 
195 #endif // KRATOS_RELAXED_DOF_UPDATER_H_INCLUDED defined
std::string Info() const override
Turn back information as a string.
Definition: periodic_interface_process.hpp:93
Dof represents a degree of freedom (DoF).
Definition: dof.h:86
Utility class to update the values of degree of freedom (Dof) variables after solving the system.
Definition: dof_updater.h:40
typename TSparseSpace::VectorType SystemVectorType
Definition: dof_updater.h:51
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
Utility class to update the values of degree of freedom (Dof) variables after solving the system.
Definition: relaxed_dof_updater.h:46
void UpdateDofs(DofsArrayType &rDofSet, const SystemVectorType &rDx) override
Calculate new values for the problem's degrees of freedom using the update vector rDx.
BaseType::UniquePointer Create() const override
Create a new instance of this class.
Definition: relaxed_dof_updater.h:92
RelaxedDofUpdater(RelaxedDofUpdater const &rOther)=delete
Deleted copy constructor.
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: relaxed_dof_updater.h:137
void Initialize(const DofsArrayType &rDofSet, const SystemVectorType &rDx) override
Initialize the RelaxedDofUpdater in preparation for a subsequent UpdateDofs call.
RelaxedDofUpdater & operator=(RelaxedDofUpdater const &rOther)=delete
Deleted assignment operator.
RelaxedDofUpdater(const double RelaxationFactor)
Default constructor.
Definition: relaxed_dof_updater.h:69
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: relaxed_dof_updater.h:131
std::string Info() const override
Turn back information as a string.
virtual ~RelaxedDofUpdater()=default
Destructor.
void Clear() override
Free internal storage to reset the instance and/or optimize memory consumption.
KRATOS_CLASS_POINTER_DEFINITION(RelaxedDofUpdater)
Pointer definition of RelaxedDofUpdater.
Vector VectorType
Definition: geometrical_transformation_utilities.h:56
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
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