13 #if !defined(KRATOS_TEMPORAL_RMS_METHOD_H_INCLUDED)
14 #define KRATOS_TEMPORAL_RMS_METHOD_H_INCLUDED
37 namespace TemporalMethods
39 template <
class TContainerType,
class TContainerItemType,
template <
class T>
class TDataRetrievalFunctor,
template <
class T>
class TDataStorageFunctor>
43 template <
class TDataType>
51 const std::string& rNormType,
56 mrInputVariable(rInputVariable),
57 mrOutputVariable(rOutputVariable)
63 TContainerType& r_container =
64 MethodUtilities::GetDataContainer<TContainerType>(this->
GetModelPart());
68 const double total_time = old_total_time +
delta_time;
70 const int number_of_items = r_container.size();
71 #pragma omp parallel for
72 for (
int i = 0;
i < number_of_items; ++
i)
74 TContainerItemType& r_item = *(r_container.begin() +
i);
75 const TDataType& r_input_value =
76 TDataRetrievalFunctor<TContainerItemType>()(r_item, mrInputVariable);
77 TDataType& r_output_value =
78 TDataStorageFunctor<TContainerItemType>()(r_item, mrOutputVariable);
81 TemporalRootMeanSquareMethod::CalculateRootMeanSquare<TDataType>(
82 r_output_value, r_input_value,
delta_time, old_total_time, total_time);
86 <<
"Calculated temporal value root mean square for "
87 << mrInputVariable.Name() <<
" input variable with "
88 << mrOutputVariable.Name() <<
" root mean square variable for "
94 TContainerType& r_container =
95 MethodUtilities::GetDataContainer<TContainerType>(this->
GetModelPart());
97 auto& initializer_method =
98 TemporalMethodUtilities::InitializeVariables<TContainerType, TContainerItemType, TDataRetrievalFunctor, TDataStorageFunctor, TDataType>;
99 initializer_method(r_container, mrOutputVariable, mrInputVariable);
102 <<
"Initialized temporal value root mean square method for "
103 << mrInputVariable.Name() <<
" input variable with "
104 << mrOutputVariable.Name() <<
" root mean square variable for "
113 template <
class TDataType>
121 const std::string& rNormType,
126 mNormType(rNormType),
127 mrInputVariable(rInputVariable),
128 mrOutputVariable(rOutputVariable)
134 TContainerType& r_container =
135 MethodUtilities::GetDataContainer<TContainerType>(this->
GetModelPart());
137 const auto& norm_method =
142 const double total_time = old_total_time +
delta_time;
144 const int number_of_items = r_container.size();
145 #pragma omp parallel for
146 for (
int i = 0;
i < number_of_items; ++
i)
148 TContainerItemType& r_item = *(r_container.begin() +
i);
149 const TDataType& r_input_value =
150 TDataRetrievalFunctor<TContainerItemType>()(r_item, mrInputVariable);
151 const double input_norm_value = norm_method(r_input_value);
152 double& r_output_value =
153 TDataStorageFunctor<TContainerItemType>()(r_item, mrOutputVariable);
155 TemporalRootMeanSquareMethod::CalculateRootMeanSquare<double>(
156 r_output_value, input_norm_value,
delta_time, old_total_time, total_time);
160 <<
"Calculated temporal norm root mean square for "
161 << mrInputVariable.Name() <<
" input variable with "
162 << mrOutputVariable.
Name() <<
" root mean square variable for "
168 TContainerType& r_container =
169 MethodUtilities::GetDataContainer<TContainerType>(this->
GetModelPart());
171 auto& initializer_method =
172 TemporalMethodUtilities::InitializeVariables<TContainerType, TContainerItemType, TDataStorageFunctor>;
173 initializer_method(r_container, mrOutputVariable, 0.0);
176 <<
"Initialized temporal norm root mean square method for "
177 << mrInputVariable.Name() <<
" input variable with "
178 << mrOutputVariable.
Name() <<
" root mean square variable for "
183 const std::string mNormType;
195 "input_variables" : [],
196 "output_variables" : []
200 const std::vector<std::string>& input_variable_names_list =
202 const std::vector<std::string>& output_variable_names_list =
205 std::vector<TemporalMethod::Pointer> method_list;
206 if (rNormType ==
"none")
209 input_variable_names_list, output_variable_names_list);
210 const int number_of_variables = input_variable_names_list.size();
211 for (
int i = 0;
i < number_of_variables; ++
i)
213 const std::string& r_variable_input_name = input_variable_names_list[
i];
214 const std::string& r_variable_output_name =
215 output_variable_names_list[
i];
217 rModelPart, rNormType, r_variable_input_name,
EchoLevel,
225 const int number_of_variables = input_variable_names_list.size();
226 for (
int i = 0;
i < number_of_variables; ++
i)
228 const std::string& r_variable_input_name = input_variable_names_list[
i];
229 const std::string& r_variable_output_name =
230 output_variable_names_list[
i];
232 rModelPart, rNormType, r_variable_input_name,
EchoLevel,
233 r_variable_output_name, method_list,
NormMethod)
243 template <
class TDataType>
244 void static CalculateRootMeanSquare(
245 TDataType& rRootMeanSquare,
246 const TDataType& rNewDataPoint,
247 const double DeltaTime,
248 const double OldTotalTime,
249 const double CurrentTotalTime)
251 rRootMeanSquare = MethodUtilities::RaiseToPower<TDataType>(
252 (MethodUtilities::RaiseToPower<TDataType>(rRootMeanSquare, 2) * OldTotalTime +
254 (1.0 / CurrentTotalTime),
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_rms_method.h:115
Definition: temporal_rms_method.h:45
Definition: temporal_rms_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
double GetTotalTime() const
Definition: temporal_method.h:87
KRATOS_CLASS_POINTER_DEFINITION(NormMethod)
static std::vector< TemporalMethod::Pointer > CreateTemporalMethodObject(ModelPart &rModelPart, const std::string &rNormType, const int EchoLevel, Parameters Params)
Definition: temporal_rms_method.h:188
ValueMethod(ModelPart &rModelPart, const std::string &rNormType, const Variable< TDataType > &rInputVariable, const int EchoLevel, const Variable< TDataType > &rOutputVariable)
Definition: temporal_rms_method.h:49
void CalculateStatistics() override
Definition: temporal_rms_method.h:132
void CalculateStatistics() override
Definition: temporal_rms_method.h:61
KRATOS_CLASS_POINTER_DEFINITION(ValueMethod)
void InitializeStatisticsVariables() override
Definition: temporal_rms_method.h:166
void InitializeStatisticsVariables() override
Definition: temporal_rms_method.h:92
NormMethod(ModelPart &rModelPart, const std::string &rNormType, const Variable< TDataType > &rInputVariable, const int EchoLevel, const Variable< double > &rOutputVariable)
Definition: temporal_rms_method.h:119
#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
TDataType RaiseToPower(const TDataType &rData, const double Power)
Definition: method_utilities.cpp:34
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