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 // KRATOS___
2 // // ) )
3 // // ___ ___
4 // // ____ //___) ) // ) )
5 // // / / // // / /
6 // ((____/ / ((____ ((___/ / MECHANICS
7 //
8 // License: geo_mechanics_application/license.txt
9 //
10 // Main authors: Ignasi de Pouplana,
11 // Vahid Galavi
12 //
13 
14 #pragma once
15 
16 #include "includes/table.h"
17 #include "includes/kratos_flags.h"
19 #include "processes/process.h"
20 
22 
23 #include <boost/algorithm/string.hpp>
24 
25 namespace Kratos
26 {
27 
28 class ApplyComponentTableProcess : public Process
29 {
30 
31 public:
32 
34 
37 
39  Parameters rParameters
41  {
43 
44  //only include validation with c++11 since raw_literals do not exist in c++03
45  Parameters default_parameters( R"(
46  {
47  "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
48  "variable_name": "PLEASE_PRESCRIBE_VARIABLE_NAME",
49  "is_fixed": false,
50  "value" : 1.0,
51  "table" : 1
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  mIsFixedProvided = rParameters.Has("is_fixed");
61  // Now validate against defaults -- this also ensures no type mismatch
62  rParameters.ValidateAndAssignDefaults(default_parameters);
63 
64  mVariableName = rParameters["variable_name"].GetString();
65  mIsFixed = rParameters["is_fixed"].GetBool();
66  mInitialValue = rParameters["value"].GetDouble();
67 
68  unsigned int TableId = rParameters["table"].GetInt();
69  mpTable = model_part.pGetTable(TableId);
70  mTimeUnitConverter = model_part.GetProcessInfo()[TIME_UNIT_CONVERTER];
71 
72  KRATOS_CATCH("")
73  }
76  ~ApplyComponentTableProcess() override = default;
77 
80  void ExecuteInitialize() override
81  {
83 
85  << "Variable '" << mVariableName << "' is not defined!" << std::endl;
86 
87  const auto& var = KratosComponents<Variable<double>>::Get(mVariableName);
88 
89  auto variable_name_1 = mpTable->NameOfX();
90  boost::to_upper(variable_name_1);
91  if (variable_name_1 == "TIME") {
92  block_for_each(mrModelPart.Nodes(), [&var, this](Node& rNode) {
93  if (mIsFixed) rNode.Fix(var);
94  else if (mIsFixedProvided) rNode.Free(var);
95  rNode.FastGetSolutionStepValue(var) = mInitialValue;
96  });
97  }
98  else if (variable_name_1 == "X") {
99  block_for_each(mrModelPart.Nodes(), [&var, this](auto& node){
100  node.FastGetSolutionStepValue(var) = mpTable->GetValue(node.X());
101  });
102  }
103  else {
104  KRATOS_ERROR << "Failed to initialize ApplyComponentTableProcess: got unknown table variable '"
105  << variable_name_1 << "'";
106  }
107 
108  KRATOS_CATCH("")
109  }
110 
113  {
114  KRATOS_TRY
115 
116  const auto& var = KratosComponents<Variable<double>>::Get(mVariableName);
117 
118  auto variable_name_1 = mpTable->NameOfX();
119  boost::to_upper(variable_name_1);
120  if (variable_name_1 == "TIME") {
121  const double Time = mrModelPart.GetProcessInfo()[TIME]/mTimeUnitConverter;
122  const double value = mpTable->GetValue(Time);
123  block_for_each(mrModelPart.Nodes(), [&var, &value](Node& rNode) {
124  rNode.FastGetSolutionStepValue(var) = value;
125  });
126  }
127 
128  KRATOS_CATCH("")
129  }
130 
132  std::string Info() const override
133  {
134  return "ApplyComponentTableProcess";
135  }
136 
137 private:
140  std::string mVariableName;
141  bool mIsFixed;
142  bool mIsFixedProvided;
143  double mInitialValue;
144  TableType::Pointer mpTable;
145  double mTimeUnitConverter;
146 };
147 
148 }
Definition: apply_component_table_process.hpp:28
ModelPart & mrModelPart
Member Variables.
Definition: apply_component_table_process.hpp:168
ApplyComponentTableProcess & operator=(const ApplyComponentTableProcess &)=delete
ApplyComponentTableProcess(ModelPart &model_part, Parameters rParameters)
Definition: apply_component_table_process.hpp:38
ApplyComponentTableProcess(const ApplyComponentTableProcess &)=delete
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:112
TableType::Pointer mpTable
Definition: apply_component_table_process.hpp:172
~ApplyComponentTableProcess() override=default
void ExecuteInitialize() override
Definition: apply_component_table_process.hpp:80
std::string Info() const override
Turn back information as a string.
Definition: apply_component_table_process.hpp:132
bool mIsFixed
Definition: apply_component_table_process.hpp:170
double mInitialValue
Definition: apply_component_table_process.hpp:171
std::string mVariableName
Definition: apply_component_table_process.hpp:169
KRATOS_CLASS_POINTER_DEFINITION(ApplyComponentTableProcess)
Definition: flags.h:58
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
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
This class defines the node.
Definition: node.h:65
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 Has(const std::string &rEntry) const
This method checks if the Parameter contains a certain entry.
Definition: kratos_parameters.cpp:520
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
#define KRATOS_ERROR
Definition: exception.h:161
#define KRATOS_ERROR_IF_NOT(conditional)
Definition: exception.h: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
void block_for_each(TIterator itBegin, TIterator itEnd, TFunction &&rFunction)
Execute a functor on all items of a range in parallel.
Definition: parallel_utilities.h:299
model_part
Definition: face_heat.py:14
Definition: mesh_converter.cpp:38