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.
dam_input_table_nodal_young_modulus_process.hpp
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Lorenzo Gracia
11 //
12 //
13 
14 #if !defined(KRATOS_DAM_INPUT_TABLE_NODAL_YOUNG_MODULUS_PROCESS )
15 #define KRATOS_DAM_INPUT_TABLE_NODAL_YOUNG_MODULUS_PROCESS
16 
17 #include <cmath>
18 
19 // Project includes
20 #include "includes/kratos_flags.h"
22 #include "processes/process.h"
23 
24 // Application include
26 
27 namespace Kratos
28 {
29 
31 {
32 
33 public:
34 
36 
38 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
39 
42  Parameters& rParameters
43  ) : Process(Flags()) , mrModelPart(rModelPart) , mrTable(Table)
44  {
46 
47  //only include validation with c++11 since raw_literals do not exist in c++03
48  Parameters default_parameters( R"(
49  {
50  "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
51  "variable_name" : "PLEASE_PRESCRIBE_VARIABLE_NAME",
52  "initial_value" : 0.0,
53  "input_file_name" : "",
54  "interval":[
55  0.0,
56  0.0
57  ]
58  } )" );
59 
60  // Some values need to be mandatorily prescribed since no meaningful default value exist. For this reason try accessing to them
61  // So that an error is thrown if they don't exist
62  rParameters["variable_name"];
63  rParameters["model_part_name"];
64 
65  // Now validate agains defaults -- this also ensures no type mismatch
66  rParameters.ValidateAndAssignDefaults(default_parameters);
67 
68  mVariableName = rParameters["variable_name"].GetString();
69  mInitialValue = rParameters["initial_value"].GetDouble();
70  mInputFile = rParameters["input_file_name"].GetString();
71 
72  KRATOS_CATCH("");
73  }
74 
76 
79 
80 
81 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
82 
83  void Execute() override
84  {
85  }
86 
87  //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
88 
89  void ExecuteInitialize() override
90  {
91 
92  KRATOS_TRY;
93 
95  const int nnodes = mrModelPart.GetMesh(0).Nodes().size();
96 
97  if(nnodes != 0)
98  {
99  ModelPart::NodesContainerType::iterator it_begin = mrModelPart.GetMesh(0).NodesBegin();
100 
101  if ((mInputFile == "") || (mInputFile == "- No file") || (mInputFile == "- Add new file"))
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  it->FastGetSolutionStepValue(var) = mInitialValue;
109 
110  }
111  }
112  else
113  {
114  #pragma omp parallel for
115  for(int i = 0; i<nnodes; i++)
116  {
117  ModelPart::NodesContainerType::iterator it = it_begin + i;
118 
119  it->FastGetSolutionStepValue(var) = mrTable.GetValue(it->Id());
120 
121  }
122  }
123  }
124 
125  KRATOS_CATCH("");
126  }
127 
128 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
129 
131  {
132 
133  KRATOS_TRY;
134 
136  const int nnodes = mrModelPart.GetMesh(0).Nodes().size();
137 
138 
139  if(nnodes != 0)
140  {
141  ModelPart::NodesContainerType::iterator it_begin = mrModelPart.GetMesh(0).NodesBegin();
142 
143  if ((mInputFile == "") || (mInputFile == "- No file") || (mInputFile == "- Add new file"))
144  {
145  #pragma omp parallel for
146  for(int i = 0; i<nnodes; i++)
147  {
148  ModelPart::NodesContainerType::iterator it = it_begin + i;
149 
150  it->FastGetSolutionStepValue(var) = mInitialValue;
151 
152  }
153  }
154  else
155  {
156  #pragma omp parallel for
157  for(int i = 0; i<nnodes; i++)
158  {
159  ModelPart::NodesContainerType::iterator it = it_begin + i;
160 
161  it->FastGetSolutionStepValue(var) = mrTable.GetValue(it->Id());
162 
163  }
164  }
165  }
166 
167  KRATOS_CATCH("");
168  }
169 
170 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
171 
173  std::string Info() const override
174  {
175  return "DamInputTableNodalYoungModulusProcess";
176  }
177 
179  void PrintInfo(std::ostream& rOStream) const override
180  {
181  rOStream << "DamInputTableNodalYoungModulusProcess";
182  }
183 
185  void PrintData(std::ostream& rOStream) const override
186  {
187  }
188 
190 
191 protected:
192 
194 
197  std::string mVariableName;
199  std::string mInputFile;
200 
201 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
202 
203 private:
204 
207 
208 };//Class
209 
210 
212 inline std::istream& operator >> (std::istream& rIStream,
214 
216 inline std::ostream& operator << (std::ostream& rOStream,
218 {
219  rThis.PrintInfo(rOStream);
220  rOStream << std::endl;
221  rThis.PrintData(rOStream);
222 
223  return rOStream;
224 }
225 
226 } /* namespace Kratos.*/
227 
228 #endif /* KRATOS_DAM_INPUT_TABLE_NODAL_YOUNG_MODULUS_PROCESS defined */
229 
Definition: dam_input_table_nodal_young_modulus_process.hpp:31
ModelPart & mrModelPart
Member Variables.
Definition: dam_input_table_nodal_young_modulus_process.hpp:195
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: dam_input_table_nodal_young_modulus_process.hpp:83
TableType & mrTable
Definition: dam_input_table_nodal_young_modulus_process.hpp:196
std::string Info() const override
Turn back information as a string.
Definition: dam_input_table_nodal_young_modulus_process.hpp:173
KRATOS_CLASS_POINTER_DEFINITION(DamInputTableNodalYoungModulusProcess)
DamInputTableNodalYoungModulusProcess(ModelPart &rModelPart, TableType &Table, Parameters &rParameters)
Constructor.
Definition: dam_input_table_nodal_young_modulus_process.hpp:41
std::string mInputFile
Definition: dam_input_table_nodal_young_modulus_process.hpp:199
void ExecuteInitialize() override
This function is designed for being called at the beginning of the computations right after reading t...
Definition: dam_input_table_nodal_young_modulus_process.hpp:89
void ExecuteInitializeSolutionStep() override
This function will be executed at every time step BEFORE performing the solve phase.
Definition: dam_input_table_nodal_young_modulus_process.hpp:130
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: dam_input_table_nodal_young_modulus_process.hpp:185
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: dam_input_table_nodal_young_modulus_process.hpp:179
std::string mVariableName
Definition: dam_input_table_nodal_young_modulus_process.hpp:197
virtual ~DamInputTableNodalYoungModulusProcess()
Destructor.
Definition: dam_input_table_nodal_young_modulus_process.hpp:78
Table< double, double > TableType
Definition: dam_input_table_nodal_young_modulus_process.hpp:37
double mInitialValue
Definition: dam_input_table_nodal_young_modulus_process.hpp:198
Definition: flags.h:58
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
NodesContainerType & Nodes()
Definition: mesh.h:346
NodeIterator NodesBegin()
Definition: mesh.h:326
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
The base class for all processes in Kratos.
Definition: process.h:49
Definition: table.h:435
TResultType GetValue(TArgumentType const &X) const
Definition: table.h:502
This class represents the value of its variable depending to other variable.
Definition: piecewize_linear_table.h:63
#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
int nnodes
Definition: sensitivityMatrix.py:24
integer i
Definition: TensorModule.f:17