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: KratosPoromechanicsApplication $
3 // Last modified by: $Author: Ignasi de Pouplana $
4 // Date: $Date: June 2016 $
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 
21 class ApplyComponentTableProcess : public Process
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  } )" );
49 
50  // Some values need to be mandatorily prescribed since no meaningful default value exist. For this reason try accessing to them
51  // So that an error is thrown if they don't exist
52  rParameters["table"];
53  rParameters["variable_name"];
54  rParameters["model_part_name"];
55 
56  // Now validate agains defaults -- this also ensures no type mismatch
57  rParameters.ValidateAndAssignDefaults(default_parameters);
58 
59  mvariable_name = rParameters["variable_name"].GetString();
60  mis_fixed = rParameters["is_fixed"].GetBool();
61  minitial_value = rParameters["value"].GetDouble();
62 
63  unsigned int TableId = rParameters["table"].GetInt();
64  mpTable = model_part.pGetTable(TableId);
65  mTimeUnitConverter = model_part.GetProcessInfo()[TIME_UNIT_CONVERTER];
66 
67  KRATOS_CATCH("");
68  }
69 
71 
74 
76 
78  void Execute() override
79  {
80  }
81 
84  void ExecuteInitialize() override
85  {
86  KRATOS_TRY;
87 
88  typedef Variable<double> component_type;
89  const component_type& var_component = KratosComponents< component_type >::Get(mvariable_name);
90 
91  const int nnodes = static_cast<int>(mr_model_part.Nodes().size());
92 
93  if(nnodes != 0)
94  {
95  ModelPart::NodesContainerType::iterator it_begin = mr_model_part.NodesBegin();
96 
97  #pragma omp parallel for
98  for(int i = 0; i<nnodes; i++)
99  {
100  ModelPart::NodesContainerType::iterator it = it_begin + i;
101 
102  if(mis_fixed)
103  {
104  it->Fix(var_component);
105  }
106 
107  it->FastGetSolutionStepValue(var_component) = minitial_value;
108  }
109  }
110 
111  KRATOS_CATCH("");
112  }
113 
116  {
117  KRATOS_TRY;
118 
119  typedef Variable<double> component_type;
120  const component_type& var_component = KratosComponents< component_type >::Get(mvariable_name);
121 
122  const double Time = mr_model_part.GetProcessInfo()[TIME]/mTimeUnitConverter;
123  double value = mpTable->GetValue(Time);
124 
125  const int nnodes = static_cast<int>(mr_model_part.Nodes().size());
126 
127  if(nnodes != 0)
128  {
129  ModelPart::NodesContainerType::iterator it_begin = mr_model_part.NodesBegin();
130 
131  #pragma omp parallel for
132  for(int i = 0; i<nnodes; i++)
133  {
134  ModelPart::NodesContainerType::iterator it = it_begin + i;
135 
136  it->FastGetSolutionStepValue(var_component) = value;
137  }
138  }
139 
140  KRATOS_CATCH("");
141  }
142 
144  std::string Info() const override
145  {
146  return "ApplyComponentTableProcess";
147  }
148 
150  void PrintInfo(std::ostream& rOStream) const override
151  {
152  rOStream << "ApplyComponentTableProcess";
153  }
154 
156  void PrintData(std::ostream& rOStream) const override
157  {
158  }
159 
161 
162 protected:
163 
165 
167  std::string mvariable_name;
168  bool mis_fixed;
170  TableType::Pointer mpTable;
171  double mTimeUnitConverter;
172 
174 
175 private:
176 
178  ApplyComponentTableProcess& operator=(ApplyComponentTableProcess const& rOther);
179 
181  //ApplyComponentTableProcess(ApplyComponentTableProcess const& rOther);
182 
183 }; // Class ApplyComponentTableProcess
184 
186 inline std::istream& operator >> (std::istream& rIStream,
188 
190 inline std::ostream& operator << (std::ostream& rOStream,
191  const ApplyComponentTableProcess& rThis)
192 {
193  rThis.PrintInfo(rOStream);
194  rOStream << std::endl;
195  rThis.PrintData(rOStream);
196 
197  return rOStream;
198 }
199 
200 } // namespace Kratos.
201 
202 #endif /* KRATOS_APPLY_COMPONENT_TABLE_PROCESS defined */
Definition: apply_component_table_process.hpp:28
std::string mvariable_name
Definition: apply_component_table_process.hpp:167
double minitial_value
Definition: apply_component_table_process.hpp:169
~ApplyComponentTableProcess() override
Destructor.
Definition: apply_component_table_process.hpp:73
Table< double, double > TableType
Defining a table with double argument and result type as table type.
Definition: apply_component_table_process.hpp:29
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: apply_component_table_process.hpp:150
void Execute() override
Execute method is used to execute the ApplyComponentTableProcess algorithms.
Definition: apply_component_table_process.hpp:78
bool mis_fixed
Definition: apply_component_table_process.hpp:168
ApplyComponentTableProcess(ModelPart &model_part, Parameters rParameters)
Constructor.
Definition: apply_component_table_process.hpp:34
double mTimeUnitConverter
Definition: apply_component_table_process.hpp:173
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: apply_component_table_process.hpp:115
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: apply_component_table_process.hpp:156
TableType::Pointer mpTable
Definition: apply_component_table_process.hpp:172
void ExecuteInitialize() override
Definition: apply_component_table_process.hpp:84
std::string Info() const override
Turn back information as a string.
Definition: apply_component_table_process.hpp:144
KRATOS_CLASS_POINTER_DEFINITION(ApplyComponentTableProcess)
ModelPart & mr_model_part
Member Variables.
Definition: apply_component_table_process.hpp:166
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