54 #if !defined(KRATOS_PRESSURE_CALCULATE_PROCESS_INCLUDED )
55 #define KRATOS_PRESSURE_CALCULATE_PROCESS_INCLUDED
73 #include "utilities/geometry_utilities.h"
156 for(ModelPart::NodesContainerType::iterator in = mr_model_part.
NodesBegin() ;
157 in != mr_model_part.
NodesEnd() ; ++in)
159 in->FastGetSolutionStepValue(PRESSURE)= 0.0;
162 for(ModelPart::ElementsContainerType::iterator
im = mr_model_part.
ElementsBegin() ;
165 im->Calculate(PRESSURE,
dummy,proc_info);
205 for(ModelPart::NodesContainerType::iterator in = mr_model_part.
NodesBegin() ;
206 in != mr_model_part.
NodesEnd() ; ++in)
208 if( (in->GetValue(NEIGHBOUR_ELEMENTS)).size() != 0)
210 const double& ar=in->FastGetSolutionStepValue(NODAL_AREA);
211 in->FastGetSolutionStepValue(PRESSURE)/=ar;
217 in->FastGetSolutionStepValue(PRESSURE)=0.0;
240 std::string
Info()
const override
242 return "PressureCalculateProcess";
248 rOStream <<
"PressureCalculateProcess";
311 unsigned int mdomain_size;
369 rOStream << std::endl;
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
Short class definition.
Definition: pressure_calculate_process.h:109
KRATOS_CLASS_POINTER_DEFINITION(PressureCalculateProcess)
Pointer definition of PressureCalculateProcess.
~PressureCalculateProcess() override
Destructor.
Definition: pressure_calculate_process.h:128
PressureCalculateProcess(ModelPart &model_part, unsigned int domain_size)
Default constructor.
Definition: pressure_calculate_process.h:122
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: pressure_calculate_process.h:252
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: pressure_calculate_process.h:246
void operator()()
Definition: pressure_calculate_process.h:137
std::string Info() const override
Turn back information as a string.
Definition: pressure_calculate_process.h:240
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: pressure_calculate_process.h:147
The base class for all processes in Kratos.
Definition: process.h:49
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_WATCH(variable)
Definition: define.h:806
#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
int domain_size
Definition: face_heat.py:4
model_part
Definition: face_heat.py:14
dummy
Definition: script.py:194