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_sum_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_SUM_METHOD_H_INCLUDED)
14 #define KRATOS_TEMPORAL_SUM_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>& rOutputVariable)
55  : TemporalMethod(rModelPart, EchoLevel),
56  mrInputVariable(rInputVariable),
57  mrOutputVariable(rOutputVariable)
58  {
59  }
60 
61  void CalculateStatistics() override
62  {
63  TContainerType& r_container =
64  MethodUtilities::GetDataContainer<TContainerType>(this->GetModelPart());
65 
66  const double delta_time = this->GetDeltaTime();
67 
68  const int number_of_items = r_container.size();
69 #pragma omp parallel for
70  for (int i = 0; i < number_of_items; ++i)
71  {
72  TContainerItemType& r_item = *(r_container.begin() + i);
73  const TDataType& r_input_value =
74  TDataRetrievalFunctor<TContainerItemType>()(r_item, mrInputVariable);
75  TDataType& r_output_value =
76  TDataStorageFunctor<TContainerItemType>()(r_item, mrOutputVariable);
77  MethodUtilities::DataTypeSizeChecker(r_input_value, r_output_value);
78 
79  TemporalSumMethod::CalculateSum<TDataType>(
80  r_output_value, r_input_value, delta_time);
81  }
82 
83  KRATOS_INFO_IF("TemporalValueSumMethod", this->GetEchoLevel() > 1)
84  << "Calculated temporal value sum for " << mrInputVariable.Name()
85  << " input variable with " << mrOutputVariable.Name()
86  << " output variable for " << this->GetModelPart().Name() << ".\n";
87  }
88 
90  {
91  TContainerType& r_container =
92  MethodUtilities::GetDataContainer<TContainerType>(this->GetModelPart());
93 
94  auto& initializer_method =
95  TemporalMethodUtilities::InitializeVariables<TContainerType, TContainerItemType, TDataRetrievalFunctor, TDataStorageFunctor, TDataType>;
96  initializer_method(r_container, mrOutputVariable, mrInputVariable);
97 
98  KRATOS_INFO_IF("TemporalValueSumMethod", this->GetEchoLevel() > 0)
99  << "Initialized temporal value sum method for "
100  << mrInputVariable.Name() << " input variable with "
101  << mrOutputVariable.Name() << " output variable for "
102  << this->GetModelPart().Name() << ".\n";
103  }
104 
105  private:
106  const Variable<TDataType>& mrInputVariable;
107  const Variable<TDataType>& mrOutputVariable;
108  };
109 
110  template <class TDataType>
111  class NormMethod : public TemporalMethod
112  {
113  public:
115 
117  ModelPart& rModelPart,
118  const std::string& rNormType,
119  const Variable<TDataType>& rInputVariable,
120  const int EchoLevel,
121  const Variable<double>& rOutputVariable)
122  : TemporalMethod(rModelPart, EchoLevel),
123  mNormType(rNormType),
124  mrInputVariable(rInputVariable),
125  mrOutputVariable(rOutputVariable)
126  {
127  }
128 
129  void CalculateStatistics() override
130  {
131  TContainerType& r_container =
132  MethodUtilities::GetDataContainer<TContainerType>(this->GetModelPart());
133 
134  const auto& norm_method =
135  MethodUtilities::GetNormMethod(mrInputVariable, mNormType);
136 
137  const double delta_time = this->GetDeltaTime();
138 
139  const int number_of_items = r_container.size();
140 #pragma omp parallel for
141  for (int i = 0; i < number_of_items; ++i)
142  {
143  TContainerItemType& r_item = *(r_container.begin() + i);
144  const TDataType& r_input_value =
145  TDataRetrievalFunctor<TContainerItemType>()(r_item, mrInputVariable);
146  const double input_norm_value = norm_method(r_input_value);
147  double& r_output_value =
148  TDataStorageFunctor<TContainerItemType>()(r_item, mrOutputVariable);
149 
150  TemporalSumMethod::CalculateSum<double>(
151  r_output_value, input_norm_value, delta_time);
152  }
153 
154  KRATOS_INFO_IF("TemporalNormSumMethod", this->GetEchoLevel() > 1)
155  << "Calculated temporal norm sum for " << mrInputVariable.Name()
156  << " input variable with " << mrOutputVariable.Name()
157  << " output variable for " << this->GetModelPart().Name() << ".\n";
158  }
159 
161  {
162  TContainerType& r_container =
163  MethodUtilities::GetDataContainer<TContainerType>(this->GetModelPart());
164 
165  auto& initializer_method =
166  TemporalMethodUtilities::InitializeVariables<TContainerType, TContainerItemType, TDataStorageFunctor>;
167  initializer_method(r_container, mrOutputVariable, 0.0);
168 
169  KRATOS_INFO_IF("TemporalNormSumMethod", this->GetEchoLevel() > 0)
170  << "Initialized temporal norm sum method for "
171  << mrInputVariable.Name() << " input variable with "
172  << mrOutputVariable.Name() << " output variable for "
173  << this->GetModelPart().Name() << ".\n";
174  }
175 
176  private:
177  const std::string mNormType;
178  const Variable<TDataType>& mrInputVariable;
179  const Variable<double>& mrOutputVariable;
180  };
181 
182  std::vector<TemporalMethod::Pointer> static CreateTemporalMethodObject(
183  ModelPart& rModelPart, const std::string& rNormType, const int EchoLevel, Parameters Params)
184  {
185  KRATOS_TRY
186 
187  Parameters default_parameters = Parameters(R"(
188  {
189  "input_variables" : [],
190  "output_variables" : []
191  })");
192  Params.RecursivelyValidateAndAssignDefaults(default_parameters);
193 
194  const std::vector<std::string>& input_variable_names_list =
195  Params["input_variables"].GetStringArray();
196  const std::vector<std::string>& output_variable_names_list =
197  Params["output_variables"].GetStringArray();
198 
199  std::vector<TemporalMethod::Pointer> method_list;
200  if (rNormType == "none") // for non norm types
201  {
203  input_variable_names_list, output_variable_names_list);
204  const int number_of_variables = input_variable_names_list.size();
205  for (int i = 0; i < number_of_variables; ++i)
206  {
207  const std::string& r_variable_input_name = input_variable_names_list[i];
208  const std::string& r_variable_output_name =
209  output_variable_names_list[i];
211  rModelPart, rNormType, r_variable_input_name, EchoLevel,
212  r_variable_output_name, method_list, ValueMethod)
213  }
214  }
215  else // for values with norms
216  {
217  MethodUtilities::CheckVariableType<double>(output_variable_names_list);
218 
219  const int number_of_variables = input_variable_names_list.size();
220  for (int i = 0; i < number_of_variables; ++i)
221  {
222  const std::string& r_variable_input_name = input_variable_names_list[i];
223  const std::string& r_variable_output_name =
224  output_variable_names_list[i];
226  rModelPart, rNormType, r_variable_input_name, EchoLevel,
227  r_variable_output_name, method_list, NormMethod)
228  }
229  }
230 
231  return method_list;
232 
233  KRATOS_CATCH("");
234  }
235 
236 private:
237  template <class TDataType>
238  void static CalculateSum(TDataType& rSum, const TDataType& rNewDataPoint, const double DeltaTime)
239  {
240  rSum = (rSum + rNewDataPoint * DeltaTime);
241  }
242 };
243 } // namespace TemporalMethods
244 } // namespace Kratos
245 
246 #endif // KRATOS_TEMPORAL_SUM_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_sum_method.h:112
Definition: temporal_sum_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_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
static std::vector< TemporalMethod::Pointer > CreateTemporalMethodObject(ModelPart &rModelPart, const std::string &rNormType, const int EchoLevel, Parameters Params)
Definition: temporal_sum_method.h:182
void CalculateStatistics() override
Definition: temporal_sum_method.h:61
void InitializeStatisticsVariables() override
Definition: temporal_sum_method.h:160
void CalculateStatistics() override
Definition: temporal_sum_method.h:129
void InitializeStatisticsVariables() override
Definition: temporal_sum_method.h:89
ValueMethod(ModelPart &rModelPart, const std::string &rNormType, const Variable< TDataType > &rInputVariable, const int EchoLevel, const Variable< TDataType > &rOutputVariable)
Definition: temporal_sum_method.h:49
NormMethod(ModelPart &rModelPart, const std::string &rNormType, const Variable< TDataType > &rInputVariable, const int EchoLevel, const Variable< double > &rOutputVariable)
Definition: temporal_sum_method.h:116
#define ADD_TEMPORAL_NORM_METHOD_ONE_OUTPUT_VARIABLE_OBJECT( model_part, norm_type, input_variable, echo_level, output_variable, object_list, method)
Definition: method_utilities.h:77
#define ADD_TEMPORAL_VALUE_METHOD_ONE_OUTPUT_VARIABLE_OBJECT( model_part, norm_type, input_variable, echo_level, output_variable, object_list, method)
Definition: method_utilities.h:27
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