14 #if !defined(KRATOS_DAM_T_SOL_AIR_HEAT_FLUX_PROCESS)
15 #define KRATOS_DAM_T_SOL_AIR_HEAT_FLUX_PROCESS
49 "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
50 "variable_name": "PLEASE_PRESCRIBE_VARIABLE_NAME",
52 "ambient_temperature" : 0.0,
53 "table_ambient_temperature" : 0,
56 "absorption_index" : 0.0,
57 "total_insolation" : 0.0,
67 rParameters[
"delta_R"];
68 rParameters[
"absorption_index"];
119 #pragma omp parallel for
122 ModelPart::NodesContainerType::iterator it = it_begin +
i;
124 const double temp_current = it->FastGetSolutionStepValue(TEMPERATURE);
125 const double heat_flux =
mH0 * (t_sol_air - temp_current);
126 it->FastGetSolutionStepValue(var) = heat_flux;
158 #pragma omp parallel for
161 ModelPart::NodesContainerType::iterator it = it_begin +
i;
163 const double temp_current = it->FastGetSolutionStepValue(TEMPERATURE);
164 const double heat_flux =
mH0 * (t_sol_air - temp_current);
165 it->FastGetSolutionStepValue(var) = heat_flux;
175 std::string
Info()
const override
177 return "DamTSolAirHeatFluxProcess";
183 rOStream <<
"DamTSolAirHeatFluxProcess";
224 rOStream << std::endl;
Definition: dam_t_sol_air_heat_flux_process.hpp:31
void ExecuteInitializeSolutionStep() override
This function will be executed at every time step BEFORE performing the solve phase.
Definition: dam_t_sol_air_heat_flux_process.hpp:135
int mTableId
Definition: dam_t_sol_air_heat_flux_process.hpp:205
std::string mVariableName
Definition: dam_t_sol_air_heat_flux_process.hpp:196
double mTotalInsolation
Definition: dam_t_sol_air_heat_flux_process.hpp:202
double mAbsorption_index
Definition: dam_t_sol_air_heat_flux_process.hpp:201
ModelPart & mrModelPart
Member Variables.
Definition: dam_t_sol_air_heat_flux_process.hpp:195
Table< double, double > TableType
Definition: dam_t_sol_air_heat_flux_process.hpp:36
std::string Info() const override
Turn back information as a string.
Definition: dam_t_sol_air_heat_flux_process.hpp:175
virtual ~DamTSolAirHeatFluxProcess()
Destructor.
Definition: dam_t_sol_air_heat_flux_process.hpp:93
double mH0
Definition: dam_t_sol_air_heat_flux_process.hpp:197
void ExecuteInitialize() override
This function is designed for being called at the beginning of the computations right after reading t...
Definition: dam_t_sol_air_heat_flux_process.hpp:104
double mEmisivity
Definition: dam_t_sol_air_heat_flux_process.hpp:199
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: dam_t_sol_air_heat_flux_process.hpp:187
KRATOS_CLASS_POINTER_DEFINITION(DamTSolAirHeatFluxProcess)
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: dam_t_sol_air_heat_flux_process.hpp:181
DamTSolAirHeatFluxProcess(ModelPart &rModelPart, Parameters &rParameters)
Constructor.
Definition: dam_t_sol_air_heat_flux_process.hpp:41
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: dam_t_sol_air_heat_flux_process.hpp:97
double mDeltaR
Definition: dam_t_sol_air_heat_flux_process.hpp:200
double mAmbientTemperature
Definition: dam_t_sol_air_heat_flux_process.hpp:198
TableType::Pointer mpTable
Definition: dam_t_sol_air_heat_flux_process.hpp:204
double mTimeUnitConverter
Definition: dam_t_sol_air_heat_flux_process.hpp:203
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
NodesContainerType & Nodes()
Definition: mesh.h:346
NodeIterator NodesBegin()
Definition: mesh.h:326
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
TableType::Pointer pGetTable(IndexType TableId)
Definition: model_part.h:595
MeshType & GetMesh(IndexType ThisIndex=0)
Definition: model_part.h:1791
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
double GetDouble() const
This method returns the double contained in the current Parameter.
Definition: kratos_parameters.cpp:657
int GetInt() const
This method returns the integer contained in the current Parameter.
Definition: kratos_parameters.cpp:666
void ValidateAndAssignDefaults(const Parameters &rDefaultParameters)
This function is designed to verify that the parameters under testing match the form prescribed by th...
Definition: kratos_parameters.cpp:1306
std::string GetString() const
This method returns the string contained in the current Parameter.
Definition: kratos_parameters.cpp:684
size_type size() const
Returns the number of elements in the container.
Definition: pointer_vector_set.h:502
The base class for all processes in Kratos.
Definition: process.h:49
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#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
time
Definition: face_heat.py:85
int nnodes
Definition: sensitivityMatrix.py:24
integer i
Definition: TensorModule.f:17