52 #if !defined(KRATOS_SET_H_MAP_PROCESS_INCLUDED )
53 #define KRATOS_SET_H_MAP_PROCESS_INCLUDED
154 double min_dist=h_min;
156 for(ModelPart::NodesContainerType::iterator in = rNodes.begin(); in!=rNodes.end(); in++)
159 dist=in->FastGetSolutionStepValue(DISTANCE);
172 double coef=(h_max-h_min)/(0.75*max_dist-h_min);
173 double c=h_min-coef*h_min;
175 for(ModelPart::NodesContainerType::iterator in = rNodes.begin(); in!=rNodes.end(); in++)
178 dist=in->FastGetSolutionStepValue(DISTANCE);
179 if (dist<0.75*max_dist && dist>=min_dist)
181 if ((coef*
dist+
c)>h_min && (coef*
dist+
c)<h_max)
183 in->FastGetSolutionStepValue(NODAL_H)=coef*
dist+
c;
207 virtual std::string
Info()
const
209 return "SetHMapProcess";
215 rOStream <<
"SetHMapProcess";
334 rOStream << std::endl;
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
The base class for all processes in Kratos.
Definition: process.h:49
virtual void Execute()
Execute method is used to execute the Process algorithms.
Definition: process.h:101
Short class definition.
Definition: set_h_map_process.h:106
virtual void PrintInfo(std::ostream &rOStream) const
Print information about this object.
Definition: set_h_map_process.h:213
virtual ~SetHMapProcess()
Destructor.
Definition: set_h_map_process.h:125
SetHMapProcess(ModelPart &model_part)
Default constructor.
Definition: set_h_map_process.h:119
KRATOS_CLASS_POINTER_DEFINITION(SetHMapProcess)
Pointer definition of SetHMapProcess.
virtual std::string Info() const
Turn back information as a string.
Definition: set_h_map_process.h:207
void operator()()
Definition: set_h_map_process.h:134
virtual void PrintData(std::ostream &rOStream) const
Print object's data.
Definition: set_h_map_process.h:219
void CalculateOptimalH(const double h_min, const double h_max, const double mmax_dist)
Definition: set_h_map_process.h:144
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_WATCH(variable)
Definition: define.h:806
#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
float dist
Definition: edgebased_PureConvection.py:89
model_part
Definition: face_heat.py:14
c
Definition: generate_weakly_compressible_navier_stokes_element.py:108