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.
model_volume_calculation_process.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosDelaunayMeshingApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: April 2018 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_VOLUME_CALCULATION_PROCESS_H_INCLUDED )
11 #define KRATOS_VOLUME_CALCULATION_PROCESS_H_INCLUDED
12 
13 
14 // External includes
15 
16 // System includes
17 
18 // Project includes
19 #include "includes/model_part.h"
20 
22 //Data:
23 //StepData:
24 //Flags: (checked)
25 // (set)
26 // (modified)
27 // (reset)
28 
29 
30 namespace Kratos
31 {
32 
35 
37  : public Process
38 {
39 public:
42 
45 
52 
54  ModelVolumeCalculationProcess(ModelPart& rModelPart, bool Axisymmetric, int EchoLevel)
55  : mrModelPart(rModelPart)
56  {
57  mAxisymmetric = Axisymmetric;
58  mEchoLevel = EchoLevel;
59  }
60 
62  : mrModelPart(rModelPart)
63  {
64  mAxisymmetric = false;
65  }
66 
69 
70 
74 
75 
79 
80 
83  {
84 
86 
87  mModelVolume = ComputeModelVolume(mMeshVolume);
88 
89  if( mEchoLevel > 0 )
90  std::cout<<" [ Model Volume : "<<mModelVolume<<" ]"<<std::endl;
91 
92  KRATOS_CATCH( "" )
93 
94  }
95 
98  {
100 
101  std::vector<double> MeshVolume;
102  double ModelVolume = ComputeModelVolume(MeshVolume);
103 
104  if( mEchoLevel > 0 )
105  std::cout<<" [ Model Volume : "<<ModelVolume<<" ] [ Step increment : "<<ModelVolume-mModelVolume<<" ] "<<std::endl;
106 
107 
108  KRATOS_CATCH( "" )
109  }
110 
114 
115 
119 
120 
124 
126  std::string Info() const override
127  {
128  return "ModelVolumeCalculationProcess";
129  }
130 
132  void PrintInfo(std::ostream& rOStream) const override
133  {
134  rOStream << "ModelVolumeCalculationProcess";
135  }
136 
138  void PrintData(std::ostream& rOStream) const override
139  {
140  }
141 
142 
146 
147 
149 
150 
151 private:
154 
158  ModelPart& mrModelPart;
159 
160  bool mAxisymmetric;
161 
162  double mModelVolume;
163 
164  std::vector<double> mMeshVolume;
165 
166  int mEchoLevel;
167 
171 
174 
175 
177 
178  double ComputeModelVolume(std::vector<double> & rMeshVolume)
179  {
180 
181  KRATOS_TRY
182 
183  unsigned int start=0;
184  unsigned int NumberOfMeshes=mrModelPart.NumberOfMeshes();
185  if(NumberOfMeshes>1)
186  start=1;
187 
188 
189  rMeshVolume.resize(NumberOfMeshes+start);
190  std::fill( rMeshVolume.begin(), rMeshVolume.end(), 0.0 );
191 
192  double ModelVolume = 0;
193 
194  if( mAxisymmetric ){
195 
196  //std::cout<<" AXISYMMETRIC MODEL "<<std::endl;
197 
198  double radius = 0;
199  double two_pi = 6.28318530717958647693;
200 
201  //By the way: set meshes options from bools
202  for(unsigned int MeshId=start; MeshId<NumberOfMeshes; ++MeshId)
203  {
204  for(ModelPart::ElementsContainerType::const_iterator ie = mrModelPart.ElementsBegin(MeshId); ie != mrModelPart.ElementsEnd(MeshId); ++ie)
205  {
206 
207  const unsigned int dimension = ie->GetGeometry().WorkingSpaceDimension();
208 
209  if( dimension > 2 )
210  std::cout<<" Axisymmetric problem with dimension: "<<dimension<<std::endl;
211 
212  radius = 0;
213  for( unsigned int i=0; i<ie->GetGeometry().size(); i++ )
214  radius += ie->GetGeometry()[i].X();
215 
216  radius/=double(ie->GetGeometry().size());
217 
218  rMeshVolume[MeshId] += ie->GetGeometry().Area() * two_pi * radius ;
219 
220 
221  }
222 
223  ModelVolume += rMeshVolume[MeshId];
224  }
225 
226  }
227  else{
228 
229  //By the way: set meshes options from bools
230  for(unsigned int MeshId=start; MeshId<NumberOfMeshes; ++MeshId)
231  {
232  for(ModelPart::ElementsContainerType::const_iterator ie = mrModelPart.ElementsBegin(MeshId); ie != mrModelPart.ElementsEnd(MeshId); ++ie)
233  {
234  const unsigned int dimension = ie->GetGeometry().WorkingSpaceDimension();
235  if( dimension == 2){
236  rMeshVolume[MeshId] += ie->GetGeometry().Area() * ie->GetProperties()[THICKNESS];
237  }
238  else if( dimension == 3){
239  rMeshVolume[MeshId] += ie->GetGeometry().Volume();
240  }
241  else{
242  //do nothing.
243  }
244 
245  }
246 
247  ModelVolume += rMeshVolume[MeshId];
248 
249  }
250 
251  }
252 
253  return ModelVolume;
254 
255  KRATOS_CATCH( "" )
256 
257  }
258 
259 
261  //Process(Process const& rOther);
262 
263 
265 
266 }; // Class Process
267 
269 
272 
273 
277 
278 
280 inline std::istream& operator >> (std::istream& rIStream,
282 
284 inline std::ostream& operator << (std::ostream& rOStream,
285  const ModelVolumeCalculationProcess& rThis)
286 {
287  rThis.PrintInfo(rOStream);
288  rOStream << std::endl;
289  rThis.PrintData(rOStream);
290 
291  return rOStream;
292 }
294 
295 
296 } // namespace Kratos.
297 
298 #endif // KRATOS_MODEL_VOLUME_CALCULATION_PROCESS_H_INCLUDED defined
Base class for all Conditions.
Definition: condition.h:59
Geometry base class.
Definition: geometry.h:71
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
Definition: model_volume_calculation_process.hpp:38
ConditionType::GeometryType GeometryType
Definition: model_volume_calculation_process.hpp:48
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: model_volume_calculation_process.hpp:82
ModelPart::PropertiesType PropertiesType
Definition: model_volume_calculation_process.hpp:47
ModelVolumeCalculationProcess(ModelPart &rModelPart, bool Axisymmetric, int EchoLevel)
Default constructor.
Definition: model_volume_calculation_process.hpp:54
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: model_volume_calculation_process.hpp:132
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: model_volume_calculation_process.hpp:97
virtual ~ModelVolumeCalculationProcess()
Destructor.
Definition: model_volume_calculation_process.hpp:68
ModelVolumeCalculationProcess(ModelPart &rModelPart)
Definition: model_volume_calculation_process.hpp:61
std::string Info() const override
Turn back information as a string.
Definition: model_volume_calculation_process.hpp:126
ModelPart::ConditionType ConditionType
Definition: model_volume_calculation_process.hpp:46
KRATOS_CLASS_POINTER_DEFINITION(ModelVolumeCalculationProcess)
Pointer definition of Process.
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: model_volume_calculation_process.hpp:138
The base class for all processes in Kratos.
Definition: process.h:49
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
start
Definition: DEM_benchmarks.py:17
static int EchoLevel
Definition: co_sim_EMPIRE_API.h:42
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
TABLE_NUMBER_ANGULAR_VELOCITY TABLE_NUMBER_MOMENT I33 BEAM_INERTIA_ROT_UNIT_LENGHT_Y KRATOS_DEFINE_APPLICATION_VARIABLE(DEM_APPLICATION, double, BEAM_INERTIA_ROT_UNIT_LENGHT_Z) typedef std double
Definition: DEM_application_variables.h:182
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
int dimension
Definition: isotropic_damage_automatic_differentiation.py:123
float radius
Definition: mesh_to_mdpa_converter.py:18
integer i
Definition: TensorModule.f:17