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_random_fields_variable_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: Joaquín Irazábal González (jirazabal@cimne.upc.edu)
11 //
12 //
13 
14 #if !defined(KRATOS_DAM_RANDOM_FIELDS_VARIABLE_PROCESS )
15 #define KRATOS_DAM_RANDOM_FIELDS_VARIABLE_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  "mean_value" : 0.0,
53  "min_value" : 0.0,
54  "max_value" : 0.0,
55  "variance" : 0.0,
56  "corr_length" : 0,
57  "interval":[
58  0.0,
59  0.0
60  ]
61  } )" );
62 
63  // Some values need to be mandatorily prescribed since no meaningful default value exist. For this reason try accessing to them
64  // So that an error is thrown if they don't exist
65  rParameters["variable_name"];
66  rParameters["model_part_name"];
67 
68  // Now validate agains defaults -- this also ensures no type mismatch
69  rParameters.ValidateAndAssignDefaults(default_parameters);
70 
71  mVariableName = rParameters["variable_name"].GetString();
72 
73  KRATOS_CATCH("");
74  }
75 
77 
80 
81 
82 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
83 
84  void Execute() override
85  {
86  }
87 
88  //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
89 
90  void ExecuteInitialize() override
91  {
92 
93  KRATOS_TRY;
94 
96  const int nnodes = mrModelPart.GetMesh(0).Nodes().size();
97 
98  if(nnodes != 0)
99  {
100  ModelPart::NodesContainerType::iterator it_begin = mrModelPart.GetMesh(0).NodesBegin();
101 
102  #pragma omp parallel for
103  for(int i = 0; i<nnodes; i++)
104  {
105  ModelPart::NodesContainerType::iterator it = it_begin + i;
106 
107  it->FastGetSolutionStepValue(var) = mrTable.GetValue(it->Id());
108 
109  }
110  }
111 
112  KRATOS_CATCH("");
113  }
114 
115 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
116 
118  {
119  }
120 
121 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
122 
124  std::string Info() const override
125  {
126  return "DamRandomFieldsVariableProcess";
127  }
128 
130  void PrintInfo(std::ostream& rOStream) const override
131  {
132  rOStream << "DamRandomFieldsVariableProcess";
133  }
134 
136  void PrintData(std::ostream& rOStream) const override
137  {
138  }
139 
141 
142 protected:
143 
145 
148  std::string mVariableName;
149 
150 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
151 
152 private:
153 
156 
157 };//Class
158 
159 
161 inline std::istream& operator >> (std::istream& rIStream,
163 
165 inline std::ostream& operator << (std::ostream& rOStream,
166  const DamRandomFieldsVariableProcess& rThis)
167 {
168  rThis.PrintInfo(rOStream);
169  rOStream << std::endl;
170  rThis.PrintData(rOStream);
171 
172  return rOStream;
173 }
174 
175 } /* namespace Kratos.*/
176 
177 #endif /* KRATOS_DAM_RANDOM_FIELDS_VARIABLE_PROCESS defined */
178 
Definition: dam_random_fields_variable_process.hpp:31
virtual ~DamRandomFieldsVariableProcess()
Destructor.
Definition: dam_random_fields_variable_process.hpp:79
ModelPart & mrModelPart
Member Variables.
Definition: dam_random_fields_variable_process.hpp:146
std::string Info() const override
Turn back information as a string.
Definition: dam_random_fields_variable_process.hpp:124
void ExecuteInitialize() override
This function is designed for being called at the beginning of the computations right after reading t...
Definition: dam_random_fields_variable_process.hpp:90
TableType & mrTable
Definition: dam_random_fields_variable_process.hpp:147
void ExecuteInitializeSolutionStep() override
This function will be executed at every time step BEFORE performing the solve phase.
Definition: dam_random_fields_variable_process.hpp:117
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: dam_random_fields_variable_process.hpp:84
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: dam_random_fields_variable_process.hpp:136
DamRandomFieldsVariableProcess(ModelPart &rModelPart, TableType &Table, Parameters &rParameters)
Constructor.
Definition: dam_random_fields_variable_process.hpp:41
KRATOS_CLASS_POINTER_DEFINITION(DamRandomFieldsVariableProcess)
std::string mVariableName
Definition: dam_random_fields_variable_process.hpp:148
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: dam_random_fields_variable_process.hpp:130
Table< double, double > TableType
Definition: dam_random_fields_variable_process.hpp:37
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
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