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_vector_variable_to_pfem_conditions_process.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosPfemFluidDynamicsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: August 2016 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_ASSIGN_VECTOR_VARIABLE_TO_PFEM_CONDITIONS_PROCESS_H_INCLUDED)
11 #define KRATOS_ASSIGN_VECTOR_VARIABLE_TO_PFEM_CONDITIONS_PROCESS_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
19 
20 namespace Kratos
21 {
22 
25 
27 
30 {
31 public:
34 
37 
39 
44  Parameters rParameters) : BaseType(rModelPart)
45  {
47 
48  Parameters default_parameters(R"(
49  {
50  "model_part_name":"MODEL_PART_NAME",
51  "variable_name": "VARIABLE_NAME",
52  "value" : [0.0, 0.0, 0.0],
53  "compound_assignment": "direct"
54  } )");
55 
56  // Validate against defaults -- this ensures no type mismatch
57  rParameters.ValidateAndAssignDefaults(default_parameters);
58 
59  mvariable_name = rParameters["variable_name"].GetString();
60 
62  {
63  KRATOS_THROW_ERROR(std::runtime_error, "trying to set a variable that is not in the model_part - variable name is ", mvariable_name);
64  }
65 
66  mvector_value[0] = rParameters["value"][0].GetDouble();
67  mvector_value[1] = rParameters["value"][1].GetDouble();
68  mvector_value[2] = rParameters["value"][2].GetDouble();
69 
70  this->SetAssignmentType(rParameters["compound_assignment"].GetString(), mAssignment);
71 
72  KRATOS_CATCH("");
73  }
74 
76  const Variable<array_1d<double, 3>> &rVariable,
77  const array_1d<double, 3> &rvector_value) : BaseType(rModelPart), mvector_value(rvector_value)
78  {
79  KRATOS_TRY;
80 
81  mvariable_name = rVariable.Name();
82 
83  if (KratosComponents<Variable<array_1d<double, 3>>>::Has(mvariable_name) == false) //case of array_1d variable
84  KRATOS_THROW_ERROR(std::runtime_error, "trying to set a variable that is not in the model_part - variable name is ", mvariable_name);
85 
86  KRATOS_CATCH("")
87  }
88 
91 
95 
97  void operator()()
98  {
99  Execute();
100  }
101 
105 
107  void Execute() override
108  {
109 
110  KRATOS_TRY
111 
112  InternalAssignValue(KratosComponents<Variable<array_1d<double, 3>>>::Get(mvariable_name), mvector_value);
113 
114  KRATOS_CATCH("")
115  }
116 
119  void ExecuteInitialize() override
120  {
121  }
122 
126  {
127  }
128 
131  {
132  }
133 
136  {
137  }
138 
140  void ExecuteBeforeOutputStep() override
141  {
142  }
143 
145  void ExecuteAfterOutputStep() override
146  {
147  }
148 
151  void ExecuteFinalize() override
152  {
153  KRATOS_TRY
154 
156  array_1d<double, 3> vector_value;
157  vector_value.clear();
158  InternalAssignValue(KratosComponents<Variable<array_1d<double, 3>>>::Get(mvariable_name), vector_value);
159 
160  KRATOS_CATCH("")
161  }
162 
166 
170 
174 
176  std::string Info() const override
177  {
178  return "AssignVectorVariableToPfemConditionsProcess";
179  }
180 
182  void PrintInfo(std::ostream &rOStream) const override
183  {
184  rOStream << "AssignVectorVariableToPfemConditionsProcess";
185  }
186 
188  void PrintData(std::ostream &rOStream) const override
189  {
190  }
191 
196 
197 protected:
203 
206 
220 
221 private:
227 
228  array_1d<double, 3> mvector_value;
229 
233 
234  void InternalAssignValue(const Variable<array_1d<double, 3>> &rVariable,
235  const array_1d<double, 3> &rvector_value)
236  {
237 
239 
240  AssignmentMethodPointer AssignmentMethod = this->GetAssignmentMethod<AssignmentMethodPointer>();
241 
242  const int nconditions = this->mrModelPart.GetMesh().Conditions().size();
243 
244  if (nconditions != 0)
245  {
246  ModelPart::ConditionsContainerType::iterator it_begin = this->mrModelPart.GetMesh().ConditionsBegin();
247 
248 #pragma omp parallel for
249  for (int i = 0; i < nconditions; i++)
250  {
251  ModelPart::ConditionsContainerType::iterator it = it_begin + i;
252 
253  (this->*AssignmentMethod)(*it, rVariable, rvector_value);
254  }
255  }
256  }
257 
264 
267 
277 
278 }; // Class AssignVectorVariableToPfemConditionsProcess
279 
281 
284 
288 
290 inline std::istream &operator>>(std::istream &rIStream,
292 
294 inline std::ostream &operator<<(std::ostream &rOStream,
296 {
297  rThis.PrintInfo(rOStream);
298  rOStream << std::endl;
299  rThis.PrintData(rOStream);
300 
301  return rOStream;
302 }
304 
305 } // namespace Kratos.
306 
307 #endif // KRATOS_ASSIGN_VECTOR_VARIABLE_TO_PFEM_CONDITIONS_PROCESS_H_INCLUDED defined
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:34
ModelPart & mrModelPart
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:250
void SetAssignmentType(std::string method, AssignmentType &rAssignment)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:272
std::string mvariable_name
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:251
AssignmentType mAssignment
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:254
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:30
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:182
AssignScalarVariableToPfemEntitiesProcess BaseType
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:38
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:130
KRATOS_CLASS_POINTER_DEFINITION(AssignVectorVariableToPfemConditionsProcess)
Pointer definition of AssignVectorVariableToPfemConditionsProcess.
std::string Info() const override
Turn back information as a string.
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:176
void ExecuteFinalize() override
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:151
void Execute() override
Execute method is used to execute the AssignVectorVariableToPfemConditionsProcess algorithms.
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:107
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:97
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:135
AssignVectorVariableToPfemConditionsProcess(ModelPart &rModelPart, const Variable< array_1d< double, 3 >> &rVariable, const array_1d< double, 3 > &rvector_value)
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:75
AssignVectorVariableToPfemConditionsProcess(AssignVectorVariableToPfemConditionsProcess const &rOther)
Copy constructor.
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:188
void ExecuteBeforeSolutionLoop() override
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:125
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:145
AssignVectorVariableToPfemConditionsProcess(ModelPart &rModelPart, Parameters rParameters)
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:43
void ExecuteInitialize() override
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:119
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:140
~AssignVectorVariableToPfemConditionsProcess() override
Destructor.
Definition: assign_vector_variable_to_pfem_conditions_process.hpp:90
Base class for all Conditions.
Definition: condition.h:59
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
ConditionsContainerType & Conditions()
Definition: mesh.h:691
ConditionIterator ConditionsBegin()
Definition: mesh.h:671
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
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
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
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
BOOST_UBLAS_INLINE void clear()
Definition: array_1d.h:325
#define KRATOS_THROW_ERROR(ExceptionType, ErrorMessage, MoreInfo)
Definition: define.h:77
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
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
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
integer i
Definition: TensorModule.f:17