10 #if !defined(KRATOS_VOLUME_CALCULATION_PROCESS_H_INCLUDED )
11 #define KRATOS_VOLUME_CALCULATION_PROCESS_H_INCLUDED
55 : mrModelPart(rModelPart)
57 mAxisymmetric = Axisymmetric;
62 : mrModelPart(rModelPart)
64 mAxisymmetric =
false;
87 mModelVolume = ComputeModelVolume(mMeshVolume);
90 std::cout<<
" [ Model Volume : "<<mModelVolume<<
" ]"<<std::endl;
101 std::vector<double> MeshVolume;
102 double ModelVolume = ComputeModelVolume(MeshVolume);
105 std::cout<<
" [ Model Volume : "<<ModelVolume<<
" ] [ Step increment : "<<ModelVolume-mModelVolume<<
" ] "<<std::endl;
126 std::string
Info()
const override
128 return "ModelVolumeCalculationProcess";
134 rOStream <<
"ModelVolumeCalculationProcess";
164 std::vector<double> mMeshVolume;
178 double ComputeModelVolume(std::vector<double> & rMeshVolume)
183 unsigned int start=0;
184 unsigned int NumberOfMeshes=mrModelPart.NumberOfMeshes();
189 rMeshVolume.resize(NumberOfMeshes+
start);
190 std::fill( rMeshVolume.begin(), rMeshVolume.end(), 0.0 );
192 double ModelVolume = 0;
199 double two_pi = 6.28318530717958647693;
202 for(
unsigned int MeshId=
start; MeshId<NumberOfMeshes; ++MeshId)
204 for(ModelPart::ElementsContainerType::const_iterator ie = mrModelPart.ElementsBegin(MeshId); ie != mrModelPart.ElementsEnd(MeshId); ++ie)
207 const unsigned int dimension = ie->GetGeometry().WorkingSpaceDimension();
210 std::cout<<
" Axisymmetric problem with dimension: "<<
dimension<<std::endl;
213 for(
unsigned int i=0;
i<ie->GetGeometry().size();
i++ )
214 radius += ie->GetGeometry()[
i].X();
218 rMeshVolume[MeshId] += ie->GetGeometry().Area() * two_pi *
radius ;
223 ModelVolume += rMeshVolume[MeshId];
230 for(
unsigned int MeshId=
start; MeshId<NumberOfMeshes; ++MeshId)
232 for(ModelPart::ElementsContainerType::const_iterator ie = mrModelPart.ElementsBegin(MeshId); ie != mrModelPart.ElementsEnd(MeshId); ++ie)
234 const unsigned int dimension = ie->GetGeometry().WorkingSpaceDimension();
236 rMeshVolume[MeshId] += ie->GetGeometry().Area() * ie->GetProperties()[THICKNESS];
239 rMeshVolume[MeshId] += ie->GetGeometry().Volume();
247 ModelVolume += rMeshVolume[MeshId];
288 rOStream << std::endl;
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