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.
mark_fluid_process.h
Go to the documentation of this file.
1 /*
2 ==============================================================================
3 KratosULFApplication
4 A library based on:
5 Kratos
6 A General Purpose Software for Multi-Physics Finite Element Analysis
7 Version 1.0 (Released on march 05, 2007).
8 
9 Copyright 2007
10 Pooyan Dadvand, Riccardo Rossi, Pawel Ryzhakov
11 pooyan@cimne.upc.edu
12 rrossi@cimne.upc.edu
13 - CIMNE (International Center for Numerical Methods in Engineering),
14 Gran Capita' s/n, 08034 Barcelona, Spain
15 
16 
17 Permission is hereby granted, free of charge, to any person obtaining
18 a copy of this software and associated documentation files (the
19 "Software"), to deal in the Software without restriction, including
20 without limitation the rights to use, copy, modify, merge, publish,
21 distribute, sublicense and/or sell copies of the Software, and to
22 permit persons to whom the Software is furnished to do so, subject to
23 the following condition:
24 
25 Distribution of this code for any commercial purpose is permissible
26 ONLY BY DIRECT ARRANGEMENT WITH THE COPYRIGHT OWNERS.
27 
28 The above copyright notice and this permission notice shall be
29 included in all copies or substantial portions of the Software.
30 
31 THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
32 EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
33 MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.
34 IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY
35 CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT,
36 TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE
37 SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
38 
39 ==============================================================================
40 */
41 
42 
43 //
44 // Project Name: Kratos
45 // Last Modified by: $Author: rrossi $
46 // Date: $Date: 2007-05-16 13:59:01 $
47 // Revision: $Revision: 1.3 $
48 //
49 // this process puts a IS_FLUID flag on the nodes of fluid elements
50 
51 #if !defined(KRATOS_MARK_FLUID_PROCESS_INCLUDED )
52 #define KRATOS_MARK_FLUID_PROCESS_INCLUDED
53 
54 
55 
56 // System includes
57 #include <string>
58 #include <iostream>
59 #include <algorithm>
60 
61 // External includes
62 
63 
64 // Project includes
65 #include "includes/define.h"
66 #include "processes/process.h"
67 #include "includes/node.h"
68 #include "includes/element.h"
69 #include "includes/model_part.h"
70 //#include "custom_utilities/geometry_utilities2D.h"
72 
73 
74 namespace Kratos
75 {
76 
79 
83 
84 
88 
92 
96 
98 
104 class MarkFluidProcess
105  : public Process
106 {
107 public:
110 
113 
117 
120  : mr_model_part(model_part)
121  {
122  }
123 
125  ~MarkFluidProcess() override
126  {
127  }
128 
129 
133 
134  void operator()()
135  {
136  Execute();
137  }
138 
139 
143 
144  void Execute() override
145  {
146  KRATOS_TRY
147 
148  for(ModelPart::NodesContainerType::const_iterator in = mr_model_part.NodesBegin(); in!=mr_model_part.NodesEnd(); in++)
149  {
150  in->FastGetSolutionStepValue(IS_FLUID) = 0.0;
151  }
152 
153  ProcessInfo& proc_info = mr_model_part.GetProcessInfo();
154  double dummy;
155  for(ModelPart::ElementsContainerType::iterator im = mr_model_part.ElementsBegin() ;
156  im != mr_model_part.ElementsEnd() ; ++im)
157  {
158  im->Calculate(IS_FLUID,dummy,proc_info);
159  }
160  /* if( (in->GetValue(NEIGHBOUR_ELEMENTS)).size() != 0)
161  {
162  in->FastGetSolutionStepValue(IS_FLUID) = 1.0;
163 
164  }
165  else
166  in->FastGetSolutionStepValue(IS_FLUID) = 0.0;
167  }*/
168  KRATOS_CATCH("")
169  }
170 
171 
175 
176 
180 
181 
185 
187  std::string Info() const override
188  {
189  return "MarkFluidProcess";
190  }
191 
193  void PrintInfo(std::ostream& rOStream) const override
194  {
195  rOStream << "MarkFluidProcess";
196  }
197 
199  void PrintData(std::ostream& rOStream) const override
200  {
201  }
202 
203 
207 
208 
210 
211 protected:
214 
215 
219 
220 
224 
225 
229 
230 
234 
235 
239 
240 
244 
245 
247 
248 private:
251 
252 
256  ModelPart& mr_model_part;
257 
261 
262 
266 
267 
271 
272 
276 
277 
281 
283 // MarkFluidProcess& operator=(MarkFluidProcess const& rOther);
284 
286 // MarkFluidProcess(MarkFluidProcess const& rOther);
287 
288 
290 
291 }; // Class MarkFluidProcess
292 
294 
297 
298 
302 
303 
305 inline std::istream& operator >> (std::istream& rIStream,
306  MarkFluidProcess& rThis);
307 
309 inline std::ostream& operator << (std::ostream& rOStream,
310  const MarkFluidProcess& rThis)
311 {
312  rThis.PrintInfo(rOStream);
313  rOStream << std::endl;
314  rThis.PrintData(rOStream);
315 
316  return rOStream;
317 }
319 
320 
321 } // namespace Kratos.
322 
323 #endif // KRATOS_MARK_FLUID_PROCESS_INCLUDED defined
324 
325 
Short class definition.
Definition: mark_fluid_process.h:74
MarkFluidProcess(ModelPart &model_part)
Default constructor.
Definition: mark_fluid_process.h:119
void operator()()
Definition: mark_fluid_process.h:134
KRATOS_CLASS_POINTER_DEFINITION(MarkFluidProcess)
Pointer definition of PushStructureProcess.
~MarkFluidProcess() override
Destructor.
Definition: mark_fluid_process.h:125
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: mark_fluid_process.h:199
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: mark_fluid_process.h:193
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: mark_fluid_process.h:144
std::string Info() const override
Turn back information as a string.
Definition: mark_fluid_process.h:187
virtual void Execute() override
Execute method is used to execute the Process algorithms.
Definition: mark_fluid_process.h:112
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
ElementIterator ElementsBegin(IndexType ThisIndex=0)
Definition: model_part.h:1169
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
ElementIterator ElementsEnd(IndexType ThisIndex=0)
Definition: model_part.h:1179
NodeIterator NodesEnd(IndexType ThisIndex=0)
Definition: model_part.h:497
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
im
Definition: GenerateCN.py:100
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
dummy
Definition: script.py:194