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.
apply_component_table_process.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosDamApplication $
3 // Last modified by: $Author: David J. Vicente $
4 // Date: $Date: September 2018 $
5 // Revision: $Revision: 1.0 $
6 //
7 
8 #if !defined(KRATOS_APPLY_COMPONENT_TABLE_PROCESS )
9 #define KRATOS_APPLY_COMPONENT_TABLE_PROCESS
10 
11 #include "includes/table.h"
12 #include "includes/kratos_flags.h"
14 #include "processes/process.h"
15 
17 
18 namespace Kratos
19 {
20 
22 {
23 
24 public:
25 
27 
30 
32 
35  Parameters rParameters
37  {
39 
40  //only include validation with c++11 since raw_literals do not exist in c++03
41  Parameters default_parameters( R"(
42  {
43  "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
44  "variable_name": "PLEASE_PRESCRIBE_VARIABLE_NAME",
45  "is_fixed": false,
46  "value" : 1.0,
47  "table" : 1,
48  "interval":[
49  0.0,
50  0.0
51  ]
52  } )" );
53 
54  // Some values need to be mandatorily prescribed since no meaningful default value exist. For this reason try accessing to them
55  // So that an error is thrown if they don't exist
56  rParameters["table"];
57  rParameters["variable_name"];
58  rParameters["model_part_name"];
59 
60  // Now validate agains defaults -- this also ensures no type mismatch
61  rParameters.ValidateAndAssignDefaults(default_parameters);
62 
63  mvariable_name = rParameters["variable_name"].GetString();
64  mis_fixed = rParameters["is_fixed"].GetBool();
65  minitial_value = rParameters["value"].GetDouble();
66 
67  unsigned int TableId = rParameters["table"].GetInt();
68  mpTable = model_part.pGetTable(TableId);
69  mTimeUnitConverter = model_part.GetProcessInfo()[TIME_UNIT_CONVERTER];
70 
71  KRATOS_CATCH("");
72  }
73 
75 
78 
79 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
80 
82  void Execute() override
83  {
84  }
85 
86 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
87 
90  void ExecuteInitialize() override
91  {
92  KRATOS_TRY;
93 
94  typedef Variable<double> component_type;
95  const component_type& var_component = KratosComponents< component_type >::Get(mvariable_name);
96 
97  const int nnodes = static_cast<int>(mr_model_part.Nodes().size());
98 
99  if(nnodes != 0)
100  {
101  ModelPart::NodesContainerType::iterator it_begin = mr_model_part.NodesBegin();
102 
103  #pragma omp parallel for
104  for(int i = 0; i<nnodes; i++)
105  {
106  ModelPart::NodesContainerType::iterator it = it_begin + i;
107 
108  if(mis_fixed)
109  {
110  it->Fix(var_component);
111  }
112 
113  it->FastGetSolutionStepValue(var_component) = minitial_value;
114  }
115  }
116 
117  KRATOS_CATCH("");
118  }
119 
122  {
123  KRATOS_TRY;
124 
125  typedef Variable<double> component_type;
126  const component_type& var_component = KratosComponents< component_type >::Get(mvariable_name);
127 
128  const double Time = mr_model_part.GetProcessInfo()[TIME]/mTimeUnitConverter;
129  double value = mpTable->GetValue(Time);
130 
131  const int nnodes = static_cast<int>(mr_model_part.Nodes().size());
132 
133  if(nnodes != 0)
134  {
135  ModelPart::NodesContainerType::iterator it_begin = mr_model_part.NodesBegin();
136 
137  #pragma omp parallel for
138  for(int i = 0; i<nnodes; i++)
139  {
140  ModelPart::NodesContainerType::iterator it = it_begin + i;
141 
142  it->FastGetSolutionStepValue(var_component) = value;
143  }
144  }
145 
146  KRATOS_CATCH("");
147  }
148 
150  std::string Info() const override
151  {
152  return "ApplyComponentTableProcessDam";
153  }
154 
156  void PrintInfo(std::ostream& rOStream) const override
157  {
158  rOStream << "ApplyComponentTableProcessDam";
159  }
160 
162  void PrintData(std::ostream& rOStream) const override
163  {
164  }
165 
167 
168 protected:
169 
171 
173  std::string mvariable_name;
174  bool mis_fixed;
176  TableType::Pointer mpTable;
178 
180 
181 private:
182 
185 
187  //ApplyComponentTableProcessDam(ApplyComponentTableProcessDam const& rOther);
188 
189 }; // Class ApplyComponentTableProcessDam
190 
192 inline std::istream& operator >> (std::istream& rIStream,
194 
196 inline std::ostream& operator << (std::ostream& rOStream,
197  const ApplyComponentTableProcessDam& rThis)
198 {
199  rThis.PrintInfo(rOStream);
200  rOStream << std::endl;
201  rThis.PrintData(rOStream);
202 
203  return rOStream;
204 }
205 
206 } // namespace Kratos.
207 
208 #endif /* KRATOS_APPLY_COMPONENT_TABLE_PROCESS defined */
Definition: apply_component_table_process.hpp:22
ApplyComponentTableProcessDam(ModelPart &model_part, Parameters rParameters)
Constructor.
Definition: apply_component_table_process.hpp:34
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: apply_component_table_process.hpp:121
KRATOS_CLASS_POINTER_DEFINITION(ApplyComponentTableProcessDam)
~ApplyComponentTableProcessDam() override
Destructor.
Definition: apply_component_table_process.hpp:77
double minitial_value
Definition: apply_component_table_process.hpp:175
TableType::Pointer mpTable
Definition: apply_component_table_process.hpp:176
std::string mvariable_name
Definition: apply_component_table_process.hpp:173
bool mis_fixed
Definition: apply_component_table_process.hpp:174
Table< double, double > TableType
Defining a table with double argument and result type as table type.
Definition: apply_component_table_process.hpp:29
ModelPart & mr_model_part
Member Variables.
Definition: apply_component_table_process.hpp:172
void Execute() override
Execute method is used to execute the ApplyComponentTableProcessDam algorithms.
Definition: apply_component_table_process.hpp:82
void ExecuteInitialize() override
Definition: apply_component_table_process.hpp:90
std::string Info() const override
Turn back information as a string.
Definition: apply_component_table_process.hpp:150
double mTimeUnitConverter
Definition: apply_component_table_process.hpp:177
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: apply_component_table_process.hpp:162
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: apply_component_table_process.hpp:156
Definition: flags.h:58
static const TComponentType & Get(const std::string &rName)
Retrieves a component with the specified name.
Definition: kratos_components.h:114
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
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
The base class for all processes in Kratos.
Definition: process.h:49
Definition: table.h:435
#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
model_part
Definition: face_heat.py:14
int nnodes
Definition: sensitivityMatrix.py:24
integer i
Definition: TensorModule.f:17