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.
assign_torque_field_about_an_axis_to_conditions_process.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosSolidMechanicsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: January 2018 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_ASSIGN_TORQUE_FIELD_ABOUT_AN_AXIS_TO_CONDITIONS_PROCESS_H_INCLUDED)
11 #define KRATOS_ASSIGN_TORQUE_FIELD_ABOUT_AN_AXIS_TO_CONDITIONS_PROCESS_H_INCLUDED
12 
13 
14 
15 // System includes
16 
17 // External includes
18 
19 // Project includes
21 
22 namespace Kratos
23 {
24 
27 
29 
32 {
33 public:
36 
39 
43 
45  pybind11::object& pPyObject,
46  const std::string& pPyMethodName,
47  const bool SpatialFieldFunction,
48  Parameters rParameters
50  {
52 
53  Parameters default_parameters( R"(
54  {
55  "model_part_name":"MODEL_PART_NAME",
56  "variable_name": "VARIABLE_NAME",
57  "modulus" : 1.0,
58  "direction" : [],
59  "center" : []
60  } )" );
61 
62 
63  // Validate against defaults -- this ensures no type mismatch
64  rParameters.ValidateAndAssignDefaults(default_parameters);
65 
66  mvariable_name = rParameters["variable_name"].GetString();
67 
68  if( KratosComponents< Variable<array_1d<double, 3> > >::Has( mvariable_name ) ) //case of array_1d variable
69  {
70 
71  mPyObject = pPyObject;
72  mPyMethodName = pPyMethodName;
73 
74  mIsSpatialField = SpatialFieldFunction;
75 
76 
77  for( unsigned int i=0; i<3; i++)
78  {
79  mdirection[i] = rParameters["direction"][i].GetDouble();
80  mcenter[i] = rParameters["center"][i].GetDouble();
81  }
82 
83  double norm = norm_2(mdirection);
84  if(norm!=0)
85  mdirection/=norm;
86  }
87  else //case of other variable type
88  {
89  KRATOS_ERROR << "trying to set a variable that is not in the model_part - variable name is " << mvariable_name <<std::endl;
90  }
91 
92  KRATOS_CATCH("");
93  }
94 
95 
98 
99 
103 
105  void operator()()
106  {
107  Execute();
108  }
109 
110 
114 
115 
117  void Execute() override
118  {
119 
120  KRATOS_TRY;
121 
122 
123  if( ! mIsSpatialField ){
124 
125  const ProcessInfo& rCurrentProcessInfo = mrModelPart.GetProcessInfo();
126  const double& rCurrentTime = rCurrentProcessInfo[TIME];
127 
128  this->CallTimeFunction(rCurrentTime, mvalue);
129 
131 
132  }
133  else //no spatial fields accepted (it have to be implemented if needed
134  {
135  KRATOS_ERROR << "trying to set an spatial field....not implemented" << mvariable_name <<std::endl;
136  }
137 
138 
139  KRATOS_CATCH("");
140 
141  }
142 
145  void ExecuteInitialize() override
146  {
147  }
148 
152  {
153  }
154 
155 
158  {
159  }
160 
163  {
164  }
165 
166 
168  void ExecuteBeforeOutputStep() override
169  {
170  }
171 
172 
174  void ExecuteAfterOutputStep() override
175  {
176  }
177 
178 
181  void ExecuteFinalize() override
182  {
184  }
185 
186 
190 
191 
195 
196 
200 
202  std::string Info() const override
203  {
204  return "AssignTorqueFieldAboutAnAxisToConditionsProcess";
205  }
206 
208  void PrintInfo(std::ostream& rOStream) const override
209  {
210  rOStream << "AssignTorqueFieldAboutAnAxisToConditionsProcess";
211  }
212 
214  void PrintData(std::ostream& rOStream) const override
215  {
216  }
217 
218 
223 
224 protected:
225 
231 
232  pybind11::object mPyObject;
233  std::string mPyMethodName;
234 
236 
240 
243 
257 
258 private:
259 
265 
266 
270 
271 
272  void CallFunction(const Node::Pointer& pNode, const double& time, double& rValue)
273  {
274  KRATOS_TRY
275 
276  if( mIsSpatialField ){
277 
278  double x = pNode->X(), y = pNode->Y(), z = pNode->Z();
279 
280  rValue = mPyObject.attr(mPyMethodName.c_str())(x,y,z,time).cast<double>();
281  }
282  else{
283 
284  rValue = mPyObject.attr(mPyMethodName.c_str())(0.0,0.0,0.0,time).cast<double>();
285  }
286 
287  KRATOS_CATCH( "" )
288 
289  }
290 
291  void CallTimeFunction(const double& time, double& rValue)
292  {
293 
294  KRATOS_TRY
295 
296  rValue = mPyObject.attr(mPyMethodName.c_str())(0.0,0.0,0.0,time).cast<double>();
297 
298  KRATOS_CATCH( "" )
299 
300  }
301 
302 
309 
312 
313 
323 
324 }; // Class AssignTorqueFieldAboutAnAxisToConditionsProcess
325 
326 
328 
331 
332 
336 
337 
339 inline std::istream& operator >> (std::istream& rIStream,
341 
343 inline std::ostream& operator << (std::ostream& rOStream,
345 {
346  rThis.PrintInfo(rOStream);
347  rOStream << std::endl;
348  rThis.PrintData(rOStream);
349 
350  return rOStream;
351 }
353 
354 
355 } // namespace Kratos.
356 
357 #endif // KRATOS_ASSIGN_TORQUE_FIELD_ABOUT_AN_AXIS_TO_CONDITIONS_PROCESS_H_INCLUDED defined
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: assign_torque_about_an_axis_to_conditions_process.hpp:34
ModelPart & mrModelPart
Definition: assign_torque_about_an_axis_to_conditions_process.hpp:221
double mvalue
Definition: assign_torque_about_an_axis_to_conditions_process.hpp:223
array_1d< double, 3 > mcenter
Definition: assign_torque_about_an_axis_to_conditions_process.hpp:225
void Execute() override
Execute method is used to execute the AssignTorqueAboutAnAxisToConditionsProcess algorithms.
Definition: assign_torque_about_an_axis_to_conditions_process.hpp:119
array_1d< double, 3 > mdirection
Definition: assign_torque_about_an_axis_to_conditions_process.hpp:224
void ExecuteFinalize() override
Definition: assign_torque_about_an_axis_to_conditions_process.hpp:168
std::string mvariable_name
Definition: assign_torque_about_an_axis_to_conditions_process.hpp:222
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:32
AssignTorqueFieldAboutAnAxisToConditionsProcess(AssignTorqueFieldAboutAnAxisToConditionsProcess const &rOther)
Copy constructor.
void ExecuteInitialize() override
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:145
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:105
void ExecuteFinalize() override
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:181
void Execute() override
Execute method is used to execute the AssignTorqueFieldAboutAnAxisToConditionsProcess algorithms.
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:117
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:214
~AssignTorqueFieldAboutAnAxisToConditionsProcess() override
Destructor.
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:97
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:157
std::string Info() const override
Turn back information as a string.
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:202
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:174
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:208
bool mIsSpatialField
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:235
AssignTorqueFieldAboutAnAxisToConditionsProcess(ModelPart &model_part, pybind11::object &pPyObject, const std::string &pPyMethodName, const bool SpatialFieldFunction, Parameters rParameters)
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:44
pybind11::object mPyObject
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:232
void ExecuteBeforeSolutionLoop() override
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:151
KRATOS_CLASS_POINTER_DEFINITION(AssignTorqueFieldAboutAnAxisToConditionsProcess)
Pointer definition of AssignTorqueFieldAboutAnAxisToConditionsProcess.
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:168
std::string mPyMethodName
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:233
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: assign_torque_field_about_an_axis_to_conditions_process.hpp:162
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
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
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
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR
Definition: exception.h:161
z
Definition: GenerateWind.py:163
bool Has(const std::string &ModelerName)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:24
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
TExpressionType::data_type norm_2(AMatrix::MatrixExpression< TExpressionType, TCategory > const &TheExpression)
Definition: amatrix_interface.h:625
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
model_part
Definition: face_heat.py:14
y
Other simbols definition.
Definition: generate_axisymmetric_navier_stokes_element.py:54
x
Definition: sensitivityMatrix.py:49
integer i
Definition: TensorModule.f:17