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.
spatial_methods.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_SPATIAL_VARIANCE_METHOD_H_INCLUDED)
14 #define KRATOS_SPATIAL_VARIANCE_METHOD_H_INCLUDED
15 
16 // System includes
17 #include <algorithm>
18 #include <cmath>
19 #include <functional>
20 #include <numeric>
21 #include <tuple>
22 #include <vector>
23 
24 // External includes
25 
26 // Project includes
27 #include "includes/communicator.h"
28 #include "includes/define.h"
29 #include "includes/model_part.h"
30 
31 // Application includes
33 
34 namespace Kratos
35 {
38 
41 
42 namespace SpatialMethods
43 {
44 template <class TContainerType, class TContainerItemType, template <class T> class TDataRetrievalFunctor>
46 {
47 public:
48  // special overloaded method for flags
49  int static CalculateSum(const ModelPart& rModelPart, const Flags& rVariable)
50  {
51  const TContainerType& r_container =
52  MethodUtilities::GetLocalDataContainer<TContainerType>(rModelPart);
53  int sum = 0;
54 
55 #pragma omp parallel for reduction(+ : sum)
56  for (int i = 0; i < static_cast<int>(r_container.size()); ++i)
57  {
58  const TContainerItemType& r_item = *(r_container.begin() + i);
59  sum += r_item.Is(rVariable);
60  }
61 
62  sum = rModelPart.GetCommunicator().GetDataCommunicator().SumAll(sum);
63  return sum;
64  }
65 
66  template <class TDataType>
67  TDataType static CalculateSum(const ModelPart& rModelPart, const Variable<TDataType>& rVariable)
68  {
70 
71  const TContainerType& r_container =
72  MethodUtilities::GetLocalDataContainer<TContainerType>(rModelPart);
73 
74  TDataType global_sum = rVariable.Zero();
75 
76  if (r_container.size() > 0)
77  {
78  const TDataType& r_initial_value =
79  TDataRetrievalFunctor<TContainerItemType>()(*r_container.begin(), rVariable);
80 
81  MethodUtilities::DataTypeSizeInitializer(global_sum, r_initial_value);
82 
83 #pragma omp parallel
84  {
85  TDataType sum = rVariable.Zero();
87 
88 #pragma omp for
89  for (int i = 0; i < static_cast<int>(r_container.size()); ++i)
90  {
91  const TContainerItemType& r_item = *(r_container.begin() + i);
92  const TDataType& current_value =
93  TDataRetrievalFunctor<TContainerItemType>()(r_item, rVariable);
95  sum += current_value;
96  }
97 #pragma omp critical
98  {
99  global_sum += sum;
100  }
101  }
102  }
103 
104  global_sum = rModelPart.GetCommunicator().GetDataCommunicator().SumAll(global_sum);
105 
106  return global_sum;
107 
108  KRATOS_CATCH("");
109  }
110 
111  template <class TDataType>
112  double static CalculateNormSum(
113  const ModelPart& rModelPart,
114  const Variable<TDataType>& rVariable,
115  const std::string& rNormType,
116  Parameters Params)
117  {
118  KRATOS_TRY
119 
120  const TContainerType& r_container =
121  MethodUtilities::GetLocalDataContainer<TContainerType>(rModelPart);
122 
123  double global_sum = 0.0;
124  const auto& norm_method =
125  MethodUtilities::GetNormMethod<TDataType>(rVariable, rNormType);
126 
127 #pragma omp parallel
128  {
129  double sum = 0.0;
130 #pragma omp for
131  for (int i = 0; i < static_cast<int>(r_container.size()); ++i)
132  {
133  const TContainerItemType& r_item = *(r_container.begin() + i);
134  const TDataType& current_value =
135  TDataRetrievalFunctor<TContainerItemType>()(r_item, rVariable);
136  sum += norm_method(current_value);
137  }
138 #pragma omp atomic
139  global_sum += sum;
140  }
141 
142  global_sum = rModelPart.GetCommunicator().GetDataCommunicator().SumAll(global_sum);
143  return global_sum;
144 
145  KRATOS_CATCH("");
146  }
147 
148  template <class TDataType>
149  TDataType static CalculateRootMeanSquare(const ModelPart& rModelPart, const Variable<TDataType>& rVariable)
150  {
151  KRATOS_TRY
152 
153  const TContainerType& r_container =
154  MethodUtilities::GetLocalDataContainer<TContainerType>(rModelPart);
155 
156  TDataType global_sum = rVariable.Zero();
157 
158  if (r_container.size() > 0)
159  {
160  const TDataType& r_initial_value =
161  TDataRetrievalFunctor<TContainerItemType>()(*r_container.begin(), rVariable);
162 
163  MethodUtilities::DataTypeSizeInitializer(global_sum, r_initial_value);
164 
165 #pragma omp parallel
166  {
167  TDataType sum = rVariable.Zero();
169 
170 #pragma omp for
171  for (int i = 0; i < static_cast<int>(r_container.size()); ++i)
172  {
173  const TContainerItemType& r_item = *(r_container.begin() + i);
174  const TDataType& current_value =
175  TDataRetrievalFunctor<TContainerItemType>()(r_item, rVariable);
177  sum += MethodUtilities::RaiseToPower<TDataType>(current_value, 2);
178  }
179 #pragma omp critical
180  {
181  global_sum += sum;
182  }
183  }
184  }
185 
186  global_sum = rModelPart.GetCommunicator().GetDataCommunicator().SumAll(global_sum);
187  const double number_of_items =
188  rModelPart.GetCommunicator().GetDataCommunicator().SumAll(
189  static_cast<double>(r_container.size()));
190  global_sum = MethodUtilities::RaiseToPower<TDataType>(
191  global_sum * (1.0 / std::max(number_of_items, 1.0)), 0.5);
192 
193  return global_sum;
194 
195  KRATOS_CATCH("");
196  }
197 
198  template <class TDataType>
200  const ModelPart& rModelPart,
201  const Variable<TDataType>& rVariable,
202  const std::string& rNormType,
203  Parameters Params)
204  {
205  KRATOS_TRY
206 
207  const TContainerType& r_container =
208  MethodUtilities::GetLocalDataContainer<TContainerType>(rModelPart);
209 
210  double global_sum = 0.0;
211  const auto& norm_method =
212  MethodUtilities::GetNormMethod<TDataType>(rVariable, rNormType);
213 
214 #pragma omp parallel
215  {
216  double sum = 0.0;
217 #pragma omp for
218  for (int i = 0; i < static_cast<int>(r_container.size()); ++i)
219  {
220  const TContainerItemType& r_item = *(r_container.begin() + i);
221  const TDataType& current_value =
222  TDataRetrievalFunctor<TContainerItemType>()(r_item, rVariable);
223  sum += std::pow(norm_method(current_value), 2);
224  }
225 #pragma omp atomic
226  global_sum += sum;
227  }
228 
229  global_sum = rModelPart.GetCommunicator().GetDataCommunicator().SumAll(global_sum);
230  const double number_of_items =
231  rModelPart.GetCommunicator().GetDataCommunicator().SumAll(
232  static_cast<double>(r_container.size()));
233  return std::sqrt(global_sum / std::max(number_of_items, 1.0));
234 
235  KRATOS_CATCH("");
236  }
237 
238  template <class TDataType>
239  TDataType static CalculateMean(const ModelPart& rModelPart, const Variable<TDataType>& rVariable)
240  {
241  const TDataType& sum = CalculateSum<TDataType>(rModelPart, rVariable);
242  const TContainerType& r_container =
243  MethodUtilities::GetLocalDataContainer<TContainerType>(rModelPart);
244 
245  const double number_of_items =
246  rModelPart.GetCommunicator().GetDataCommunicator().SumAll(
247  static_cast<double>(r_container.size()));
248 
249  return sum * (1.0 / std::max(number_of_items, 1.0));
250  }
251 
252  template <class TDataType>
253  double static CalculateNormMean(
254  const ModelPart& rModelPart,
255  const Variable<TDataType>& rVariable,
256  const std::string& rNormType,
257  Parameters Params)
258  {
259  const double sum =
260  CalculateNormSum<TDataType>(rModelPart, rVariable, rNormType, Params);
261  const TContainerType& r_container =
262  MethodUtilities::GetLocalDataContainer<TContainerType>(rModelPart);
263 
264  const double number_of_items =
265  rModelPart.GetCommunicator().GetDataCommunicator().SumAll(
266  static_cast<double>(r_container.size()));
267 
268  if (number_of_items > 0)
269  {
270  return sum * (1.0 / number_of_items);
271  }
272 
273  return 0.0;
274  }
275 
276  template <class TDataType>
277  std::tuple<TDataType, TDataType> static CalculateVariance(
278  const ModelPart& rModelPart, const Variable<TDataType>& rVariable)
279  {
280  TDataType mean = CalculateMean<TDataType>(rModelPart, rVariable);
281  TDataType global_variance = rVariable.Zero();
282 
283  const TContainerType& r_container =
284  MethodUtilities::GetLocalDataContainer<TContainerType>(rModelPart);
285 
286  if (r_container.size() > 0)
287  {
288  const TDataType& r_initial_value =
289  TDataRetrievalFunctor<TContainerItemType>()(*r_container.begin(), rVariable);
290 
291  MethodUtilities::DataTypeSizeInitializer(global_variance, r_initial_value);
292 
293 #pragma omp parallel
294  {
295  TDataType variance = rVariable.Zero();
296  MethodUtilities::DataTypeSizeInitializer(variance, r_initial_value);
297 
298 #pragma omp for
299  for (int i = 0; i < static_cast<int>(r_container.size()); ++i)
300  {
301  const TContainerItemType& r_item = *(r_container.begin() + i);
302  const TDataType& current_value =
303  TDataRetrievalFunctor<TContainerItemType>()(r_item, rVariable);
304  MethodUtilities::DataTypeSizeChecker(current_value, variance);
305 
306  variance += MethodUtilities::RaiseToPower(current_value, 2);
307  }
308 #pragma omp critical
309  {
310  global_variance += variance;
311  }
312  }
313  }
314 
315  global_variance =
316  rModelPart.GetCommunicator().GetDataCommunicator().SumAll(global_variance);
317  const double number_of_items =
318  rModelPart.GetCommunicator().GetDataCommunicator().SumAll(
319  static_cast<double>(r_container.size()));
320 
321  if (number_of_items > 0)
322  {
323  global_variance *= (1.0 / number_of_items);
324  global_variance -= MethodUtilities::RaiseToPower(mean, 2);
325  }
326 
327  return std::make_tuple<TDataType, TDataType>(
328  std::forward<TDataType>(mean), std::forward<TDataType>(global_variance));
329  }
330 
331  template <class TDataType>
332  std::tuple<double, double> static CalculateNormVariance(
333  const ModelPart& rModelPart,
334  const Variable<TDataType>& rVariable,
335  const std::string& rNormType,
336  Parameters Params)
337  {
338  double mean = CalculateNormMean<TDataType>(rModelPart, rVariable, rNormType, Params);
339 
340  const TContainerType& r_container =
341  MethodUtilities::GetLocalDataContainer<TContainerType>(rModelPart);
342 
343  double global_variance = 0.0;
344  const auto& norm_method =
345  MethodUtilities::GetNormMethod<TDataType>(rVariable, rNormType);
346 
347 #pragma omp parallel
348  {
349  double variance = 0.0;
350 
351 #pragma omp for
352  for (int i = 0; i < static_cast<int>(r_container.size()); ++i)
353  {
354  const TContainerItemType& r_item = *(r_container.begin() + i);
355  const TDataType& current_value =
356  TDataRetrievalFunctor<TContainerItemType>()(r_item, rVariable);
357 
358  variance += std::pow(norm_method(current_value), 2);
359  }
360 #pragma omp atomic
361  global_variance += variance;
362  }
363  global_variance =
364  rModelPart.GetCommunicator().GetDataCommunicator().SumAll(global_variance);
365  const double number_of_items =
366  rModelPart.GetCommunicator().GetDataCommunicator().SumAll(
367  static_cast<double>(r_container.size()));
368 
369  if (number_of_items > 0)
370  {
371  global_variance *= (1.0 / number_of_items);
372  global_variance -= MethodUtilities::RaiseToPower(mean, 2);
373  }
374 
375  return std::make_tuple<double, double>(
376  std::forward<double>(mean), std::forward<double>(global_variance));
377  }
378 
379  template <class TDataType>
380  std::tuple<double, std::size_t> static GetNormMax(
381  const ModelPart& rModelPart,
382  const Variable<TDataType>& rVariable,
383  const std::string& rNormType,
384  Parameters Params)
385  {
386  KRATOS_TRY
387 
388  const TContainerType& r_container =
389  MethodUtilities::GetLocalDataContainer<TContainerType>(rModelPart);
390 
391  double global_max = std::numeric_limits<double>::lowest();
392  unsigned int global_id = 0;
393  const auto& norm_method =
394  MethodUtilities::GetNormMethod<TDataType>(rVariable, rNormType);
395 
396 #pragma omp parallel
397  {
398  double current_max = std::numeric_limits<double>::lowest();
399  unsigned int current_id = 0;
400 #pragma omp for
401  for (int i = 0; i < static_cast<int>(r_container.size()); ++i)
402  {
403  const TContainerItemType& r_item = *(r_container.begin() + i);
404  const TDataType& current_value =
405  TDataRetrievalFunctor<TContainerItemType>()(r_item, rVariable);
406  const double value_norm = norm_method(current_value);
407  if (value_norm > current_max)
408  {
409  current_max = value_norm;
410  current_id = r_item.Id();
411  }
412  }
413 #pragma omp critical
414  {
415  if (current_max > global_max)
416  {
417  global_max = current_max;
418  global_id = current_id;
419  }
420  }
421  }
422 
423  const DataCommunicator& r_data_communicator =
424  rModelPart.GetCommunicator().GetDataCommunicator();
425  const auto& global_max_value_array =
426  r_data_communicator.AllGather(std::vector<double>{global_max});
427  const auto& global_max_id_array =
428  r_data_communicator.AllGather(std::vector<unsigned int>{global_id});
429 
430  for (std::size_t i = 0; i < global_max_value_array.size(); ++i)
431  {
432  if (global_max_value_array[i] > global_max)
433  {
434  global_max = global_max_value_array[i];
435  global_id = global_max_id_array[i];
436  }
437  }
438 
439  return std::make_tuple<double, unsigned int>(
440  std::forward<double>(global_max), std::forward<unsigned int>(global_id));
441 
442  KRATOS_CATCH("");
443  }
444 
445  template <class TDataType>
446  std::tuple<double, std::size_t> static GetNormMin(
447  const ModelPart& rModelPart,
448  const Variable<TDataType>& rVariable,
449  const std::string& rNormType,
450  Parameters Params)
451  {
452  KRATOS_TRY
453 
454  const TContainerType& r_container =
455  MethodUtilities::GetLocalDataContainer<TContainerType>(rModelPart);
456 
457  double global_min = std::numeric_limits<double>::max();
458  unsigned int global_id = 0;
459  const auto& norm_method =
460  MethodUtilities::GetNormMethod<TDataType>(rVariable, rNormType);
461 
462 #pragma omp parallel
463  {
464  double current_min = std::numeric_limits<double>::max();
465  unsigned int current_id = 0;
466 #pragma omp for
467  for (int i = 0; i < static_cast<int>(r_container.size()); ++i)
468  {
469  const TContainerItemType& r_item = *(r_container.begin() + i);
470  const TDataType& current_value =
471  TDataRetrievalFunctor<TContainerItemType>()(r_item, rVariable);
472  const double value_norm = norm_method(current_value);
473  if (value_norm < current_min)
474  {
475  current_min = value_norm;
476  current_id = r_item.Id();
477  }
478  }
479 #pragma omp critical
480  {
481  if (current_min < global_min)
482  {
483  global_min = current_min;
484  global_id = current_id;
485  }
486  }
487  }
488 
489  const DataCommunicator& r_data_communicator =
490  rModelPart.GetCommunicator().GetDataCommunicator();
491  const auto& global_min_value_array =
492  r_data_communicator.AllGather(std::vector<double>{global_min});
493  const auto& global_min_id_array =
494  r_data_communicator.AllGather(std::vector<unsigned int>{global_id});
495 
496  for (std::size_t i = 0; i < global_min_value_array.size(); ++i)
497  {
498  if (global_min_value_array[i] < global_min)
499  {
500  global_min = global_min_value_array[i];
501  global_id = global_min_id_array[i];
502  }
503  }
504 
505  return std::make_tuple<double, unsigned int>(
506  std::forward<double>(global_min), std::forward<unsigned int>(global_id));
507 
508  KRATOS_CATCH("");
509  }
510 
511  template <class TDataType>
512  double static GetNormMedian(
513  const ModelPart& rModelPart,
514  const Variable<TDataType>& rVariable,
515  const std::string& rNormType,
516  Parameters Params)
517  {
518  KRATOS_TRY
519 
520  const TContainerType& r_container =
521  MethodUtilities::GetLocalDataContainer<TContainerType>(rModelPart);
522 
523  const auto& norm_method =
524  MethodUtilities::GetNormMethod<TDataType>(rVariable, rNormType);
525 
526  std::vector<double> local_values;
527  local_values.resize(r_container.size());
528  local_values.shrink_to_fit();
529 
530 #pragma omp parallel for
531  for (int i = 0; i < static_cast<int>(r_container.size()); ++i)
532  {
533  const TContainerItemType& r_item = *(r_container.begin() + i);
534  const TDataType& current_value =
535  TDataRetrievalFunctor<TContainerItemType>()(r_item, rVariable);
536  local_values[i] = norm_method(current_value);
537  }
538 
539  std::sort(local_values.begin(), local_values.end());
540 
541  const DataCommunicator& r_data_communicator =
542  rModelPart.GetCommunicator().GetDataCommunicator();
543  const std::vector<std::vector<double>>& global_values =
544  r_data_communicator.Gatherv(local_values, 0);
545 
546  double median = 0.0;
547  if (r_data_communicator.Rank() == 0)
548  {
549  const std::vector<double>& sorted_values_list =
551  const int number_of_values = sorted_values_list.size();
552 
553  if (number_of_values > 0)
554  {
555  if (number_of_values % 2 != 0)
556  {
557  median = sorted_values_list[number_of_values / 2];
558  }
559  else
560  {
561  median = (sorted_values_list[(number_of_values - 1) / 2] +
562  sorted_values_list[number_of_values / 2]) *
563  0.5;
564  }
565  }
566  }
567 
568  r_data_communicator.Broadcast(median, 0);
569  return median;
570 
571  KRATOS_CATCH("");
572  }
573 
574  template <class TDataType>
575  std::tuple<double, double, std::vector<double>, std::vector<int>, std::vector<double>, std::vector<double>, std::vector<double>> static GetNormDistribution(
576 
577  const ModelPart& rModelPart,
578  const Variable<TDataType>& rVariable,
579  const std::string& rNormType,
580  Parameters Params)
581  {
582  KRATOS_TRY
583 
584  Parameters default_parameters = Parameters(R"(
585  {
586  "number_of_value_groups" : 10,
587  "min_value" : "min",
588  "max_value" : "max"
589  })");
590 
591  if (Params.Has("min_value") && Params["min_value"].IsDouble())
592  {
593  default_parameters["min_value"].SetDouble(0.0);
594  }
595  if (Params.Has("max_value") && Params["max_value"].IsDouble())
596  {
597  default_parameters["max_value"].SetDouble(0.0);
598  }
599  Params.RecursivelyValidateAndAssignDefaults(default_parameters);
600 
601  double min_value{0.0};
602  if (Params["min_value"].IsDouble())
603  {
604  min_value = Params["min_value"].GetDouble();
605  }
606  else if (
607  Params["min_value"].IsString() &&
608  Params["min_value"].GetString() == "min")
609  {
610  const auto& min_data =
611  GetNormMin<TDataType>(rModelPart, rVariable, rNormType, Params);
612  min_value = std::get<0>(min_data);
613  }
614  else
615  {
616  KRATOS_ERROR << "Unknown min_value. Allowed only double or \"min\" "
617  "string as a value. [ min_value = "
618  << Params["min_value"] << " ]\n.";
619  }
620 
621  double max_value{0.0};
622  if (Params["max_value"].IsDouble())
623  {
624  max_value = Params["max_value"].GetDouble();
625  }
626  else if (
627  Params["max_value"].IsString() &&
628  Params["max_value"].GetString() == "max")
629  {
630  const auto& max_data =
631  GetNormMax<TDataType>(rModelPart, rVariable, rNormType, Params);
632  max_value = std::get<0>(max_data);
633  }
634  else
635  {
636  KRATOS_ERROR << "Unknown max_value. Allowed only double or \"max\" "
637  "string as a value. [ max_value = "
638  << Params["max_value"] << " ]\n.";
639  }
640 
641  const int number_of_groups = Params["number_of_value_groups"].GetInt();
642 
643  const TContainerType& r_container =
644  MethodUtilities::GetLocalDataContainer<TContainerType>(rModelPart);
645 
646  const auto& norm_method =
647  MethodUtilities::GetNormMethod<TDataType>(rVariable, rNormType);
648 
649  std::vector<double> group_limits;
650  for (int i = 0; i < number_of_groups + 1; ++i)
651  {
652  group_limits.push_back(
653  min_value + (max_value - min_value) * static_cast<double>(i) /
654  static_cast<double>(number_of_groups));
655  }
656 
657  // final group limit is extended by a small amount. epsilon in numeric
658  // limits cannot be used since testing also need to have the same
659  // extending value in python. Therefore hard coded value is used
660  group_limits[group_limits.size() - 1] += 1e-16;
661  group_limits.push_back(std::numeric_limits<double>::max());
662 
663  group_limits.shrink_to_fit();
664  const int number_of_limits = group_limits.size();
665 
666  std::vector<int> distribution;
667  std::vector<double> group_means, group_variances;
668  for (int i = 0; i < number_of_limits; ++i)
669  {
670  distribution.push_back(0);
671  group_means.push_back(0.0);
672  group_variances.push_back(0.0);
673  }
674  distribution.shrink_to_fit();
675  group_means.shrink_to_fit();
676  group_variances.shrink_to_fit();
677 
678 #pragma omp parallel
679  {
680  std::vector<int> local_distribution;
681  std::vector<double> local_means, local_variances;
682  for (int i = 0; i < number_of_limits; ++i)
683  {
684  local_distribution.push_back(0);
685  local_means.push_back(0.0);
686  local_variances.push_back(0.0);
687  }
688  local_distribution.shrink_to_fit();
689  local_means.shrink_to_fit();
690  local_variances.shrink_to_fit();
691 
692 #pragma omp for
693  for (int i = 0; i < static_cast<int>(r_container.size()); ++i)
694  {
695  const TContainerItemType& r_item = *(r_container.begin() + i);
696  const TDataType& current_value =
697  TDataRetrievalFunctor<TContainerItemType>()(r_item, rVariable);
698  const double value_norm = norm_method(current_value);
699  for (int i = 0; i < number_of_limits; ++i)
700  {
701  if (value_norm < group_limits[i])
702  {
703  ++local_distribution[i];
704  local_means[i] += value_norm;
705  local_variances[i] += std::pow(value_norm, 2);
706  break;
707  }
708  }
709  }
710 #pragma omp critical
711  {
712  for (int i = 0; i < number_of_limits; ++i)
713  {
714  distribution[i] += local_distribution[i];
715  group_means[i] += local_means[i];
716  group_variances[i] += local_variances[i];
717  }
718  }
719  }
720 
721  std::vector<int> global_distribution =
722  rModelPart.GetCommunicator().GetDataCommunicator().SumAll(distribution);
723  std::vector<double> global_mean_distribution =
724  rModelPart.GetCommunicator().GetDataCommunicator().SumAll(group_means);
725  std::vector<double> global_variance_distribution =
726  rModelPart.GetCommunicator().GetDataCommunicator().SumAll(group_variances);
727 
728  const double number_of_items = static_cast<double>(std::max(
729  std::accumulate(global_distribution.begin(), global_distribution.end(), 0), 1));
730  std::vector<double> global_percentage_distributions;
731  for (int i = 0; i < number_of_limits; ++i)
732  {
733  const double number_of_values_in_group =
734  static_cast<double>(global_distribution[i]);
735  global_percentage_distributions.push_back(number_of_values_in_group / number_of_items);
736  if (number_of_values_in_group > 0.0)
737  {
738  global_mean_distribution[i] /= number_of_values_in_group;
739  global_variance_distribution[i] /= number_of_values_in_group;
740  global_variance_distribution[i] -=
741  std::pow(global_mean_distribution[i], 2);
742  }
743  }
744 
745  // reversing group limit is extention
746  group_limits[group_limits.size() - 2] -= 1e-16;
747  group_limits[group_limits.size() - 1] = max_value;
748 
749  return std::make_tuple<
750  double, double, std::vector<double>, std::vector<int>,
751  std::vector<double>, std::vector<double>, std::vector<double>>(
752  std::forward<double>(min_value), std::forward<double>(max_value),
753  std::forward<std::vector<double>>(group_limits),
754  std::forward<std::vector<int>>(global_distribution),
755  std::forward<std::vector<double>>(global_percentage_distributions),
756  std::forward<std::vector<double>>(global_mean_distribution),
757  std::forward<std::vector<double>>(global_variance_distribution));
758 
759  KRATOS_CATCH("");
760  }
761 };
762 
766 
770 
772  : public SpatialMethods::ContainerSpatialMethods<NodesContainerType, NodeType, MethodUtilities::HistoricalDataValueRetrievalFunctor>
773 {
774 };
775 
777  : public SpatialMethods::ContainerSpatialMethods<NodesContainerType, NodeType, MethodUtilities::NonHistoricalDataValueRetrievalFunctor>
778 {
779 };
780 
782  : public SpatialMethods::ContainerSpatialMethods<ConditionsContainerType, ConditionType, MethodUtilities::NonHistoricalDataValueRetrievalFunctor>
783 {
784 };
785 
787  : public SpatialMethods::ContainerSpatialMethods<ElementsContainerType, ElementType, MethodUtilities::NonHistoricalDataValueRetrievalFunctor>
788 {
789 };
790 
791 } // namespace SpatialMethods
792 
796 
800 
804 
808 
810 
813 
817 
819 
821 
822 } // namespace Kratos.
823 
824 #endif // KRATOS_SPATIAL_VARIANCE_METHOD_H_INCLUDED defined
virtual const DataCommunicator & GetDataCommunicator() const
Definition: communicator.cpp:340
Base class for all Conditions.
Definition: condition.h:59
Serial (do-nothing) version of a wrapper class for MPI communication.
Definition: data_communicator.h:318
virtual int Rank() const
Get the parallel rank for this DataCommunicator.
Definition: data_communicator.h:587
void Broadcast(TObject &rBroadcastObject, const int SourceRank) const
Synchronize a buffer to the value held by the broadcasting rank.
Definition: data_communicator.h:473
Base class for all Elements.
Definition: element.h:60
Definition: flags.h:58
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType::ConditionsContainerType ConditionsContainerType
Condintions container. A vector set of Conditions with their Id's as key.
Definition: model_part.h:183
MeshType::ElementsContainerType ElementsContainerType
Element container. A vector set of Elements with their Id's as key.
Definition: model_part.h:168
Element ElementType
Definition: model_part.h:120
Communicator & GetCommunicator()
Definition: model_part.h:1821
Node NodeType
Definition: model_part.h:117
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
Condition ConditionType
Definition: model_part.h:121
This class defines the node.
Definition: node.h:65
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
int GetInt() const
This method returns the integer contained in the current Parameter.
Definition: kratos_parameters.cpp:666
bool IsDouble() const
This method checks if the parameter is a double.
Definition: kratos_parameters.cpp:544
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
bool Has(const std::string &rEntry) const
This method checks if the Parameter contains a certain entry.
Definition: kratos_parameters.cpp:520
Definition: spatial_methods.h:46
static double CalculateNormRootMeanSquare(const ModelPart &rModelPart, const Variable< TDataType > &rVariable, const std::string &rNormType, Parameters Params)
Definition: spatial_methods.h:199
static TDataType CalculateRootMeanSquare(const ModelPart &rModelPart, const Variable< TDataType > &rVariable)
Definition: spatial_methods.h:149
static int CalculateSum(const ModelPart &rModelPart, const Flags &rVariable)
Definition: spatial_methods.h:49
static double CalculateNormMean(const ModelPart &rModelPart, const Variable< TDataType > &rVariable, const std::string &rNormType, Parameters Params)
Definition: spatial_methods.h:253
static std::tuple< double, std::size_t > GetNormMax(const ModelPart &rModelPart, const Variable< TDataType > &rVariable, const std::string &rNormType, Parameters Params)
Definition: spatial_methods.h:380
static double CalculateNormSum(const ModelPart &rModelPart, const Variable< TDataType > &rVariable, const std::string &rNormType, Parameters Params)
Definition: spatial_methods.h:112
static std::tuple< double, double, std::vector< double >, std::vector< int >, std::vector< double >, std::vector< double >, std::vector< double > > GetNormDistribution(const ModelPart &rModelPart, const Variable< TDataType > &rVariable, const std::string &rNormType, Parameters Params)
Definition: spatial_methods.h:575
static std::tuple< double, double > CalculateNormVariance(const ModelPart &rModelPart, const Variable< TDataType > &rVariable, const std::string &rNormType, Parameters Params)
Definition: spatial_methods.h:332
static TDataType CalculateSum(const ModelPart &rModelPart, const Variable< TDataType > &rVariable)
Definition: spatial_methods.h:67
static std::tuple< double, std::size_t > GetNormMin(const ModelPart &rModelPart, const Variable< TDataType > &rVariable, const std::string &rNormType, Parameters Params)
Definition: spatial_methods.h:446
static std::tuple< TDataType, TDataType > CalculateVariance(const ModelPart &rModelPart, const Variable< TDataType > &rVariable)
Definition: spatial_methods.h:277
static double GetNormMedian(const ModelPart &rModelPart, const Variable< TDataType > &rVariable, const std::string &rNormType, Parameters Params)
Definition: spatial_methods.h:512
static TDataType CalculateMean(const ModelPart &rModelPart, const Variable< TDataType > &rVariable)
Definition: spatial_methods.h:239
Definition: spatial_methods.h:773
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
const TDataType & Zero() const
This method returns the zero value of the variable type.
Definition: variable.h:346
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR
Definition: exception.h:161
static double max(double a, double b)
Definition: GeometryFunctions.h:79
std::vector< double > SortSortedValuesList(const std::vector< std::vector< double >> &rValues)
Definition: method_utilities.cpp:725
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 DataTypeSizeInitializer(TDataType &rData, const TDataType &rReferenceData)
Definition: method_utilities.cpp:82
ModelPart::ElementsContainerType ElementsContainerType
Definition: spatial_methods.h:768
ModelPart::ConditionsContainerType ConditionsContainerType
Definition: spatial_methods.h:769
ModelPart::NodesContainerType NodesContainerType
Definition: spatial_methods.h:767
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
TExpressionType::data_type sum(AMatrix::MatrixExpression< TExpressionType, TCategory > const &TheExpression)
Definition: amatrix_interface.h:675
TABLE_NUMBER_ANGULAR_VELOCITY TABLE_NUMBER_MOMENT I33 BEAM_INERTIA_ROT_UNIT_LENGHT_Y KRATOS_DEFINE_APPLICATION_VARIABLE(DEM_APPLICATION, double, BEAM_INERTIA_ROT_UNIT_LENGHT_Z) typedef std double
Definition: DEM_application_variables.h:182
int current_id
Output settings end ####.
Definition: script.py:194
integer i
Definition: TensorModule.f:17
e
Definition: run_cpp_mpi_tests.py:31