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.
temporal_variance_method.h
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: Suneth Warnakulasuriya (https://github.com/sunethwarna)
11 //
12 
13 #if !defined(KRATOS_TEMPORAL_VARIANCE_METHOD_H_INCLUDED)
14 #define KRATOS_TEMPORAL_VARIANCE_METHOD_H_INCLUDED
15 
16 // System includes
17 
18 // External includes
19 
20 // Project includes
21 #include "includes/define.h"
22 #include "includes/model_part.h"
23 
24 // Application includes
28 
29 namespace Kratos
30 {
33 
36 
37 namespace TemporalMethods
38 {
39 template <class TContainerType, class TContainerItemType, template <class T> class TDataRetrievalFunctor, template <class T> class TDataStorageFunctor>
41 {
42 public:
43  template <class TDataType>
44  class ValueMethod : public TemporalMethod
45  {
46  public:
48 
50  ModelPart& rModelPart,
51  const std::string& rNormType,
52  const Variable<TDataType>& rInputVariable,
53  const int EchoLevel,
54  const Variable<TDataType>& rOutputMeanVariable,
55  const Variable<TDataType>& rOutputVarianceVariable)
56  : TemporalMethod(rModelPart, EchoLevel),
57  mrInputVariable(rInputVariable),
58  mrOutputMeanVariable(rOutputMeanVariable),
59  mrOutputVarianceVariable(rOutputVarianceVariable)
60  {
62 
63  KRATOS_ERROR_IF(rOutputMeanVariable == rOutputVarianceVariable) << "Same variable is given for mean and variance in value variance method with input variable "
64  << rInputVariable
65  .Name()
66  << ". Please provide two different variables. [ variable = "
67  << rOutputMeanVariable
68  .Name()
69  << " ].\n";
70 
71  KRATOS_CATCH("");
72  }
73 
74  void CalculateStatistics() override
75  {
76  TContainerType& r_container =
77  MethodUtilities::GetDataContainer<TContainerType>(this->GetModelPart());
78 
79  const double delta_time = this->GetDeltaTime();
80  const double old_total_time = this->GetTotalTime();
81  const double total_time = old_total_time + delta_time;
82 
83  const int number_of_items = r_container.size();
84 #pragma omp parallel for
85  for (int i = 0; i < number_of_items; ++i)
86  {
87  TContainerItemType& r_item = *(r_container.begin() + i);
88  const TDataType& r_input_value =
89  TDataRetrievalFunctor<TContainerItemType>()(r_item, mrInputVariable);
90  TDataType& r_output_mean_value =
91  TDataStorageFunctor<TContainerItemType>()(r_item, mrOutputMeanVariable);
92  TDataType& r_output_variance_value =
93  TDataStorageFunctor<TContainerItemType>()(r_item, mrOutputVarianceVariable);
94 
95  MethodUtilities::DataTypeSizeChecker(r_input_value, r_output_mean_value);
96  MethodUtilities::DataTypeSizeChecker(r_input_value, r_output_variance_value);
97 
98  TemporalVarianceMethod::CalculateMeanAndVariance<TDataType>(
99  r_output_mean_value, r_output_variance_value, r_input_value,
100  delta_time, old_total_time, total_time);
101  }
102 
103  KRATOS_INFO_IF("TemporalValueVarianceMethod", this->GetEchoLevel() > 1)
104  << "Calculated temporal value variance for "
105  << mrInputVariable.Name() << " input variable with "
106  << mrOutputMeanVariable.Name() << " mean variable and "
107  << mrOutputVarianceVariable.Name() << " variance variable for "
108  << this->GetModelPart().Name() << ".\n";
109  }
110 
112  {
113  TContainerType& r_container =
114  MethodUtilities::GetDataContainer<TContainerType>(this->GetModelPart());
115 
116  auto& initializer_method =
117  TemporalMethodUtilities::InitializeVariables<TContainerType, TContainerItemType, TDataRetrievalFunctor, TDataStorageFunctor, TDataType>;
118  initializer_method(r_container, mrOutputMeanVariable, mrInputVariable);
119  initializer_method(r_container, mrOutputVarianceVariable, mrInputVariable);
120 
121  KRATOS_INFO_IF("TemporalValueVarianceMethod", this->GetEchoLevel() > 0)
122  << "Initialized temporal value variance method for "
123  << mrInputVariable.Name() << " input variable with "
124  << mrOutputMeanVariable.Name() << " mean variable and "
125  << mrOutputVarianceVariable.Name() << " variance variable for "
126  << this->GetModelPart().Name() << ".\n";
127  }
128 
129  private:
130  const Variable<TDataType>& mrInputVariable;
131  const Variable<TDataType>& mrOutputMeanVariable;
132  const Variable<TDataType>& mrOutputVarianceVariable;
133  };
134 
135  template <class TDataType>
136  class NormMethod : public TemporalMethod
137  {
138  public:
140 
142  ModelPart& rModelPart,
143  const std::string& rNormType,
144  const Variable<TDataType>& rInputVariable,
145  const int EchoLevel,
146  const Variable<double>& rOutputMeanVariable,
147  const Variable<double>& rOutputVarianceVariable)
148  : TemporalMethod(rModelPart, EchoLevel),
149  mNormType(rNormType),
150  mrInputVariable(rInputVariable),
151  mrOutputMeanVariable(rOutputMeanVariable),
152  mrOutputVarianceVariable(rOutputVarianceVariable)
153  {
154  KRATOS_TRY
155 
156  KRATOS_ERROR_IF(rOutputMeanVariable == rOutputVarianceVariable) << "Same variable is given for mean and variance in norm variance method with input variable "
157  << rInputVariable
158  .Name()
159  << ". Please provide two different variables. [ variable = "
160  << rOutputMeanVariable
161  .Name()
162  << " ].\n";
163 
164  KRATOS_CATCH("");
165  }
166 
167  void CalculateStatistics() override
168  {
169  TContainerType& r_container =
170  MethodUtilities::GetDataContainer<TContainerType>(this->GetModelPart());
171 
172  const auto& norm_method =
173  MethodUtilities::GetNormMethod(mrInputVariable, mNormType);
174 
175  const double delta_time = this->GetDeltaTime();
176  const double old_total_time = this->GetTotalTime();
177  const double total_time = old_total_time + delta_time;
178 
179  const int number_of_items = r_container.size();
180 #pragma omp parallel for
181  for (int i = 0; i < number_of_items; ++i)
182  {
183  TContainerItemType& r_item = *(r_container.begin() + i);
184  const TDataType& r_input_value =
185  TDataRetrievalFunctor<TContainerItemType>()(r_item, mrInputVariable);
186  const double input_norm_value = norm_method(r_input_value);
187  double& r_output_mean_value =
188  TDataStorageFunctor<TContainerItemType>()(r_item, mrOutputMeanVariable);
189  double& r_output_variance_value =
190  TDataStorageFunctor<TContainerItemType>()(r_item, mrOutputVarianceVariable);
191 
192  TemporalVarianceMethod::CalculateMeanAndVariance<double>(
193  r_output_mean_value, r_output_variance_value,
194  input_norm_value, delta_time, old_total_time, total_time);
195  }
196 
197  KRATOS_INFO_IF("TemporalNormVarianceMethod", this->GetEchoLevel() > 1)
198  << "Calculated temporal norm variance for " << mrInputVariable.Name()
199  << " input variable with " << mrOutputMeanVariable.Name()
200  << " mean variable and " << mrOutputVarianceVariable.Name()
201  << " variance variable for " << this->GetModelPart().Name() << ".\n";
202  }
203 
204  // norm output variable initialization
206  {
207  TContainerType& r_container =
208  MethodUtilities::GetDataContainer<TContainerType>(this->GetModelPart());
209 
210  auto& initializer_method =
211  TemporalMethodUtilities::InitializeVariables<TContainerType, TContainerItemType, TDataStorageFunctor>;
212  initializer_method(r_container, mrOutputMeanVariable, 0.0);
213  initializer_method(r_container, mrOutputVarianceVariable, 0.0);
214 
215  KRATOS_INFO_IF("TemporalNormVarianceMethod", this->GetEchoLevel() > 0)
216  << "Initialized temporal norm variance method for "
217  << mrInputVariable.Name() << " input variable with "
218  << mrOutputMeanVariable.Name() << " mean variable and "
219  << mrOutputVarianceVariable.Name() << " variance variable for "
220  << this->GetModelPart().Name() << ".\n";
221  }
222 
223  private:
224  const std::string mNormType;
225  const Variable<TDataType>& mrInputVariable;
226  const Variable<double>& mrOutputMeanVariable;
227  const Variable<double>& mrOutputVarianceVariable;
228  };
229 
230  std::vector<TemporalMethod::Pointer> static CreateTemporalMethodObject(
231  ModelPart& rModelPart, const std::string& rNormType, const int EchoLevel, Parameters Params)
232  {
233  KRATOS_TRY
234 
235  Parameters default_parameters = Parameters(R"(
236  {
237  "input_variables" : [],
238  "output_mean_variables" : [],
239  "output_variance_variables" : []
240  })");
241  Params.RecursivelyValidateAndAssignDefaults(default_parameters);
242 
243  const std::vector<std::string>& input_variable_names_list =
244  Params["input_variables"].GetStringArray();
245  const std::vector<std::string>& output_variable_1_names_list =
246  Params["output_mean_variables"].GetStringArray();
247  const std::vector<std::string>& output_variable_2_names_list =
248  Params["output_variance_variables"].GetStringArray();
249 
250  std::vector<TemporalMethod::Pointer> method_list;
251  if (rNormType == "none") // for non norm types
252  {
254  input_variable_names_list, output_variable_1_names_list);
256  input_variable_names_list, output_variable_2_names_list);
257 
258  const int number_of_variables = input_variable_names_list.size();
259  for (int i = 0; i < number_of_variables; ++i)
260  {
261  const std::string& r_variable_input_name = input_variable_names_list[i];
262  const std::string& r_variable_1_output_name =
263  output_variable_1_names_list[i];
264  const std::string& r_variable_2_output_name =
265  output_variable_2_names_list[i];
267  rModelPart, rNormType, r_variable_input_name, EchoLevel,
268  r_variable_1_output_name, r_variable_2_output_name, method_list, ValueMethod)
269  }
270  }
271  else // for values with norms
272  {
273  MethodUtilities::CheckVariableType<double>(output_variable_1_names_list);
274  MethodUtilities::CheckVariableType<double>(output_variable_2_names_list);
275 
276  const int number_of_variables = input_variable_names_list.size();
277  for (int i = 0; i < number_of_variables; ++i)
278  {
279  const std::string& r_variable_input_name = input_variable_names_list[i];
280  const std::string& r_variable_1_output_name =
281  output_variable_1_names_list[i];
282  const std::string& r_variable_2_output_name =
283  output_variable_2_names_list[i];
285  rModelPart, rNormType, r_variable_input_name, EchoLevel,
286  r_variable_1_output_name, r_variable_2_output_name, method_list, NormMethod)
287  }
288  }
289 
290  return method_list;
291 
292  KRATOS_CATCH("");
293  }
294 
295 private:
296  template <class TDataType>
297  void static CalculateMeanAndVariance(
298  TDataType& rMean,
299  TDataType& rVariance,
300  const TDataType& rNewDataPoint,
301  const double DeltaTime,
302  const double OldTotalTime,
303  const double CurrentTotalTime)
304  {
305  const TDataType new_mean =
306  (rMean * OldTotalTime + rNewDataPoint * DeltaTime) * (1.0 / CurrentTotalTime);
307  rVariance =
308  ((rVariance + MethodUtilities::RaiseToPower<TDataType>(rMean, 2)) * OldTotalTime +
309  MethodUtilities::RaiseToPower<TDataType>(rNewDataPoint, 2) * DeltaTime) *
310  (1 / CurrentTotalTime) -
311  MethodUtilities::RaiseToPower<TDataType>(new_mean, 2);
312  rMean = new_mean;
313  }
314 };
315 } // namespace TemporalMethods
316 } // namespace Kratos
317 
318 #endif // KRATOS_TEMPORAL_VARIANCE_METHOD_H_INCLUDED
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
std::string & Name()
Definition: model_part.h:1811
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
std::vector< std::string > GetStringArray() const
This method returns the array of strings in the current Parameter.
Definition: kratos_parameters.cpp:693
void RecursivelyValidateAndAssignDefaults(const Parameters &rDefaultParameters)
This function is designed to verify that the parameters under testing match the form prescribed by th...
Definition: kratos_parameters.cpp:1389
Definition: temporal_method.h:40
Definition: temporal_variance_method.h:137
Definition: temporal_variance_method.h:45
Definition: temporal_variance_method.h:41
const std::string & Name() const
Definition: variable_data.h:201
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR_IF(conditional)
Definition: exception.h:162
#define KRATOS_INFO_IF(label, conditional)
Definition: logger.h:251
double GetDeltaTime() const
Definition: temporal_method.h:106
int GetEchoLevel() const
Definition: temporal_method.h:112
ModelPart & GetModelPart() const
Definition: temporal_method.h:82
double GetTotalTime() const
Definition: temporal_method.h:87
void InitializeStatisticsVariables() override
Definition: temporal_variance_method.h:205
void CalculateStatistics() override
Definition: temporal_variance_method.h:167
NormMethod(ModelPart &rModelPart, const std::string &rNormType, const Variable< TDataType > &rInputVariable, const int EchoLevel, const Variable< double > &rOutputMeanVariable, const Variable< double > &rOutputVarianceVariable)
Definition: temporal_variance_method.h:141
void CalculateStatistics() override
Definition: temporal_variance_method.h:74
static std::vector< TemporalMethod::Pointer > CreateTemporalMethodObject(ModelPart &rModelPart, const std::string &rNormType, const int EchoLevel, Parameters Params)
Definition: temporal_variance_method.h:230
ValueMethod(ModelPart &rModelPart, const std::string &rNormType, const Variable< TDataType > &rInputVariable, const int EchoLevel, const Variable< TDataType > &rOutputMeanVariable, const Variable< TDataType > &rOutputVarianceVariable)
Definition: temporal_variance_method.h:49
void InitializeStatisticsVariables() override
Definition: temporal_variance_method.h:111
#define ADD_TEMPORAL_NORM_METHOD_TWO_OUTPUT_VARIABLE_OBJECT( model_part, norm_type, input_variable, echo_level, output_variable_1, output_variable_2, object_list, method)
Definition: method_utilities.h:190
#define ADD_TEMPORAL_VALUE_METHOD_TWO_OUTPUT_VARIABLE_OBJECT( model_part, norm_type, input_variable, echo_level, output_variable_1, output_variable_2, object_list, method)
Definition: method_utilities.h:127
static int EchoLevel
Definition: co_sim_EMPIRE_API.h:42
const std::function< double(const TDataType &)> GetNormMethod(const Variable< TDataType > &rVariable, const std::string &rNormType)
Definition: method_utilities.cpp:313
void DataTypeSizeChecker(const TDataType &rData, const TDataType &rReferenceData)
Definition: method_utilities.cpp:131
void CheckInputOutputVariables(const std::vector< std::string > &rInputVariableNamesList, const std::vector< std::string > &rOutputVariableNamesList)
Definition: method_utilities.cpp:678
template void CheckVariableType< double >(const std::vector< std::string > &rVariableNamesList)
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
delta_time
Definition: generate_frictional_mortar_condition.py:130
integer i
Definition: TensorModule.f:17