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.
inlet_mesher_process.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosPfemApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: July 2018 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_INLET_MESHER_PROCESS_H_INCLUDED )
11 #define KRATOS_INLET_MESHER_PROCESS_H_INCLUDED
12 
13 
14 // External includes
15 
16 // System includes
17 
18 // Project includes
19 #include "includes/model_part.h"
22 
23 
24 namespace Kratos
25 {
26 
29 
31 
35  : public MesherProcess
36 {
37  public:
40 
43 
48 
52 
55  MesherUtilities::MeshingParameters& rRemeshingParameters,
56  int EchoLevel)
57  : mrModelPart(rModelPart),
58  mrRemesh(rRemeshingParameters)
59  {
60  mEchoLevel = EchoLevel;
61  }
62 
63 
65  virtual ~InletMesherProcess() {}
66 
67 
71 
73  void operator()()
74  {
75  Execute();
76  }
77 
78 
82 
84  void Execute() override
85  {
87 
88  if( mEchoLevel > 1 )
89  std::cout<<" [ INLET MANAGEMENT PROCESS: "<<std::endl;
90 
91  if( mrModelPart.Name() != mrRemesh.SubModelPartName )
92  std::cout<<" ModelPart Supplied do not corresponds to the Meshing Domain: ("<<mrModelPart.Name()<<" != "<<mrRemesh.SubModelPartName<<")"<<std::endl;
93 
94  this->SetInletLayer();
95 
96  if( mEchoLevel > 1 )
97  std::cout<<" INLET MANAGEMENT PROCESS ]; "<<std::endl;
98 
99  KRATOS_CATCH(" ")
100  }
101 
102 
106 
107 
111 
112 
116 
118  std::string Info() const override
119  {
120  return "InletMesherProcess";
121  }
122 
124  void PrintInfo(std::ostream& rOStream) const override
125  {
126  rOStream << "InletMesherProcess";
127  }
128 
130  void PrintData(std::ostream& rOStream) const override
131  {
132  }
133 
134 
138 
140 
141 
142  private:
145 
149  ModelPart& mrModelPart;
150 
152 
153  int mEchoLevel;
154 
158 
159 
163 
164 
165  void SetInletLayer()
166 
167  {
168  KRATOS_TRY
169 
170  double critical_distance = 4.0*mrRemesh.Refine->CriticalRadius;
171 
172  unsigned int NodeId = MesherUtilities::GetMaxNodeId( mrModelPart.GetParentModelPart() ) + 1;
173 
175 
176  for(ModelPart::NodesContainerType::iterator i_node = mrModelPart.NodesBegin() ; i_node != mrModelPart.NodesEnd() ; ++i_node)
177  {
178  if(i_node->Is(INLET) ){
179 
180  double distance = norm_2(i_node->FastGetSolutionStepValue(DISPLACEMENT));
181 
182  if( distance > critical_distance ){
183 
184  // add new inlet node
185  NodeType::Pointer pnode = i_node->Clone();
186 
187  pnode->SetId(NodeId);
188  pnode->Coordinates() = pnode->GetInitialPosition();
189 
190  pnode->Set(INLET,true);
191  // release old inlet node
192  i_node->Set(INLET,false);
193 
194  Node::DofsContainerType& Dofs = i_node->GetDofs();
195  //free dofs
196  for(Node::DofsContainerType::iterator i_dof = Dofs.begin(); i_dof != Dofs.end(); ++i_dof)
197  {
198  (*i_dof)->FreeDof();
199  }
200 
201  noalias(pnode->FastGetSolutionStepValue(DISPLACEMENT)) = ZeroVector(3);
202  noalias(pnode->FastGetSolutionStepValue(DISPLACEMENT,1)) = ZeroVector(3);
203 
204  noalias(i_node->FastGetSolutionStepValue(VELOCITY,1)) = i_node->FastGetSolutionStepValue(VELOCITY);
205  noalias(i_node->FastGetSolutionStepValue(ACCELERATION)) = ZeroVector(3);
206  noalias(i_node->FastGetSolutionStepValue(ACCELERATION,1)) = ZeroVector(3);
207 
208  InletNodes.push_back(pnode);
209 
210  ++NodeId;
211 
212  }
213 
214  }
215 
216  }
217 
218  if( InletNodes.size() !=0 )
219  mrModelPart.AddNodes(InletNodes.begin(), InletNodes.end());
220 
221  KRATOS_CATCH("")
222  }
223 
227 
231 
235 
236 
238  InletMesherProcess& operator=(InletMesherProcess const& rOther);
239 
240 
242 
243 
245  //Process(Process const& rOther);
246 
247 
249 
250 }; // Class Process
251 
253 
256 
257 
261 
262 
264 inline std::istream& operator >> (std::istream& rIStream,
265  InletMesherProcess& rThis);
266 
268 inline std::ostream& operator << (std::ostream& rOStream,
269  const InletMesherProcess& rThis)
270 {
271  rThis.PrintInfo(rOStream);
272  rOStream << std::endl;
273  rThis.PrintData(rOStream);
274 
275  return rOStream;
276 }
278 
279 
280 } // namespace Kratos.
281 
282 #endif // KRATOS_INLET_MESHER_PROCESS_H_INCLUDED defined
Base class for all Conditions.
Definition: condition.h:59
Geometry base class.
Definition: geometry.h:71
Insert a new Inlet Layer when the previous one has gone away.
Definition: inlet_mesher_process.hpp:36
virtual ~InletMesherProcess()
Destructor.
Definition: inlet_mesher_process.hpp:65
void PrintData(std::ostream &rOStream) const override
Print object's data.s.
Definition: inlet_mesher_process.hpp:130
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: inlet_mesher_process.hpp:124
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: inlet_mesher_process.hpp:84
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: inlet_mesher_process.hpp:73
std::string Info() const override
Turn back information as a string.
Definition: inlet_mesher_process.hpp:118
InletMesherProcess(ModelPart &rModelPart, MesherUtilities::MeshingParameters &rRemeshingParameters, int EchoLevel)
Default constructor.
Definition: inlet_mesher_process.hpp:54
ModelPart::PropertiesType PropertiesType
Definition: inlet_mesher_process.hpp:46
KRATOS_CLASS_POINTER_DEFINITION(InletMesherProcess)
Pointer definition of Process.
ModelPart::NodeType NodeType
Definition: inlet_mesher_process.hpp:44
ConditionType::GeometryType GeometryType
Definition: inlet_mesher_process.hpp:47
ModelPart::ConditionType ConditionType
Definition: inlet_mesher_process.hpp:45
The base class for processes passed to the solution scheme.
Definition: mesher_process.hpp:37
static unsigned int GetMaxNodeId(ModelPart &rModelPart)
Definition: mesher_utilities.hpp:1408
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
std::string & Name()
Definition: model_part.h:1811
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
This class defines the node.
Definition: node.h:65
std::vector< std::unique_ptr< Dof< double > >> DofsContainerType
The DoF container type definition.
Definition: node.h:92
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
static int EchoLevel
Definition: co_sim_EMPIRE_API.h:42
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
TExpressionType::data_type norm_2(AMatrix::MatrixExpression< TExpressionType, TCategory > const &TheExpression)
Definition: amatrix_interface.h:625
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
Definition: mesher_utilities.hpp:631
RefiningParameters::Pointer Refine
Definition: mesher_utilities.hpp:684
std::string SubModelPartName
Definition: mesher_utilities.hpp:642