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.
ulf_apply_bc_process.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: Author Julio Marti.
11 //
12 
13 
14 //
15 // Project Name: Kratos
16 // Last Modified by: $Author: anonymous $
17 // Date: $Date: 2008-10-23 12:50:01 $
18 // Revision: $Revision: 1.4 $
19 //
20 // this process save structural elements in a separate list
21 
22 #if !defined(KRATOS_ULF_APPLY_BC_PROCESS_INCLUDED )
23 #define KRATOS_ULF_APPLY_BC_PROCESS_INCLUDED
24 
25 
26 
27 // System includes
28 #include <string>
29 #include <iostream>
30 #include <algorithm>
31 
32 // External includes
33 
34 
35 // Project includes
36 #include "includes/define.h"
37 #include "processes/process.h"
38 #include "includes/node.h"
39 #include "includes/element.h"
40 #include "includes/model_part.h"
41 //#include "custom_utilities/geometry_utilities2D.h"
42 //#include "custom_elements/updated_lagrangian_fluid.h"
43 
44 
45 namespace Kratos
46 {
47 
50 
54 
55 
59 
63 
67 
69 
76  : public Process
77  {
78  public:
81 
84 
88 
91  : mr_model_part(model_part)
92  {
93  }
94 
97  {
98  }
99 
100 
104 
105  void operator()()
106  {
107  Execute();
108  }
109 
110 
114 
115  virtual void Execute() override
116  {
117  KRATOS_TRY
118  for(ModelPart::NodesContainerType::const_iterator in = mr_model_part.NodesBegin(); in!=mr_model_part.NodesEnd(); in++)
119  {
120  //marking wet nodes
121  /* if(in->FastGetSolutionStepValue(IS_STRUCTURE) )
122  if( (in->GetValue(NEIGHBOUR_ELEMENTS)).size() != 0)
123  in->FastGetSolutionStepValue(IS_FLUID) = 1.0;
124  else //it is not anymore of fluid
125  in->FastGetSolutionStepValue(IS_FLUID) = 0.0;
126  //marking as free surface the lonely nodes
127  else
128  */ if( (in->GetValue(NEIGHBOUR_ELEMENTS)).size() == 0)
129  in->FastGetSolutionStepValue(IS_BOUNDARY) = 1.0;
130  }
131 
132  //identify the free surface
133  for(ModelPart::NodeIterator i = mr_model_part.NodesBegin() ;
134  i != mr_model_part.NodesEnd() ; ++i)
135  {
136  //reset the free surface
137  i->FastGetSolutionStepValue(IS_FREE_SURFACE) = 0;
138  //identify the free surface and fix the pressure accordingly
139  if( i->FastGetSolutionStepValue(IS_BOUNDARY) != 0
140  &&
141  i->FastGetSolutionStepValue(IS_STRUCTURE) == 0)
142  {
143  i->FastGetSolutionStepValue(IS_FREE_SURFACE) = 1;
144  }
145  }
146 
147 
148  KRATOS_CATCH("")
149  }
150 
151 
155 
156 
160 
161 
165 
167  virtual std::string Info() const override
168  {
169  return "Pfem2ApplyBCProcess";
170  }
171 
173  virtual void PrintInfo(std::ostream& rOStream) const override
174  {
175  rOStream << "Pfem2ApplyBCProcess";
176  }
177 
179  virtual void PrintData(std::ostream& rOStream) const override
180  {
181  }
182 
183 
187 
188 
190 
191  protected:
194 
195 
199 
200 
204 
205 
209 
210 
214 
215 
219 
220 
224 
225 
227 
228  private:
231 
232 
236  ModelPart& mr_model_part;
237 
241 
242 
246 
247 
251 
252 
256 
257 
261 
263  // UlfApplyBCProcess& operator=(UlfApplyBCProcess const& rOther);
264 
266  // UlfApplyBCProcess(UlfApplyBCProcess const& rOther);
267 
268 
270 
271 }; // Class UlfApplyBCProcess
272 
274 
277 
278 
282 
283 
285  inline std::istream& operator >> (std::istream& rIStream,
286  Pfem2ApplyBCProcess& rThis);
287 
289  inline std::ostream& operator << (std::ostream& rOStream,
290  const Pfem2ApplyBCProcess& rThis)
291  {
292  rThis.PrintInfo(rOStream);
293  rOStream << std::endl;
294  rThis.PrintData(rOStream);
295 
296  return rOStream;
297  }
299 
300 
301 } // namespace Kratos.
302 
303 #endif // KRATOS_ULF_APPLY_BC_PROCESS_INCLUDED defined
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
MeshType::NodeIterator NodeIterator
Definition: model_part.h:134
NodeIterator NodesEnd(IndexType ThisIndex=0)
Definition: model_part.h:497
Short class definition.
Definition: ulf_apply_bc_process.h:77
Pfem2ApplyBCProcess(ModelPart &model_part)
Default constructor.
Definition: ulf_apply_bc_process.h:90
KRATOS_CLASS_POINTER_DEFINITION(Pfem2ApplyBCProcess)
Pointer definition of PushStructureProcess.
virtual std::string Info() const override
Turn back information as a string.
Definition: ulf_apply_bc_process.h:167
virtual void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: ulf_apply_bc_process.h:173
virtual void Execute() override
Execute method is used to execute the Process algorithms.
Definition: ulf_apply_bc_process.h:115
virtual void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: ulf_apply_bc_process.h:179
virtual ~Pfem2ApplyBCProcess()
Destructor.
Definition: ulf_apply_bc_process.h:96
void operator()()
Definition: ulf_apply_bc_process.h:105
The base class for all processes in Kratos.
Definition: process.h:49
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
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
model_part
Definition: face_heat.py:14
integer i
Definition: TensorModule.f:17