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.
face_heat_utilities.h
Go to the documentation of this file.
1 // KRATOS ___ ___ _ ___ __ ___ ___ ___ ___
2 // / __/ _ \| \| \ \ / /__| \_ _| __| __|
3 // | (_| (_) | .` |\ V /___| |) | || _|| _|
4 // \___\___/|_|\_| \_/ |___/___|_| |_| APPLICATION
5 //
6 // License: BSD License
7 // Kratos default license: kratos/license.txt
8 //
9 // Main authors: Riccardo Rossi
10 //
11 
12 #if !defined(KRATOS_FACE_HEAT_UTILITIES_INCLUDED )
13 #define KRATOS_FACE_HEAT_UTILITIES_INCLUDED
14 
15 
16 
17 // System includes
18 #include <string>
19 #include <iostream>
20 #include <algorithm>
21 
22 // External includes
23 
24 
25 // Project includes
26 #include "includes/define.h"
27 #include "includes/model_part.h"
29 #include "includes/node.h"
30 #include "utilities/geometry_utilities.h"
33 
34 namespace Kratos
35 {
37 {
38 public:
39 
40  //**********************************************************************************************
41  //**********************************************************************************************
44  double face_heat_source
45  )
46  {
48  for(ModelPart::ConditionsContainerType::iterator iii = conditions.begin(); iii != conditions.end(); iii++)
49  {
50  Geometry< Node >& geom = iii->GetGeometry();
51 
52  for(unsigned int k = 0; k<geom.size(); k++)
53  geom[k].FastGetSolutionStepValue(FACE_HEAT_FLUX) = face_heat_source;
54 
55  }
56  std::cout << "Conditions are generated" << std::endl;
57 
58  KRATOS_CATCH("");
59  }
60 
61  void GenerateModelPart( ModelPart& OriginModelPart, ModelPart& DestinationModelPart, Element const& rReferenceElement, Condition const& rReferenceBoundaryCondition)
62  {
63  KRATOS_TRY;
64 
65  //assigning the nodes to the new model part
66  DestinationModelPart.rProperties() = OriginModelPart.rProperties();
67  DestinationModelPart.Nodes().clear();
68  DestinationModelPart.Nodes() = OriginModelPart.Nodes();
69 
70  //generating the elements
71  int id = 1;
72  Properties::Pointer properties = OriginModelPart.GetMesh().pGetProperties(1);
73  for(ModelPart::ElementsContainerType::iterator iii = OriginModelPart.ElementsBegin(); iii != OriginModelPart.ElementsEnd(); iii++)
74  {
75  Geometry< Node >& geom = iii->GetGeometry();
76  Properties::Pointer properties = iii->pGetProperties();
77  Element::Pointer p_element = rReferenceElement.Create(id, geom ,properties);
78  DestinationModelPart.Elements().push_back(p_element);
79  id = id + 1;
80  }
81  std::cout << "Elements are generated" << std::endl;
82 
83  //generating the conditions
84  id = 1;
85  for(ModelPart::ConditionsContainerType::iterator iii = OriginModelPart.ConditionsBegin(); iii != OriginModelPart.ConditionsEnd(); iii++)
86  {
87  Geometry< Node >& geom = iii->GetGeometry();
88  double nfree_surf = 0;
89  for(unsigned int k = 0; k<geom.size(); k++)
90  nfree_surf += geom[k].FastGetSolutionStepValue(IS_FREE_SURFACE);
91 
92  if(nfree_surf > 1)
93  {
94  Properties::Pointer properties = iii->pGetProperties();
95  Condition::Pointer p_condition = rReferenceBoundaryCondition.Create(id, geom,properties);
96  DestinationModelPart.Conditions().push_back(p_condition);
97  id = id + 1;
98  }
99  }
100  std::cout << "Conditions are generated" << std::endl;
101 
102  KRATOS_CATCH("");
103  }
104 
105 
106 
107 
108 
109  void ConditionModelPart(ModelPart& temperature_model_part, ModelPart& full_model_part, const int TDim)
110  {
111  KRATOS_TRY;
112 
113 
114  temperature_model_part.Conditions().clear();
115  int nd=TDim;
116  Properties::Pointer properties = full_model_part.GetMesh().pGetProperties(1);
117 
118 
119  //int id=full_model_part.Conditions().size();
120  //int n_int=0.0;
121 
122 
123  for(ModelPart::ElementsContainerType::iterator im = full_model_part.ElementsBegin() ; im != full_model_part.ElementsEnd() ; ++im)
124  {
125  if(nd==2)
126  {
127 
128  int n_int=im->GetGeometry()[0].FastGetSolutionStepValue(MATERIAL_VARIABLE);
129  for(int j=1; j<nd+1; j++) n_int+= im->GetGeometry()[j].FastGetSolutionStepValue(MATERIAL_VARIABLE);
130 
131  if (n_int==3)
132  {
133  }
134  else
135  {
136 
137 
138  n_int=im->GetGeometry()[1].FastGetSolutionStepValue(MATERIAL_VARIABLE) + im->GetGeometry()[2].FastGetSolutionStepValue(MATERIAL_VARIABLE);
139 
140  if(n_int==2)
141  {
142 
144  temp1.reserve(2);
145  temp1.push_back(im->GetGeometry()(1));
146  temp1.push_back(im->GetGeometry()(2));
148  int id = (im->Id()-1)*3;
149  Condition::Pointer p_cond = Kratos::make_intrusive<ThermalFace>(id, cond, properties);
150  temperature_model_part.Conditions().push_back(p_cond);
151 
152  }
153 
154  n_int= im->GetGeometry()[2].FastGetSolutionStepValue(MATERIAL_VARIABLE) + im->GetGeometry()[0].FastGetSolutionStepValue(MATERIAL_VARIABLE);
155 
156 
157  if(n_int==2)
158  {
159 
161  temp1.reserve(2);
162  temp1.push_back(im->GetGeometry()(2));
163  temp1.push_back(im->GetGeometry()(0));
165  int id = (im->Id()-1)*3+1;
166  Condition::Pointer p_cond = Kratos::make_intrusive<ThermalFace>(id, cond, properties);
167  temperature_model_part.Conditions().push_back(p_cond);
168  }
169 
170 
171  n_int= im->GetGeometry()[0].FastGetSolutionStepValue(MATERIAL_VARIABLE) + im->GetGeometry()[1].FastGetSolutionStepValue(MATERIAL_VARIABLE) ;
172 
173  if(n_int==2)
174  {
175 
177  temp1.reserve(2);
178  temp1.push_back(im->GetGeometry()(0));
179  temp1.push_back(im->GetGeometry()(1));
181  int id = (im->Id()-1)*3+2;
182 
183  Condition::Pointer p_cond = Kratos::make_intrusive<ThermalFace>(id, cond, properties);
184  temperature_model_part.Conditions().push_back(p_cond);
185 
186  }
187 
188  }
189 
190  }
191 
192  else
193  {
194 
195  int n_int=im->GetGeometry()[0].FastGetSolutionStepValue(MATERIAL_VARIABLE);
196  for(int j=1; j<nd+1; j++) n_int+= im->GetGeometry()[j].FastGetSolutionStepValue(MATERIAL_VARIABLE);
197  if (n_int==4)
198 
199  {
200 
201  }
202  else
203  {
204  n_int=0.0;
205  n_int=im->GetGeometry()[1].FastGetSolutionStepValue(MATERIAL_VARIABLE);
206  n_int+=im->GetGeometry()[2].FastGetSolutionStepValue(MATERIAL_VARIABLE);
207  n_int+=im->GetGeometry()[3].FastGetSolutionStepValue(MATERIAL_VARIABLE);
208  if(n_int==3)
209  {
210 
212  temp.reserve(3);
213  temp.push_back(im->GetGeometry()(1));
214  temp.push_back(im->GetGeometry()(2));
215  temp.push_back(im->GetGeometry()(3));
217  int id = (im->Id()-1)*4;
218  Condition::Pointer p_cond = Kratos::make_intrusive<ThermalFace>(id, cond, properties);
219  temperature_model_part.Conditions().push_back(p_cond);
220  }
221  n_int=0.0;
222  n_int=im->GetGeometry()[0].FastGetSolutionStepValue(MATERIAL_VARIABLE);
223  n_int+=im->GetGeometry()[3].FastGetSolutionStepValue(MATERIAL_VARIABLE);
224  n_int+=im->GetGeometry()[2].FastGetSolutionStepValue(MATERIAL_VARIABLE);
225  if(n_int==3)
226  {
227 
229  temp.reserve(3);
230  temp.push_back(im->GetGeometry()(0));
231  temp.push_back(im->GetGeometry()(3));
232  temp.push_back(im->GetGeometry()(2));
234  int id = (im->Id()-1)*4;
235  Condition::Pointer p_cond = Kratos::make_intrusive<ThermalFace>(id, cond, properties);
236  temperature_model_part.Conditions().push_back(p_cond);
237  }
238  n_int=0.0;
239  n_int=im->GetGeometry()[0].FastGetSolutionStepValue(MATERIAL_VARIABLE);
240  n_int+=im->GetGeometry()[1].FastGetSolutionStepValue(MATERIAL_VARIABLE);
241  n_int+=im->GetGeometry()[3].FastGetSolutionStepValue(MATERIAL_VARIABLE);
242  if(n_int==3)
243  {
244 
246  temp.reserve(3);
247  temp.push_back(im->GetGeometry()(0));
248  temp.push_back(im->GetGeometry()(1));
249  temp.push_back(im->GetGeometry()(3));
251  int id = (im->Id()-1)*4;
252  Condition::Pointer p_cond = Kratos::make_intrusive<ThermalFace>(id, cond, properties);
253  temperature_model_part.Conditions().push_back(p_cond);
254  }
255 
256 
257  n_int=0.0;
258  n_int=im->GetGeometry()[0].FastGetSolutionStepValue(MATERIAL_VARIABLE);
259  n_int+=im->GetGeometry()[2].FastGetSolutionStepValue(MATERIAL_VARIABLE);
260  n_int+=im->GetGeometry()[1].FastGetSolutionStepValue(MATERIAL_VARIABLE);
261  if(n_int==3)
262  {
263 
265  temp.reserve(3);
266  temp.push_back(im->GetGeometry()(0));
267  temp.push_back(im->GetGeometry()(2));
268  temp.push_back(im->GetGeometry()(1));
270  int id = (im->Id()-1)*4;
271  Condition::Pointer p_cond = Kratos::make_intrusive<ThermalFace>(id, cond, properties);
272  temperature_model_part.Conditions().push_back(p_cond);
273 
274  }
275 
276  }
277 
278  }
279  }
280 
281 
282  KRATOS_CATCH("");
283  }
284 
285 
286 private:
287 
288 };
289 
290 } // namespace Kratos.
291 
292 #endif // KRATOS_FACE_HEAT_UTILITIES_INCLUDED defined
Base class for all Conditions.
Definition: condition.h:59
virtual Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const
It creates a new condition pointer.
Definition: condition.h:205
Base class for all Elements.
Definition: element.h:60
virtual Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const
It creates a new element pointer.
Definition: element.h:202
Definition: face_heat_utilities.h:37
void GenerateModelPart(ModelPart &OriginModelPart, ModelPart &DestinationModelPart, Element const &rReferenceElement, Condition const &rReferenceBoundaryCondition)
Definition: face_heat_utilities.h:61
void ApplyFaceHeat(ModelPart::ConditionsContainerType &conditions, double face_heat_source)
Definition: face_heat_utilities.h:42
void ConditionModelPart(ModelPart &temperature_model_part, ModelPart &full_model_part, const int TDim)
Definition: face_heat_utilities.h:109
Geometry base class.
Definition: geometry.h:71
SizeType size() const
Definition: geometry.h:518
PropertiesType::Pointer pGetProperties(IndexType PropertiesId)
Definition: mesh.h:394
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
MeshType::ConditionsContainerType ConditionsContainerType
Condintions container. A vector set of Conditions with their Id's as key.
Definition: model_part.h:183
ConditionIterator ConditionsBegin(IndexType ThisIndex=0)
Definition: model_part.h:1361
ConditionsContainerType & Conditions(IndexType ThisIndex=0)
Definition: model_part.h:1381
ElementsContainerType & Elements(IndexType ThisIndex=0)
Definition: model_part.h:1189
PropertiesContainerType & rProperties(IndexType ThisIndex=0)
Definition: model_part.h:1003
ElementIterator ElementsEnd(IndexType ThisIndex=0)
Definition: model_part.h:1179
MeshType & GetMesh(IndexType ThisIndex=0)
Definition: model_part.h:1791
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
ConditionIterator ConditionsEnd(IndexType ThisIndex=0)
Definition: model_part.h:1371
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
void reserve(size_type dim)
Definition: pointer_vector.h:319
void push_back(const TPointerType &x)
Definition: pointer_vector.h:270
A three node 3D triangle geometry with linear shape functions.
Definition: triangle_3d_3.h:77
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#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
int k
Definition: quadrature.py:595
int j
Definition: quadrature.py:648
float temp
Definition: rotating_cone.py:85
temperature_model_part
Definition: script_THERMAL_CORRECT.py:32