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.
dam_temperature_by_device_process.hpp
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Lorenzo Gracia
11 //
12 //
13 
14 #if !defined(KRATOS_DAM_TEMPERATURE_BY_DEVICE_PROCESS)
15 #define KRATOS_DAM_TEMPERATURE_BY_DEVICE_PROCESS
16 
17 #include <cmath>
18 
19 // Project includes
20 #include "includes/kratos_flags.h"
22 #include "processes/process.h"
23 
24 // Application include
26 
27 namespace Kratos
28 {
29 
31 {
32 
33  public:
35 
37 
38  //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
39 
42  Parameters &rParameters) : Process(Flags()), mrModelPart(rModelPart)
43  {
45 
46  //only include validation with c++11 since raw_literals do not exist in c++03
47  Parameters default_parameters(R"(
48  {
49  "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
50  "variable_name": "PLEASE_PRESCRIBE_VARIABLE_NAME",
51  "is_fixed" : false,
52  "value" : 0.0,
53  "table" : 0,
54  "position" : [0.0,0.0,0.0],
55  "interval":[
56  0.0,
57  0.0
58  ]
59  } )");
60 
61  // Some values need to be mandatorily prescribed since no meaningful default value exist. For this reason try accessing to them
62  // So that an error is thrown if they don't exist
63  rParameters["value"];
64  rParameters["variable_name"];
65  rParameters["model_part_name"];
66 
67  // Now validate agains defaults -- this also ensures no type mismatch
68  rParameters.ValidateAndAssignDefaults(default_parameters);
69 
70  mVariableName = rParameters["variable_name"].GetString();
71  mIsFixed = rParameters["is_fixed"].GetBool();
72  mValue = rParameters["value"].GetDouble();
73 
74  unsigned int Dim = rModelPart.GetProcessInfo()[DOMAIN_SIZE];
75 
76  // Getting the values of the device coordinates
77  mDeviceCoordinates.resize(Dim, false);
78  mDeviceCoordinates[0] = rParameters["position"][0].GetDouble();
79  mDeviceCoordinates[1] = rParameters["position"][1].GetDouble();
80  if (Dim ==3) mDeviceCoordinates[2] = rParameters["position"][2].GetDouble();
81 
82  mTimeUnitConverter = mrModelPart.GetProcessInfo()[TIME_UNIT_CONVERTER];
83  mTableId = rParameters["table"].GetInt();
84 
85  if (mTableId != 0)
87 
88  KRATOS_CATCH("");
89  }
90 
92 
95 
96  //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
97 
99  {
100 
101  KRATOS_TRY;
102 
103  const int nelements = mrModelPart.GetMesh(0).Elements().size();
105  bool IsInside = false;
106  array_1d<double, 3> LocalCoordinates;
107  Element::Pointer pSelectedElement;
108 
109  // Getting the values of table in case that it exist
110  if (mTableId != 0)
111  {
112  double time = mrModelPart.GetProcessInfo()[TIME];
114  mValue = mpTable->GetValue(time);
115  }
116 
117  if (nelements != 0)
118  {
119  ModelPart::ElementsContainerType::iterator el_begin = mrModelPart.ElementsBegin();
120  int PointsNumber = 0;
121 
122  for (int k = 0; k < nelements; k++)
123  {
124  ModelPart::ElementsContainerType::iterator it = el_begin + k;
125  pSelectedElement = (*(it.base()));
126  IsInside = pSelectedElement->GetGeometry().IsInside(mDeviceCoordinates, LocalCoordinates);
127  if (IsInside)
128  break;
129  }
130  if (IsInside == false)
131  {
132  for (int k = 0; k < nelements; k++)
133  {
134  ModelPart::ElementsContainerType::iterator it = el_begin + k;
135  pSelectedElement = (*(it.base()));
136  IsInside = pSelectedElement->GetGeometry().IsInside(mDeviceCoordinates, LocalCoordinates, 1.0e-5);
137  if (IsInside)
138  break;
139  }
140  }
141  if (IsInside == false)
142  {
143  KRATOS_ERROR << "ERROR!!, PLEASE REPEAT THE SEARCH " << std::endl;
144  }
145 
146  PointsNumber = pSelectedElement->GetGeometry().PointsNumber();
147 
148  for (int j = 0; j < PointsNumber; j++)
149  {
150  pSelectedElement->GetGeometry().GetPoint(j).FastGetSolutionStepValue(var) = mValue;
151  pSelectedElement->GetGeometry().GetPoint(j).Fix(var);
152  }
153  }
154 
155  KRATOS_CATCH("");
156  }
157 
159 
161  std::string Info() const override
162  {
163  return "DamTemperaturebyDeviceProcess";
164  }
165 
167  void PrintInfo(std::ostream &rOStream) const override
168  {
169  rOStream << "DamTemperaturebyDeviceProcess";
170  }
171 
173  void PrintData(std::ostream &rOStream) const override
174  {
175  }
176 
178 
179  protected:
181 
183  std::string mVariableName;
184  bool mIsFixed;
185  double mValue;
188  TableType::Pointer mpTable;
189  int mTableId;
190 
191  //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
192 
193  private:
196 
197 }; //Class
198 
200 inline std::istream &operator>>(std::istream &rIStream,
202 
204 inline std::ostream &operator<<(std::ostream &rOStream,
205  const DamTemperaturebyDeviceProcess &rThis)
206 {
207  rThis.PrintInfo(rOStream);
208  rOStream << std::endl;
209  rThis.PrintData(rOStream);
210 
211  return rOStream;
212 }
213 
214 } /* namespace Kratos.*/
215 
216 #endif /* KRATOS_DAM_TEMPERATURE_BY_DEVICE_PROCESS defined */
Definition: dam_temperature_by_device_process.hpp:31
KRATOS_CLASS_POINTER_DEFINITION(DamTemperaturebyDeviceProcess)
DamTemperaturebyDeviceProcess(ModelPart &rModelPart, Parameters &rParameters)
Constructor.
Definition: dam_temperature_by_device_process.hpp:41
int mTableId
Definition: dam_temperature_by_device_process.hpp:189
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: dam_temperature_by_device_process.hpp:167
virtual ~DamTemperaturebyDeviceProcess()
Destructor.
Definition: dam_temperature_by_device_process.hpp:94
double mValue
Definition: dam_temperature_by_device_process.hpp:185
std::string mVariableName
Definition: dam_temperature_by_device_process.hpp:183
double mTimeUnitConverter
Definition: dam_temperature_by_device_process.hpp:187
void ExecuteInitializeSolutionStep() override
This function will be executed at every time step BEFORE performing the solve phase.
Definition: dam_temperature_by_device_process.hpp:98
bool mIsFixed
Definition: dam_temperature_by_device_process.hpp:184
array_1d< double, 3 > mDeviceCoordinates
Definition: dam_temperature_by_device_process.hpp:186
std::string Info() const override
Turn back information as a string.
Definition: dam_temperature_by_device_process.hpp:161
ModelPart & mrModelPart
Member Variables.
Definition: dam_temperature_by_device_process.hpp:182
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: dam_temperature_by_device_process.hpp:173
Table< double, double > TableType
Definition: dam_temperature_by_device_process.hpp:36
TableType::Pointer mpTable
Definition: dam_temperature_by_device_process.hpp:188
Definition: flags.h:58
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
ElementsContainerType & Elements()
Definition: mesh.h:568
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
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
bool GetBool() const
This method returns the boolean contained in the current Parameter.
Definition: kratos_parameters.cpp:675
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
Definition: table.h:435
BOOST_UBLAS_INLINE void resize(size_type array_size, bool preserve=true)
Definition: array_1d.h:242
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR
Definition: exception.h:161
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 k
Definition: quadrature.py:595
int j
Definition: quadrature.py:648