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.
variable_utils.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: Riccardo Rossi
11 // Ruben Zorrilla
12 // Vicente Mataix Ferrandiz
13 //
14 
15 #pragma once
16 
17 // System includes
18 #include <type_traits>
19 
20 // External includes
21 
22 // Project includes
23 #include "includes/define.h"
24 #include "includes/model_part.h"
25 #include "includes/checks.h"
30 
31 namespace Kratos
32 {
35 
39 
43 
47 
51 
62 class KRATOS_API(KRATOS_CORE) VariableUtils
63 {
64 public:
67 
70 
73 
76 
79 
82 
85 
88 
91 
94 
98 
108 
112 
125  template <class TVarType>
127  const TVarType &rVariable,
128  const TVarType &rDestinationVariable,
129  const ModelPart &rOriginModelPart,
130  ModelPart &rDestinationModelPart,
131  const unsigned int ReadBufferStep,
132  const unsigned int WriteBufferStep )
133  {
134  const int n_orig_nodes = rOriginModelPart.NumberOfNodes();
135  const int n_dest_nodes = rDestinationModelPart.NumberOfNodes();
136 
137  KRATOS_ERROR_IF_NOT(n_orig_nodes == n_dest_nodes)
138  << "Origin and destination model parts have different number of nodes."
139  << "\n\t- Number of origin nodes: " << n_orig_nodes
140  << "\n\t- Number of destination nodes: " << n_dest_nodes << std::endl;
141 
142  IndexPartition<std::size_t>(n_orig_nodes).for_each([&](std::size_t index)
143  {
144  auto it_dest_node = rDestinationModelPart.NodesBegin() + index;
145  const auto it_orig_node = rOriginModelPart.NodesBegin() + index;
146  const auto &r_value = it_orig_node->GetSolutionStepValue(rVariable, ReadBufferStep);
147  it_dest_node->FastGetSolutionStepValue(rDestinationVariable, WriteBufferStep) = r_value; });
148 
149  rDestinationModelPart.GetCommunicator().SynchronizeVariable(rDestinationVariable);
150  }
151 
162  template <class TVarType>
164  const TVarType &rVariable,
165  const TVarType &rDestinationVariable,
166  const ModelPart &rOriginModelPart,
167  ModelPart &rDestinationModelPart,
168  const unsigned int BuffStep = 0)
169  {
170  this->CopyModelPartNodalVar(rVariable, rDestinationVariable, rOriginModelPart, rDestinationModelPart, BuffStep, BuffStep);
171  }
172 
182  template< class TVarType >
184  const TVarType& rVariable,
185  const ModelPart& rOriginModelPart,
186  ModelPart& rDestinationModelPart,
187  const unsigned int BuffStep = 0)
188  {
189  this->CopyModelPartNodalVar(rVariable, rVariable, rOriginModelPart, rDestinationModelPart, BuffStep);
190  }
191 
192  template< class TVarType >
194  const TVarType &rVariable,
195  const TVarType &rDestinationVariable,
196  const ModelPart &rOriginModelPart,
197  ModelPart &rDestinationModelPart,
198  const unsigned int BuffStep = 0)
199  {
200  const int n_orig_nodes = rOriginModelPart.NumberOfNodes();
201  const int n_dest_nodes = rDestinationModelPart.NumberOfNodes();
202 
203  KRATOS_ERROR_IF_NOT(n_orig_nodes == n_dest_nodes) <<
204  "Origin and destination model parts have different number of nodes." <<
205  "\n\t- Number of origin nodes: " << n_orig_nodes <<
206  "\n\t- Number of destination nodes: " << n_dest_nodes << std::endl;
207 
208  IndexPartition<std::size_t>(n_orig_nodes).for_each([&](std::size_t index){
209  auto it_dest_node = rDestinationModelPart.NodesBegin() + index;
210  const auto it_orig_node = rOriginModelPart.NodesBegin() + index;
211  const auto& r_value = it_orig_node->GetSolutionStepValue(rVariable, BuffStep);
212  it_dest_node->GetValue(rDestinationVariable) = r_value;
213  });
214 
215  rDestinationModelPart.GetCommunicator().SynchronizeNonHistoricalVariable(rDestinationVariable);
216  }
217 
218  template< class TVarType >
220  const TVarType &rVariable,
221  const ModelPart &rOriginModelPart,
222  ModelPart &rDestinationModelPart,
223  const unsigned int BuffStep = 0)
224  {
225  this->CopyModelPartNodalVarToNonHistoricalVar(rVariable, rVariable, rOriginModelPart, rDestinationModelPart, BuffStep);
226  }
227 
228  template <class TDataType>
230  const Variable<TDataType>& rOriginVariable,
231  const Variable<TDataType>& rDestinationVariable,
232  const ModelPart& rOriginModelPart,
233  ModelPart& rDestinationModelPart,
234  const Flags& rFlag,
235  const bool CheckValue = true,
236  const unsigned int ReadBufferStep = 0,
237  const unsigned int WriteBufferStep = 0)
238  {
239  KRATOS_TRY
240 
242  rOriginModelPart.FullName() == rDestinationModelPart.FullName() &&
243  rOriginVariable == rDestinationVariable &&
244  ReadBufferStep == WriteBufferStep)
245  << "Trying to copy flagged nodal solution step values with the same origin and destination model parts/variables/buffer steps. This is not permitted ( Origin model part: "
246  << rOriginModelPart.Name() << ", destination model part: " << rDestinationModelPart.Name()
247  << ", variable: " << rOriginVariable.Name() << ", buffer step: " << ReadBufferStep << " ) !";
248 
249  KRATOS_ERROR_IF_NOT(rOriginModelPart.HasNodalSolutionStepVariable(rOriginVariable))
250  << rOriginVariable.Name() << " is not found in nodal solution step variables list in origin model part ( "
251  << rOriginModelPart.Name() << " ).";
252 
253  KRATOS_ERROR_IF_NOT(rDestinationModelPart.HasNodalSolutionStepVariable(rDestinationVariable))
254  << rDestinationVariable.Name() << " is not found in nodal solution step variables list in destination model part ( "
255  << rDestinationModelPart.Name() << " ).";
256 
257  KRATOS_ERROR_IF(ReadBufferStep >= rOriginModelPart.GetBufferSize())
258  << "Origin model part ( " << rOriginModelPart.Name()
259  << " ) buffer size is smaller or equal than read buffer size [ "
260  << rOriginModelPart.GetBufferSize() << " <= " << ReadBufferStep << " ].";
261 
262  KRATOS_ERROR_IF(WriteBufferStep >= rDestinationModelPart.GetBufferSize())
263  << "Destination model part ( " << rDestinationModelPart.Name()
264  << " ) buffer size is smaller or equal than read buffer size [ "
265  << rDestinationModelPart.GetBufferSize() << " <= " << WriteBufferStep << " ].";
266 
267  CopyModelPartFlaggedVariable<NodesContainerType>(
268  rOriginModelPart, rDestinationModelPart, rFlag, CheckValue,
269  [&](NodeType& rDestNode, const TDataType& rValue) {
270  rDestNode.FastGetSolutionStepValue(
271  rDestinationVariable, WriteBufferStep) = rValue;
272  },
273  [&](const NodeType& rOriginNode) -> const TDataType& {
274  return rOriginNode.FastGetSolutionStepValue(rOriginVariable, ReadBufferStep);
275  });
276 
277  rDestinationModelPart.GetCommunicator().SynchronizeVariable(rDestinationVariable);
278 
279  KRATOS_CATCH("");
280  }
281 
282  template <class TDataType>
284  const Variable<TDataType>& rOriginVariable,
285  const Variable<TDataType>& rDestinationVariable,
286  ModelPart& rModelPart,
287  const Flags& rFlag,
288  const bool CheckValue = true,
289  const unsigned int ReadBufferStep = 0,
290  const unsigned int WriteBufferStep = 0)
291  {
292  KRATOS_TRY
293 
294  CopyModelPartFlaggedNodalHistoricalVarToHistoricalVar(
295  rOriginVariable, rDestinationVariable, rModelPart, rModelPart,
296  rFlag, CheckValue, ReadBufferStep, WriteBufferStep);
297 
298  KRATOS_CATCH("");
299  }
300 
301  template <class TDataType>
303  const Variable<TDataType>& rVariable,
304  const ModelPart& rOriginModelPart,
305  ModelPart& rDestinationModelPart,
306  const Flags& rFlag,
307  const bool CheckValue = true,
308  const unsigned int ReadBufferStep = 0,
309  const unsigned int WriteBufferStep = 0)
310  {
311  KRATOS_TRY
312 
313  CopyModelPartFlaggedNodalHistoricalVarToHistoricalVar(
314  rVariable, rVariable, rOriginModelPart, rDestinationModelPart,
315  rFlag, CheckValue, ReadBufferStep, WriteBufferStep);
316 
317  KRATOS_CATCH("");
318  }
319 
320  template <class TDataType>
322  const Variable<TDataType>& rOriginVariable,
323  const Variable<TDataType>& rDestinationVariable,
324  const ModelPart& rOriginModelPart,
325  ModelPart& rDestinationModelPart,
326  const Flags& rFlag,
327  const bool CheckValue = true,
328  const unsigned int ReadBufferStep = 0)
329  {
330  KRATOS_TRY
331 
332  KRATOS_ERROR_IF_NOT(rOriginModelPart.HasNodalSolutionStepVariable(rOriginVariable))
333  << rOriginVariable.Name() << " is not found in nodal solution step variables list in origin model part ( "
334  << rOriginModelPart.Name() << " ).";
335 
336  KRATOS_ERROR_IF(ReadBufferStep >= rOriginModelPart.GetBufferSize())
337  << "Origin model part ( " << rOriginModelPart.Name()
338  << " ) buffer size is smaller or equal than read buffer size [ "
339  << rOriginModelPart.GetBufferSize() << " <= " << ReadBufferStep << " ].";
340 
341 
342  CopyModelPartFlaggedVariable<NodesContainerType>(
343  rOriginModelPart, rDestinationModelPart, rFlag, CheckValue,
344  [&](NodeType& rDestNode, const TDataType& rValue) {
345  rDestNode.SetValue(rDestinationVariable, rValue);
346  },
347  [&](const NodeType& rOriginNode) -> const TDataType& {
348  return rOriginNode.FastGetSolutionStepValue(rOriginVariable, ReadBufferStep);
349  });
350 
351  rDestinationModelPart.GetCommunicator().SynchronizeNonHistoricalVariable(rDestinationVariable);
352 
353  KRATOS_CATCH("");
354  }
355 
356  template <class TDataType>
358  const Variable<TDataType>& rOriginVariable,
359  const Variable<TDataType>& rDestinationVariable,
360  ModelPart& rModelPart,
361  const Flags& rFlag,
362  const bool CheckValue = true,
363  const unsigned int ReadBufferStep = 0)
364  {
365  CopyModelPartFlaggedNodalHistoricalVarToNonHistoricalVar(
366  rOriginVariable, rDestinationVariable, rModelPart, rModelPart,
367  rFlag, CheckValue, ReadBufferStep);
368  }
369 
370  template <class TDataType>
372  const Variable<TDataType>& rVariable,
373  const ModelPart& rOriginModelPart,
374  ModelPart& rDestinationModelPart,
375  const Flags& rFlag,
376  const bool CheckValue = true,
377  const unsigned int ReadBufferStep = 0)
378  {
379  CopyModelPartFlaggedNodalHistoricalVarToNonHistoricalVar(
380  rVariable, rVariable, rOriginModelPart, rDestinationModelPart,
381  rFlag, CheckValue, ReadBufferStep);
382  }
383 
384  template <class TDataType>
386  const Variable<TDataType>& rVariable,
387  ModelPart& rModelPart,
388  const Flags& rFlag,
389  const bool CheckValue = true,
390  const unsigned int ReadBufferStep = 0)
391  {
392  CopyModelPartFlaggedNodalHistoricalVarToNonHistoricalVar(
393  rVariable, rVariable, rModelPart, rModelPart,
394  rFlag, CheckValue, ReadBufferStep);
395  }
396 
397  template <class TDataType>
399  const Variable<TDataType>& rOriginVariable,
400  const Variable<TDataType>& rDestinationVariable,
401  const ModelPart& rOriginModelPart,
402  ModelPart& rDestinationModelPart,
403  const Flags& rFlag,
404  const bool CheckValue = true,
405  const unsigned int WriteBufferStep = 0)
406  {
407  KRATOS_TRY
408 
409  KRATOS_ERROR_IF_NOT(rDestinationModelPart.HasNodalSolutionStepVariable(rDestinationVariable))
410  << rDestinationVariable.Name() << " is not found in nodal solution step variables list in destination model part ( "
411  << rDestinationModelPart.Name() << " ).";
412 
413  KRATOS_ERROR_IF(WriteBufferStep >= rDestinationModelPart.GetBufferSize())
414  << "Destination model part ( " << rDestinationModelPart.Name()
415  << " ) buffer size is smaller or equal than read buffer size [ "
416  << rDestinationModelPart.GetBufferSize() << " <= " << WriteBufferStep << " ].";
417 
418  CopyModelPartFlaggedVariable<NodesContainerType>(
419  rOriginModelPart, rDestinationModelPart, rFlag, CheckValue,
420  [&](NodeType& rDestNode, const TDataType& rValue) {
421  rDestNode.FastGetSolutionStepValue(
422  rDestinationVariable, WriteBufferStep) = rValue;
423  },
424  [&](const NodeType& rOriginNode) -> const TDataType& {
425  return rOriginNode.GetValue(rOriginVariable);
426  });
427 
428  rDestinationModelPart.GetCommunicator().SynchronizeVariable(rDestinationVariable);
429 
430  KRATOS_CATCH("");
431  }
432 
433  template <class TDataType>
435  const Variable<TDataType>& rOriginVariable,
436  const Variable<TDataType>& rDestinationVariable,
437  ModelPart& rModelPart,
438  const Flags& rFlag,
439  const bool CheckValue = true,
440  const unsigned int WriteBufferStep = 0)
441  {
442  CopyModelPartFlaggedNodalNonHistoricalVarToHistoricalVar(
443  rOriginVariable, rDestinationVariable, rModelPart, rModelPart,
444  rFlag, CheckValue, WriteBufferStep);
445  }
446 
447  template <class TDataType>
449  const Variable<TDataType>& rVariable,
450  const ModelPart& rOriginModelPart,
451  ModelPart& rDestinationModelPart,
452  const Flags& rFlag,
453  const bool CheckValue = true,
454  const unsigned int WriteBufferStep = 0)
455  {
456  CopyModelPartFlaggedNodalNonHistoricalVarToHistoricalVar(
457  rVariable, rVariable, rOriginModelPart, rDestinationModelPart,
458  rFlag, CheckValue, WriteBufferStep);
459  }
460 
461  template <class TDataType>
463  const Variable<TDataType>& rVariable,
464  ModelPart& rModelPart,
465  const Flags& rFlag,
466  const bool CheckValue = true,
467  const unsigned int WriteBufferStep = 0)
468  {
469  CopyModelPartFlaggedNodalNonHistoricalVarToHistoricalVar(
470  rVariable, rVariable, rModelPart, rModelPart,
471  rFlag, CheckValue, WriteBufferStep);
472  }
473 
474  template <class TDataType>
476  const Variable<TDataType>& rOriginVariable,
477  const Variable<TDataType>& rDestinationVariable,
478  const ModelPart& rOriginModelPart,
479  ModelPart& rDestinationModelPart,
480  const Flags& rFlag,
481  const bool CheckValue = true)
482  {
483  KRATOS_TRY
484 
486  rOriginModelPart.FullName() == rDestinationModelPart.FullName() &&
487  rOriginVariable == rDestinationVariable
488  ) << "Trying to copy flagged nodal non-historical values with the same model parts/variables. This is not permitted ( Origin model part: "
489  << rOriginModelPart.Name() << ", destination model part: " << rDestinationModelPart.Name()
490  << ", variable: " << rOriginVariable.Name() << " ) !";
491 
492  CopyModelPartFlaggedVariable<NodesContainerType>(
493  rOriginModelPart, rDestinationModelPart, rFlag, CheckValue,
494  [&](NodeType& rDestNode, const TDataType& rValue) {
495  rDestNode.SetValue(rDestinationVariable, rValue);
496  },
497  [&](const NodeType& rOriginNode) -> const TDataType& {
498  return rOriginNode.GetValue(rOriginVariable);
499  });
500 
501  rDestinationModelPart.GetCommunicator().SynchronizeNonHistoricalVariable(rDestinationVariable);
502 
503  KRATOS_CATCH("");
504  }
505 
506  template <class TDataType>
508  const Variable<TDataType>& rOriginVariable,
509  const Variable<TDataType>& rDestinationVariable,
510  ModelPart& rModelPart,
511  const Flags& rFlag,
512  const bool CheckValue = true)
513  {
514  CopyModelPartFlaggedNodalNonHistoricalVarToNonHistoricalVar(
515  rOriginVariable, rDestinationVariable, rModelPart, rModelPart, rFlag, CheckValue);
516  }
517 
518  template <class TDataType>
520  const Variable<TDataType>& rVariable,
521  const ModelPart& rOriginModelPart,
522  ModelPart& rDestinationModelPart,
523  const Flags& rFlag,
524  const bool CheckValue = true)
525  {
526  CopyModelPartFlaggedNodalNonHistoricalVarToNonHistoricalVar(
527  rVariable, rVariable, rOriginModelPart, rDestinationModelPart, rFlag, CheckValue);
528  }
529 
530  template <class TDataType>
532  const Variable<TDataType>& rOriginVariable,
533  const Variable<TDataType>& rDestinationVariable,
534  const ModelPart& rOriginModelPart,
535  ModelPart& rDestinationModelPart,
536  const Flags& rFlag,
537  const bool CheckValue = true)
538  {
539  KRATOS_TRY
540 
541  KRATOS_ERROR_IF(rOriginModelPart.FullName() == rDestinationModelPart.FullName() && rOriginVariable == rDestinationVariable)
542  << "Trying to copy flagged elemental variable data with the same model "
543  "parts/variables. This is not permitted ( Origin model part: "
544  << rOriginModelPart.Name() << ", destination model part: " << rDestinationModelPart.Name()
545  << ", variable: " << rOriginVariable.Name() << " ) !";
546 
547  CopyModelPartFlaggedVariable<ElementsContainerType>(
548  rOriginModelPart, rDestinationModelPart, rFlag, CheckValue,
549  [&](ElementType& rDestElement, const TDataType& rValue) {
550  rDestElement.SetValue(rDestinationVariable, rValue);
551  },
552  [&](const ElementType& rOriginElement) -> const TDataType& {
553  return rOriginElement.GetValue(rOriginVariable);
554  });
555 
556  KRATOS_CATCH("");
557  }
558 
559  template <class TDataType>
561  const Variable<TDataType>& rOriginVariable,
562  const Variable<TDataType>& rDestinationVariable,
563  ModelPart& rModelPart,
564  const Flags& rFlag,
565  const bool CheckValue = true)
566  {
567  CopyModelPartFlaggedElementVar(
568  rOriginVariable, rDestinationVariable, rModelPart, rModelPart, rFlag, CheckValue);
569  }
570 
571  template <class TDataType>
573  const Variable<TDataType>& rVariable,
574  const ModelPart& rOriginModelPart,
575  ModelPart& rDestinationModelPart,
576  const Flags& rFlag,
577  const bool CheckValue = true)
578  {
579  CopyModelPartFlaggedElementVar(
580  rVariable, rVariable, rOriginModelPart, rDestinationModelPart, rFlag, CheckValue);
581  }
582 
583  template <class TDataType>
585  const Variable<TDataType>& rOriginVariable,
586  const Variable<TDataType>& rDestinationVariable,
587  const ModelPart& rOriginModelPart,
588  ModelPart& rDestinationModelPart,
589  const Flags& rFlag,
590  const bool CheckValue = true)
591  {
592  KRATOS_TRY
593 
594  KRATOS_ERROR_IF(rOriginModelPart.FullName() == rDestinationModelPart.FullName() && rOriginVariable == rDestinationVariable)
595  << "Trying to copy flagged condition variable data with the same model "
596  "parts/variables. This is not permitted ( Origin model part: "
597  << rOriginModelPart.Name() << ", destination model part: " << rDestinationModelPart.Name()
598  << ", variable: " << rOriginVariable.Name() << " ) !";
599 
600  CopyModelPartFlaggedVariable<ConditionsContainerType>(
601  rOriginModelPart, rDestinationModelPart, rFlag, CheckValue,
602  [&](ConditionType& rDestCondition, const TDataType& rValue) {
603  rDestCondition.SetValue(rDestinationVariable, rValue);
604  },
605  [&](const ConditionType& rOriginCondition) -> const TDataType& {
606  return rOriginCondition.GetValue(rOriginVariable);
607  });
608 
609  KRATOS_CATCH("");
610  }
611 
612  template <class TDataType>
614  const Variable<TDataType>& rOriginVariable,
615  const Variable<TDataType>& rDestinationVariable,
616  ModelPart& rModelPart,
617  const Flags& rFlag,
618  const bool CheckValue = true)
619  {
620  CopyModelPartFlaggedConditionVar(
621  rOriginVariable, rDestinationVariable, rModelPart, rModelPart, rFlag, CheckValue);
622  }
623 
624  template <class TDataType>
626  const Variable<TDataType>& rVariable,
627  const ModelPart& rOriginModelPart,
628  ModelPart& rDestinationModelPart,
629  const Flags& rFlag,
630  const bool CheckValue = true)
631  {
632  CopyModelPartFlaggedConditionVar(
633  rVariable, rVariable, rOriginModelPart, rDestinationModelPart, rFlag, CheckValue);
634  }
635 
645  template< class TVarType >
647  const TVarType& rVariable,
648  const ModelPart& rOriginModelPart,
649  ModelPart& rDestinationModelPart){
650 
651  const int n_orig_elems = rOriginModelPart.NumberOfElements();
652  const int n_dest_elems = rDestinationModelPart.NumberOfElements();
653 
654  KRATOS_ERROR_IF_NOT(n_orig_elems == n_dest_elems) << "Origin and destination model parts have different number of elements."
655  << "\n\t- Number of origin elements: " << n_orig_elems
656  << "\n\t- Number of destination elements: " << n_dest_elems << std::endl;
657 
658  IndexPartition<std::size_t>(n_orig_elems).for_each([&](std::size_t index){
659  auto it_dest_elems = rDestinationModelPart.ElementsBegin() + index;
660  const auto it_orig_elems = rOriginModelPart.ElementsBegin() + index;
661  const auto& r_value = it_orig_elems->GetValue(rVariable);
662  it_dest_elems->SetValue(rVariable,r_value);
663  });
664  }
665 
674  template<class TDataType, class TVarType = Variable<TDataType> >
676  const TVarType& rVariable,
677  const TDataType& rValue,
678  NodesContainerType& rNodes,
679  const unsigned int Step = 0)
680  {
681  KRATOS_TRY
682 
683  block_for_each(rNodes, [&](Node& rNode) {
684  rNode.FastGetSolutionStepValue(rVariable, Step) = rValue;
685  });
686 
687  KRATOS_CATCH("")
688  }
689 
700  template <class TDataType, class TVarType = Variable<TDataType>>
702  const TVarType &rVariable,
703  const TDataType &rValue,
704  NodesContainerType &rNodes,
705  const Flags Flag,
706  const bool CheckValue = true)
707  {
708  KRATOS_TRY
709 
710  block_for_each(rNodes, [&](Node& rNode){
711  if(rNode.Is(Flag) == CheckValue){
712  rNode.FastGetSolutionStepValue(rVariable) = rValue;}
713  });
714 
715  KRATOS_CATCH("")
716  }
717 
723  template< class TType , class TContainerType>
725  const Variable< TType >& rVariable,
726  TContainerType& rContainer)
727  {
728  KRATOS_TRY
729  this->SetNonHistoricalVariable(rVariable, rVariable.Zero(), rContainer);
730  KRATOS_CATCH("")
731  }
732 
741  template<class TContainerType, class... TVariableArgs>
743  TContainerType& rContainer,
744  const TVariableArgs&... rVariableArgs)
745  {
746  block_for_each(rContainer, [&](auto& rEntity){
747  (rEntity.SetValue(rVariableArgs, rVariableArgs.Zero()), ...);
748  });
749  }
750 
756  template< class TType >
758  const Variable< TType >& rVariable,
759  NodesContainerType& rNodes)
760  {
761  KRATOS_TRY
762  this->SetVariable(rVariable, rVariable.Zero(), rNodes);
763  KRATOS_CATCH("")
764  }
765 
773  template<class... TVariableArgs>
775  NodesContainerType& rNodes,
776  const TVariableArgs&... rVariableArgs)
777  {
778  block_for_each(rNodes, [&](NodeType& rNode){(
779  AuxiliaryHistoricalValueSetter<typename TVariableArgs::Type>(rVariableArgs, rVariableArgs.Zero(), rNode), ...);
780  });
781  }
782 
789  template< class TType, class TContainerType, class TVarType = Variable< TType >>
791  const TVarType& rVariable,
792  const TType& Value,
793  TContainerType& rContainer
794  )
795  {
796  KRATOS_TRY
797 
798  block_for_each(rContainer, [&](typename TContainerType::value_type& rEntity){
799  rEntity.SetValue(rVariable, Value);
800  });
801 
802  KRATOS_CATCH("")
803  }
804 
813  template< class TType, class TContainerType, class TVarType = Variable< TType >>
815  const TVarType& rVariable,
816  const TType& rValue,
817  TContainerType& rContainer,
818  const Flags Flag,
819  const bool Check = true
820  )
821  {
822  KRATOS_TRY
823 
824  block_for_each(rContainer, [&](typename TContainerType::value_type& rEntity){
825  if(rEntity.Is(Flag) == Check){
826  rEntity.SetValue(rVariable, rValue);}
827  });
828 
829  KRATOS_CATCH("")
830  }
831 
837  template< class TContainerType, class TVarType>
839  const TVarType& rVariable,
840  TContainerType& rContainer
841  )
842  {
843  KRATOS_TRY
844 
845  block_for_each(rContainer, [&rVariable](auto& rEntity){
846  rEntity.GetData().Erase(rVariable);
847  });
848 
849  KRATOS_CATCH("")
850  }
851 
856  template< class TContainerType>
857  void ClearNonHistoricalData(TContainerType& rContainer)
858  {
859  KRATOS_TRY
860 
861  block_for_each(rContainer, [&](typename TContainerType::value_type& rEntity){
862  rEntity.GetData().Clear();
863  });
864 
865  KRATOS_CATCH("")
866  }
867 
886  template <class TDataType, class TContainerType, class TWeightDataType>
887  void WeightedAccumulateVariableOnNodes(
888  ModelPart& rModelPart,
889  const Variable<TDataType>& rVariable,
890  const Variable<TWeightDataType>& rWeightVariable,
891  const bool IsInverseWeightProvided = false);
892 
899  template< class TContainerType >
900  void SetFlag(
901  const Flags& rFlag,
902  const bool FlagValue,
903  TContainerType& rContainer
904  )
905  {
906  KRATOS_TRY
907 
908  block_for_each(rContainer, [&](typename TContainerType::value_type& rEntity){
909  rEntity.Set(rFlag, FlagValue);
910  });
911 
912  KRATOS_CATCH("")
913 
914  }
915 
921  template< class TContainerType >
922  void ResetFlag(
923  const Flags& rFlag,
924  TContainerType& rContainer
925  )
926  {
927  KRATOS_TRY
928 
929  block_for_each(rContainer, [&](typename TContainerType::value_type& rEntity){
930  rEntity.Reset(rFlag);
931  });
932 
933  KRATOS_CATCH("")
934  }
935 
941  template< class TContainerType >
942  void FlipFlag(
943  const Flags& rFlag,
944  TContainerType& rContainer
945  )
946  {
947  KRATOS_TRY
948 
949  block_for_each(rContainer, [&](typename TContainerType::value_type& rEntity){
950  rEntity.Flip(rFlag);
951  });
952 
953  KRATOS_CATCH("")
954  }
955 
965  template< class TDataType, class TVariableType = Variable<TDataType> >
967  const TVariableType &rOriginVariable,
968  const TVariableType &rSavedVariable,
969  NodesContainerType &rNodesContainer)
970  {
971  KRATOS_TRY
972 
973  block_for_each(rNodesContainer, [&](Node& rNode){
974  rNode.SetValue(rSavedVariable, rNode.FastGetSolutionStepValue(rOriginVariable));
975  });
976 
977  KRATOS_CATCH("")
978  }
979 
990  template< class TDataType, class TContainerType, class TVariableType = Variable<TDataType> >
992  const TVariableType &rOriginVariable,
993  const TVariableType &rSavedVariable,
994  TContainerType &rContainer
995  )
996  {
997  KRATOS_TRY
998 
999  block_for_each(rContainer, [&](typename TContainerType::value_type& rEntity){
1000  rEntity.SetValue(rSavedVariable, rEntity.GetValue(rOriginVariable));
1001  });
1002 
1003  KRATOS_CATCH("")
1004  }
1005 
1016  template< class TDataType, class TVariableType = Variable<TDataType> >
1018  const TVariableType &rOriginVariable,
1019  const TVariableType &rDestinationVariable,
1020  NodesContainerType &rNodesContainer)
1021  {
1022  KRATOS_TRY
1023 
1024  block_for_each(rNodesContainer, [&](Node& rNode){
1025  rNode.FastGetSolutionStepValue(rDestinationVariable) = rNode.FastGetSolutionStepValue(rOriginVariable);
1026  });
1027 
1028  KRATOS_CATCH("")
1029  }
1030 
1038  [[nodiscard]] NodesContainerType SelectNodeList(
1039  const DoubleVarType& Variable,
1040  const double Value,
1041  const NodesContainerType& rOriginNodes
1042  );
1043 
1050  template<class TVarType>
1052  const TVarType& rVariable,
1053  const NodesContainerType& rNodes
1054  )
1055  {
1056  KRATOS_TRY
1057 
1058  for (auto& i_node : rNodes)
1059  KRATOS_CHECK_VARIABLE_IN_NODAL_DATA(rVariable, i_node);
1060 
1061  return 0;
1062 
1063  KRATOS_CATCH("");
1064  }
1065 
1072  template< class TVarType >
1074  const TVarType& rVar,
1075  const bool IsFixed,
1076  NodesContainerType& rNodes
1077  )
1078  {
1079  KRATOS_TRY
1080 
1081  if (rNodes.size() != 0) {
1082  // checking the first node to avoid error being thrown in parallel region
1083  KRATOS_ERROR_IF_NOT(rNodes.begin()->HasDofFor(rVar)) << "Trying to fix/free dof of variable " << rVar.Name() << " but this dof does not exist in node #" << rNodes.begin()->Id() << "!" << std::endl;
1084 
1085 #ifdef KRATOS_DEBUG
1086  for (const auto& r_node : rNodes) {
1087  KRATOS_ERROR_IF_NOT(r_node.HasDofFor(rVar)) << "Trying to fix/free dof of variable " << rVar.Name() << " but this dof does not exist in node #" << r_node.Id() << "!" << std::endl;
1088  }
1089 #endif
1090 
1091  CheckVariableExists(rVar, rNodes);
1092 
1093  if (IsFixed) {
1094  block_for_each(rNodes,[&](Node& rNode){
1095  rNode.pGetDof(rVar)->FixDof();
1096  });
1097  } else {
1098  block_for_each(rNodes,[&](Node& rNode){
1099  rNode.pGetDof(rVar)->FreeDof();
1100  });
1101  }
1102  }
1103 
1104  KRATOS_CATCH("")
1105  }
1106 
1120  template< class TVarType >
1122  const TVarType& rVariable,
1123  const bool IsFixed,
1124  NodesContainerType& rNodes,
1125  const Flags& rFlag,
1126  const bool CheckValue = true)
1127  {
1128  KRATOS_TRY
1129 
1130  if (rNodes.size() != 0) {
1131  // checking the first node to avoid error being thrown in parallel region
1132  KRATOS_ERROR_IF_NOT(rNodes.begin()->HasDofFor(rVariable))
1133  << "Trying to fix/free dof of variable " << rVariable.Name()
1134  << " but this dof does not exist in node #"
1135  << rNodes.begin()->Id() << "!" << std::endl;
1136 
1137 #ifdef KRATOS_DEBUG
1138  for (const auto& r_node : rNodes) {
1139  KRATOS_ERROR_IF_NOT(r_node.HasDofFor(rVariable))
1140  << "Trying to fix/free dof of variable " << rVariable.Name()
1141  << " but this dof does not exist in node #" << r_node.Id()
1142  << "!" << std::endl;
1143  }
1144 #endif
1145 
1146  CheckVariableExists(rVariable, rNodes);
1147 
1148  if (IsFixed) {
1149  block_for_each(rNodes, [&rVariable, &rFlag, CheckValue](NodeType& rNode) {
1150  if (rNode.Is(rFlag) == CheckValue) {
1151  rNode.pGetDof(rVariable)->FixDof();
1152  }
1153  });
1154  }
1155  else {
1156  block_for_each(rNodes, [&rVariable, &rFlag, CheckValue](NodeType& rNode) {
1157  if (rNode.Is(rFlag) == CheckValue) {
1158  rNode.pGetDof(rVariable)->FreeDof();
1159  }
1160  });
1161  }
1162  }
1163 
1164  KRATOS_CATCH("");
1165  }
1166 
1167 
1177  template< class TVarType >
1179  const TVarType& rVar,
1180  const Vector& rData,
1181  NodesContainerType& rNodes
1182  )
1183  {
1184  KRATOS_TRY
1185 
1186  if(rNodes.size() != 0 && rNodes.size() == rData.size()) {
1187  // First we do a check
1188  CheckVariableExists(rVar, rNodes);
1189 
1190  IndexPartition<std::size_t>(rNodes.size()).for_each([&](std::size_t index){
1191  NodesContainerType::iterator it_node = rNodes.begin() + index;
1192  it_node->FastGetSolutionStepValue(rVar) = rData[index];
1193  });
1194  } else
1195  KRATOS_ERROR << "There is a mismatch between the size of data array and the number of nodes ";
1196 
1197  KRATOS_CATCH("")
1198  }
1199 
1206  [[nodiscard]] array_1d<double, 3> SumNonHistoricalNodeVectorVariable(
1207  const ArrayVarType& rVar,
1208  const ModelPart& rModelPart
1209  );
1210 
1217  template< class TVarType >
1219  const TVarType& rVar,
1220  const ModelPart& rModelPart
1221  )
1222  {
1223  KRATOS_TRY
1224 
1225  double sum_value = 0.0;
1226 
1227  // Getting info
1228  const auto& r_communicator = rModelPart.GetCommunicator();
1229  const auto& r_local_mesh = r_communicator.LocalMesh();
1230  const auto& r_nodes_array = r_local_mesh.Nodes();
1231 
1232  sum_value = block_for_each<SumReduction<double>>(r_nodes_array, [&](Node& rNode){
1233  return rNode.GetValue(rVar);
1234  });
1235 
1236  return r_communicator.GetDataCommunicator().SumAll(sum_value);
1237 
1238  KRATOS_CATCH("")
1239  }
1240 
1252  template< class TDataType, class TVarType = Variable<TDataType> >
1253  [[nodiscard]] TDataType SumHistoricalVariable(
1254  const TVarType &rVariable,
1255  const ModelPart &rModelPart,
1256  const unsigned int BuffStep = 0
1257  )
1258  {
1259  KRATOS_TRY
1260 
1261  const auto &r_communicator = rModelPart.GetCommunicator();
1262 
1263  TDataType sum_value = block_for_each<SumReduction<TDataType>>(r_communicator.LocalMesh().Nodes(),[&](Node& rNode){
1264  return rNode.GetSolutionStepValue(rVariable, BuffStep);
1265  });
1266 
1267  return r_communicator.GetDataCommunicator().SumAll(sum_value);
1268 
1269  KRATOS_CATCH("")
1270 
1271  }
1272 
1279  [[nodiscard]] array_1d<double, 3> SumConditionVectorVariable(
1280  const ArrayVarType& rVar,
1281  const ModelPart& rModelPart
1282  );
1283 
1290  template< class TVarType >
1291  [[nodiscard]] double SumConditionScalarVariable(
1292  const TVarType& rVar,
1293  const ModelPart& rModelPart
1294  )
1295  {
1296  KRATOS_TRY
1297 
1298  double sum_value = 0.0;
1299 
1300  // Getting info
1301  const auto& r_communicator = rModelPart.GetCommunicator();
1302  const auto& r_local_mesh = r_communicator.LocalMesh();
1303  const auto& r_conditions_array = r_local_mesh.Conditions();
1304 
1305  sum_value = block_for_each<SumReduction<double>>(r_conditions_array, [&](ConditionType& rCond){
1306  return rCond.GetValue(rVar);
1307  });
1308 
1309  return r_communicator.GetDataCommunicator().SumAll(sum_value);
1310 
1311  KRATOS_CATCH("")
1312  }
1313 
1320  array_1d<double, 3> SumElementVectorVariable(
1321  const ArrayVarType& rVar,
1322  const ModelPart& rModelPart
1323  );
1324 
1331  template< class TVarType >
1332  [[nodiscard]] double SumElementScalarVariable(
1333  const TVarType& rVar,
1334  const ModelPart& rModelPart
1335  )
1336  {
1337  KRATOS_TRY
1338 
1339  double sum_value = 0.0;
1340 
1341  // Getting info
1342  const auto& r_communicator = rModelPart.GetCommunicator();
1343  const auto& r_local_mesh = r_communicator.LocalMesh();
1344  const auto& r_elements_array = r_local_mesh.Elements();
1345 
1346  sum_value = block_for_each<SumReduction<double>>(r_elements_array, [&](ElementType& rElem){
1347  return rElem.GetValue(rVar);
1348  });
1349 
1350  return r_communicator.GetDataCommunicator().SumAll(sum_value);
1351 
1352  KRATOS_CATCH("")
1353  }
1354 
1360  template< class TVarType >
1361  void AddDof(
1362  const TVarType& rVar,
1363  ModelPart& rModelPart
1364  )
1365  {
1366  KRATOS_TRY
1367 
1368  // First we do a chek
1369  if(rModelPart.NumberOfNodes() != 0)
1370  KRATOS_ERROR_IF_NOT(rModelPart.NodesBegin()->SolutionStepsDataHas(rVar)) << "ERROR:: Variable : " << rVar << "not included in the Solution step data ";
1371 
1372  rModelPart.GetNodalSolutionStepVariablesList().AddDof(&rVar);
1373 
1374  block_for_each(rModelPart.Nodes(),[&](Node& rNode){
1375  rNode.AddDof(rVar);
1376  });
1377 
1378  KRATOS_CATCH("")
1379  }
1380 
1387  template< class TVarType >
1389  const TVarType& rVar,
1390  const TVarType& rReactionVar,
1391  ModelPart& rModelPart
1392  )
1393  {
1394  KRATOS_TRY
1395 
1396  if(rModelPart.NumberOfNodes() != 0) {
1397  KRATOS_ERROR_IF_NOT(rModelPart.NodesBegin()->SolutionStepsDataHas(rVar)) << "ERROR:: DoF Variable : " << rVar << "not included in the Soluttion step data ";
1398  KRATOS_ERROR_IF_NOT(rModelPart.NodesBegin()->SolutionStepsDataHas(rReactionVar)) << "ERROR:: Reaction Variable : " << rReactionVar << "not included in the Soluttion step data ";
1399  }
1400 
1401  // If in debug we do a check for all nodes
1402  #ifdef KRATOS_DEBUG
1403  CheckVariableExists(rVar, rModelPart.Nodes());
1404  CheckVariableExists(rReactionVar, rModelPart.Nodes());
1405  #endif
1406 
1407  rModelPart.GetNodalSolutionStepVariablesList().AddDof(&rVar, &rReactionVar);
1408 
1409  block_for_each(rModelPart.Nodes(),[&](Node& rNode){
1410  rNode.AddDof(rVar,rReactionVar);
1411  });
1412 
1413  KRATOS_CATCH("")
1414  }
1415 
1423  static void AddDofsList(
1424  const std::vector<std::string>& rDofsVarNamesList,
1425  ModelPart& rModelPart);
1426 
1434  static void AddDofsWithReactionsList(
1435  const std::vector<std::array<std::string,2>>& rDofsAndReactionsNamesList,
1436  ModelPart& rModelPart);
1437 
1442  bool CheckVariableKeys();
1443 
1448  void UpdateCurrentToInitialConfiguration(const ModelPart::NodesContainerType& rNodes);
1449 
1454  void UpdateInitialToCurrentConfiguration(const ModelPart::NodesContainerType& rNodes);
1455 
1463  void UpdateCurrentPosition(
1464  const ModelPart::NodesContainerType& rNodes,
1465  const ArrayVarType& rUpdateVariable = DISPLACEMENT,
1466  const IndexType BufferPosition = 0
1467  );
1468 
1480  template<class TVectorType=Vector>
1481  [[nodiscard]] TVectorType GetCurrentPositionsVector(
1482  const ModelPart::NodesContainerType& rNodes,
1483  const unsigned int Dimension
1484  );
1485 
1497  template<class TVectorType=Vector>
1498  [[nodiscard]] TVectorType GetInitialPositionsVector(
1499  const ModelPart::NodesContainerType& rNodes,
1500  const unsigned int Dimension
1501  );
1502 
1515  void SetCurrentPositionsVector(
1517  const Vector& rPositions
1518  );
1519 
1532  void SetInitialPositionsVector(
1534  const Vector& rPositions
1535  );
1536 
1540  template <Globals::DataLocation TLocation, class TEntity, class TValue>
1541  static bool HasValue(const TEntity& rEntity, const Variable<TValue>& rVariable)
1542  {
1543  if constexpr (TLocation == Globals::DataLocation::NodeHistorical) {
1544  static_assert(std::is_same_v<TEntity,Node>);
1545  return rEntity.SolutionStepsDataHas(rVariable);
1546  } else {
1547  static_assert(std::is_same_v<TEntity,Node>
1548  || std::is_same_v<TEntity,Element>
1549  || std::is_same_v<TEntity,Condition>
1550  || std::is_same_v<TEntity,ProcessInfo>
1551  || std::is_same_v<TEntity,ModelPart>);
1552  return rEntity.Has(rVariable);
1553  }
1554  }
1555 
1559  template <Globals::DataLocation TLocation, class TEntity, class TValue>
1560  static std::conditional_t<std::is_arithmetic_v<TValue>,
1561  TValue, // <== return by value if scalar type
1562  const TValue&> // <== return by reference if non-scalar type
1563  GetValue(const TEntity& rEntity, const Variable<TValue>& rVariable)
1564  {
1565  if constexpr (TLocation == Globals::DataLocation::NodeHistorical) {
1566  static_assert(std::is_same_v<TEntity,Node>);
1567  return rEntity.FastGetSolutionStepValue(rVariable);
1568  } else if constexpr (TLocation == Globals::DataLocation::NodeNonHistorical) {
1569  static_assert(std::is_same_v<TEntity,Node>);
1570  return rEntity.GetValue(rVariable);
1571  } else if constexpr (TLocation == Globals::DataLocation::Element) {
1572  static_assert(std::is_same_v<TEntity,Element>);
1573  return rEntity.GetValue(rVariable);
1574  } else if constexpr (TLocation == Globals::DataLocation::Condition) {
1575  static_assert(std::is_same_v<TEntity,Condition>);
1576  return rEntity.GetValue(rVariable);
1577  } else if constexpr (TLocation == Globals::DataLocation::ProcessInfo) {
1578  static_assert(std::is_same_v<TEntity,ProcessInfo>);
1579  return rEntity.GetValue(rVariable);
1580  } else if constexpr (TLocation == Globals::DataLocation::ModelPart) {
1581  static_assert(std::is_same_v<TEntity,ModelPart>);
1582  return rEntity.GetValue(rVariable);
1583  } else {
1584  static_assert(std::is_same_v<TEntity,void>, "Unsupported DataLocation");
1585  }
1586  }
1587 
1591  template <Globals::DataLocation TLocation, class TEntity, class TValue>
1592  static TValue& GetValue(TEntity& rEntity, const Variable<TValue>& rVariable)
1593  {
1594  if constexpr (TLocation == Globals::DataLocation::NodeHistorical) {
1595  static_assert(std::is_same_v<TEntity,Node>);
1596  return rEntity.FastGetSolutionStepValue(rVariable);
1597  } else if constexpr (TLocation == Globals::DataLocation::NodeNonHistorical) {
1598  static_assert(std::is_same_v<TEntity,Node>);
1599  return rEntity.GetValue(rVariable);
1600  } else if constexpr (TLocation == Globals::DataLocation::Element) {
1601  static_assert(std::is_same_v<TEntity,Element>);
1602  return rEntity.GetValue(rVariable);
1603  } else if constexpr (TLocation == Globals::DataLocation::Condition) {
1604  static_assert(std::is_same_v<TEntity,Condition>);
1605  return rEntity.GetValue(rVariable);
1606  } else if constexpr (TLocation == Globals::DataLocation::ProcessInfo) {
1607  static_assert(std::is_same_v<TEntity,ProcessInfo>);
1608  return rEntity.GetValue(rVariable);
1609  } else if constexpr (TLocation == Globals::DataLocation::ModelPart) {
1610  static_assert(std::is_same_v<TEntity,ModelPart>);
1611  return rEntity.GetValue(rVariable);
1612  } else {
1613  static_assert(std::is_same_v<TEntity,void>, "Unsupported DataLocation");
1614  }
1615  }
1616 
1621  template <Globals::DataLocation TLocation, class TEntity, class TValue>
1622  static void SetValue(TEntity& rEntity,
1623  const Variable<TValue>& rVariable,
1624  std::conditional_t<std::is_arithmetic_v<TValue>,
1625  TValue, /*pass scalar types by value*/
1626  const TValue&> /*pass non-scalar types by reference*/ Value)
1627  {
1628  if constexpr (TLocation == Globals::DataLocation::NodeHistorical) {
1629  static_assert(std::is_same_v<TEntity,Node>);
1630  rEntity.FastGetSolutionStepValue(rVariable) = Value;
1631  } else if constexpr (TLocation == Globals::DataLocation::NodeNonHistorical) {
1632  static_assert(std::is_same_v<TEntity,Node>);
1633  rEntity.SetValue(rVariable, Value);
1634  } else if constexpr (TLocation == Globals::DataLocation::Element) {
1635  static_assert(std::is_same_v<TEntity,Element>);
1636  rEntity.SetValue(rVariable, Value);
1637  } else if constexpr (TLocation == Globals::DataLocation::Condition) {
1638  static_assert(std::is_same_v<TEntity,Condition>);
1639  rEntity.SetValue(rVariable, Value);
1640  } else if constexpr (TLocation == Globals::DataLocation::ProcessInfo) {
1641  static_assert(std::is_same_v<TEntity,ProcessInfo>);
1642  rEntity.SetValue(rVariable, Value);
1643  } else if constexpr (TLocation == Globals::DataLocation::ModelPart) {
1644  static_assert(std::is_same_v<TEntity,ModelPart>);
1645  rEntity.SetValue(rVariable, Value);
1646  } else {
1647  static_assert(std::is_same_v<TEntity,void>, "Unsupported DataLocation");
1648  }
1649  }
1650 
1662  [[nodiscard]] Vector GetSolutionStepValuesVector(
1663  const ModelPart::NodesContainerType& rNodes,
1664  const Variable<array_1d<double,3>>& rVar,
1665  const unsigned int Step,
1666  const unsigned int Dimension=3
1667  );
1668 
1679  [[nodiscard]] Vector GetSolutionStepValuesVector(
1680  const ModelPart::NodesContainerType& rNodes,
1681  const Variable<double>& rVar,
1682  const unsigned int Step
1683  );
1684 
1694  void SetSolutionStepValuesVector(
1696  const Variable<array_1d<double,3>>& rVar,
1697  const Vector& rData,
1698  const unsigned int Step
1699  );
1700 
1710  void SetSolutionStepValuesVector(
1712  const Variable<double>& rVar,
1713  const Vector& rData,
1714  const unsigned int Step
1715  );
1716 
1729  [[nodiscard]] Vector GetValuesVector(
1730  const ModelPart::NodesContainerType& rNodes,
1731  const Variable<array_1d<double,3>>& rVar,
1732  const unsigned int Dimension=3
1733  );
1734 
1746  [[nodiscard]] Vector GetValuesVector(
1747  const ModelPart::NodesContainerType& rNodes,
1748  const Variable<double>& rVar
1749  );
1750 
1761  void SetValuesVector(
1763  const Variable<array_1d<double,3>>& rVar,
1764  const Vector& rData
1765  );
1766 
1777  void SetValuesVector(
1779  const Variable<double>& rVar,
1780  const Vector& rData
1781  );
1782 
1783 
1787 
1788 
1792 
1793 
1797 
1798 
1800 
1801 private:
1804 
1808 
1812 
1813 
1817 
1822  template< class TVarType >
1823  bool CheckVariableKeysHelper()
1824  {
1825  KRATOS_TRY
1826 
1827  for (const auto& var : KratosComponents< TVarType >::GetComponents()) {
1828  if (var.first == "NONE" || var.first == "")
1829  std::cout << " var first is NONE or empty " << var.first << var.second << std::endl;
1830  if (var.second->Name() == "NONE" || var.second->Name() == "")
1831  std::cout << var.first << var.second << std::endl;
1832  if (var.first != var.second->Name()) //name of registration does not correspond to the var name
1833  std::cout << "Registration Name = " << var.first << " Variable Name = " << std::endl;
1834  }
1835 
1836  return true;
1837  KRATOS_CATCH("")
1838  }
1839 
1840  template <class TContainerType>
1841  [[nodiscard]] TContainerType& GetContainer(ModelPart& rModelPart);
1842 
1843  template <class TContainerType>
1844  [[nodiscard]] const TContainerType& GetContainer(const ModelPart& rModelPart);
1845 
1846  template<class TDataType>
1847  static void AuxiliaryHistoricalValueSetter(
1848  const Variable<TDataType>& rVariable,
1849  const TDataType& rValue,
1850  NodeType& rNode);
1851 
1852  template <class TContainerType, class TSetterFunction, class TGetterFunction>
1853  void CopyModelPartFlaggedVariable(
1854  const ModelPart& rOriginModelPart,
1855  ModelPart& rDestinationModelPart,
1856  const Flags& rFlag,
1857  const bool CheckValue,
1858  TSetterFunction&& rSetterFunction,
1859  TGetterFunction&& rGetterFunction)
1860  {
1861  KRATOS_TRY
1862 
1863  const auto& r_origin_container = GetContainer<TContainerType>(rOriginModelPart);
1864  auto& r_destination_container = GetContainer<TContainerType>(rDestinationModelPart);
1865 
1866  const int number_of_origin_items = r_origin_container.size();
1867  const int number_of_destination_items = r_destination_container.size();
1868 
1869  KRATOS_ERROR_IF_NOT(number_of_origin_items == number_of_destination_items)
1870  << "Origin ( " << rOriginModelPart.Name() << " ) and destination ( "
1871  << rDestinationModelPart.Name() << " ) model parts have different number of items."
1872  << "\n\t- Number of origin items: " << number_of_origin_items
1873  << "\n\t- Number of destination items: " << number_of_destination_items
1874  << std::endl;
1875 
1876  IndexPartition<int>(number_of_origin_items).for_each([&](int i_node) {
1877  const auto& r_orig_item = *(r_origin_container.begin() + i_node);
1878  auto& r_dest_item = *(r_destination_container.begin() + i_node);
1879  if (r_orig_item.Is(rFlag) == CheckValue) {
1880  rSetterFunction(r_dest_item, rGetterFunction(r_orig_item));
1881  }
1882  });
1883 
1884  KRATOS_CATCH("");
1885  }
1886 
1890 
1891 
1895 
1896 
1900 
1901 
1903 
1904 }; /* Class VariableUtils */
1905 
1907 
1910 
1911 
1913 
1914 } /* namespace Kratos.*/
virtual bool SynchronizeNonHistoricalVariable(Variable< int > const &rThisVariable)
Definition: communicator.cpp:407
virtual bool SynchronizeVariable(Variable< int > const &rThisVariable)
Definition: communicator.cpp:357
virtual const DataCommunicator & GetDataCommunicator() const
Definition: communicator.cpp:340
MeshType & LocalMesh()
Returns the reference to the mesh storing all local entities.
Definition: communicator.cpp:245
Base class for all Conditions.
Definition: condition.h:59
Base class for all Elements.
Definition: element.h:60
Definition: flags.h:58
bool Is(Flags const &rOther) const
Definition: flags.h:274
void SetValue(const TVariableType &rThisVariable, typename TVariableType::Type const &rValue)
Definition: geometrical_object.h:238
This class is useful for index iteration over containers.
Definition: parallel_utilities.h:451
void for_each(TUnaryFunction &&f)
Definition: parallel_utilities.h:514
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
ConditionsContainerType & Conditions()
Definition: mesh.h:691
NodesContainerType & Nodes()
Definition: mesh.h:346
ElementsContainerType & Elements()
Definition: mesh.h:568
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
ElementIterator ElementsBegin(IndexType ThisIndex=0)
Definition: model_part.h:1169
std::string FullName() const
This method returns the full name of the model part (including the parents model parts)
Definition: model_part.h:1850
Communicator & GetCommunicator()
Definition: model_part.h:1821
std::string & Name()
Definition: model_part.h:1811
IndexType GetBufferSize() const
This method gets the suffer size of the model part database.
Definition: model_part.h:1876
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
bool HasNodalSolutionStepVariable(VariableData const &ThisVariable) const
Definition: model_part.h:544
SizeType NumberOfElements(IndexType ThisIndex=0) const
Definition: model_part.h:1027
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
SizeType NumberOfNodes(IndexType ThisIndex=0) const
Definition: model_part.h:341
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
VariablesList & GetNodalSolutionStepVariablesList()
Definition: model_part.h:549
This class defines the node.
Definition: node.h:65
TVariableType::Type & FastGetSolutionStepValue(const TVariableType &rThisVariable)
Definition: node.h:435
void SetValue(const TVariableType &rThisVariable, typename TVariableType::Type const &rValue)
Definition: node.h:493
DofType::Pointer pGetDof(TVariableType const &rDofVariable) const
Get DoF counted pointer for a given variable.
Definition: node.h:711
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
boost::indirect_iterator< typename TContainerType::iterator > iterator
Definition: pointer_vector_set.h:95
size_type size() const
Returns the number of elements in the container.
Definition: pointer_vector_set.h:502
iterator begin()
Returns an iterator pointing to the beginning of the container.
Definition: pointer_vector_set.h:278
const std::string & Name() const
Definition: variable_data.h:201
const TDataType & Zero() const
This method returns the zero value of the variable type.
Definition: variable.h:346
This class implements a set of auxiliar, already parallelized, methods to perform some common tasks r...
Definition: variable_utils.h:63
ModelPart::ElementType ElementType
The element type.
Definition: variable_utils.h:75
void CopyVariable(const TVariableType &rOriginVariable, const TVariableType &rDestinationVariable, NodesContainerType &rNodesContainer)
Takes the value of an historical variable and sets it in another variable This function takes the val...
Definition: variable_utils.h:1017
void CopyModelPartFlaggedNodalHistoricalVarToNonHistoricalVar(const Variable< TDataType > &rVariable, ModelPart &rModelPart, const Flags &rFlag, const bool CheckValue=true, const unsigned int ReadBufferStep=0)
Definition: variable_utils.h:385
void CopyModelPartFlaggedConditionVar(const Variable< TDataType > &rVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const Flags &rFlag, const bool CheckValue=true)
Definition: variable_utils.h:625
double SumNonHistoricalNodeScalarVariable(const TVarType &rVar, const ModelPart &rModelPart)
Returns the nodal value summation of a non-historical scalar variable.
Definition: variable_utils.h:1218
void CopyModelPartNodalVar(const TVarType &rVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const unsigned int BuffStep=0)
Copies the nodal value of a variable from an origin model part nodes to the nodes in a destination mo...
Definition: variable_utils.h:183
void AddDofWithReaction(const TVarType &rVar, const TVarType &rReactionVar, ModelPart &rModelPart)
This function add dofs to the nodes in a model part. It is useful since addition is done in parallel.
Definition: variable_utils.h:1388
void CopyModelPartFlaggedNodalHistoricalVarToNonHistoricalVar(const Variable< TDataType > &rVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const Flags &rFlag, const bool CheckValue=true, const unsigned int ReadBufferStep=0)
Definition: variable_utils.h:371
void SetHistoricalVariableToZero(const Variable< TType > &rVariable, NodesContainerType &rNodes)
Sets the nodal value of any variable to zero.
Definition: variable_utils.h:757
void SetNonHistoricalVariable(const TVarType &rVariable, const TType &rValue, TContainerType &rContainer, const Flags Flag, const bool Check=true)
Sets the container value of any type of non historical variable (considering flag)
Definition: variable_utils.h:814
void CopyModelPartFlaggedElementVar(const Variable< TDataType > &rOriginVariable, const Variable< TDataType > &rDestinationVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const Flags &rFlag, const bool CheckValue=true)
Definition: variable_utils.h:531
void ApplyFixity(const TVarType &rVar, const bool IsFixed, NodesContainerType &rNodes)
Fixes or frees a variable for all of the nodes in the list. The dof has to exist.
Definition: variable_utils.h:1073
void CopyModelPartFlaggedNodalNonHistoricalVarToHistoricalVar(const Variable< TDataType > &rVariable, ModelPart &rModelPart, const Flags &rFlag, const bool CheckValue=true, const unsigned int WriteBufferStep=0)
Definition: variable_utils.h:462
static bool HasValue(const TEntity &rEntity, const Variable< TValue > &rVariable)
Check whether a Node, Element, Condition, ProcessInfo, or ModelPart stores a value for the provided V...
Definition: variable_utils.h:1541
ModelPart::ElementsContainerType ElementsContainerType
The elements container.
Definition: variable_utils.h:87
void CopyModelPartFlaggedNodalNonHistoricalVarToNonHistoricalVar(const Variable< TDataType > &rOriginVariable, const Variable< TDataType > &rDestinationVariable, ModelPart &rModelPart, const Flags &rFlag, const bool CheckValue=true)
Definition: variable_utils.h:507
ModelPart::NodesContainerType NodesContainerType
The nodes container.
Definition: variable_utils.h:81
Variable< double > DoubleVarType
A definition of the double variable.
Definition: variable_utils.h:90
void CopyModelPartFlaggedNodalHistoricalVarToHistoricalVar(const Variable< TDataType > &rVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const Flags &rFlag, const bool CheckValue=true, const unsigned int ReadBufferStep=0, const unsigned int WriteBufferStep=0)
Definition: variable_utils.h:302
int CheckVariableExists(const TVarType &rVariable, const NodesContainerType &rNodes)
Checks if all the nodes of a node set has the specified variable.
Definition: variable_utils.h:1051
void ClearNonHistoricalData(TContainerType &rContainer)
Clears the container data value container.
Definition: variable_utils.h:857
void ApplyVector(const TVarType &rVar, const Vector &rData, NodesContainerType &rNodes)
Loops along a vector data to set its values to the nodes contained in a node set.
Definition: variable_utils.h:1178
ModelPart::NodeType NodeType
The node type.
Definition: variable_utils.h:69
static TValue & GetValue(TEntity &rEntity, const Variable< TValue > &rVariable)
Fetch the value of a variable stored in an entity.
Definition: variable_utils.h:1592
void CopyModelPartFlaggedNodalHistoricalVarToNonHistoricalVar(const Variable< TDataType > &rOriginVariable, const Variable< TDataType > &rDestinationVariable, ModelPart &rModelPart, const Flags &rFlag, const bool CheckValue=true, const unsigned int ReadBufferStep=0)
Definition: variable_utils.h:357
static void SetNonHistoricalVariablesToZero(TContainerType &rContainer, const TVariableArgs &... rVariableArgs)
Set the Non Historical Variables To Zero This method sets the provided list of variables to zero in t...
Definition: variable_utils.h:742
void SaveVariable(const TVariableType &rOriginVariable, const TVariableType &rSavedVariable, NodesContainerType &rNodesContainer)
Takes the value of a non-historical variable and saves it in another variable For a nodal container,...
Definition: variable_utils.h:966
void SetNonHistoricalVariable(const TVarType &rVariable, const TType &Value, TContainerType &rContainer)
Sets the container value of any type of non historical variable.
Definition: variable_utils.h:790
void EraseNonHistoricalVariable(const TVarType &rVariable, TContainerType &rContainer)
Erases the container non historical variable.
Definition: variable_utils.h:838
void CopyModelPartFlaggedElementVar(const Variable< TDataType > &rVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const Flags &rFlag, const bool CheckValue=true)
Definition: variable_utils.h:572
void SaveNonHistoricalVariable(const TVariableType &rOriginVariable, const TVariableType &rSavedVariable, TContainerType &rContainer)
Takes the value of a non-historical variable and saves it in another historical variable For a non-no...
Definition: variable_utils.h:991
void SetNonHistoricalVariableToZero(const Variable< TType > &rVariable, TContainerType &rContainer)
Sets the nodal value of any variable to zero.
Definition: variable_utils.h:724
void AddDof(const TVarType &rVar, ModelPart &rModelPart)
This function add dofs to the nodes in a model part. It is useful since addition is done in parallel.
Definition: variable_utils.h:1361
static std::conditional_t< std::is_arithmetic_v< TValue >, TValue, const TValue & > GetValue(const TEntity &rEntity, const Variable< TValue > &rVariable)
Fetch the value of a variable stored in an entity.
Definition: variable_utils.h:1563
void SetVariable(const TVarType &rVariable, const TDataType &rValue, NodesContainerType &rNodes, const unsigned int Step=0)
Sets the nodal value of a scalar variable.
Definition: variable_utils.h:675
ModelPart::ConditionsContainerType ConditionsContainerType
The conditions container.
Definition: variable_utils.h:84
void CopyModelPartFlaggedNodalHistoricalVarToHistoricalVar(const Variable< TDataType > &rOriginVariable, const Variable< TDataType > &rDestinationVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const Flags &rFlag, const bool CheckValue=true, const unsigned int ReadBufferStep=0, const unsigned int WriteBufferStep=0)
Definition: variable_utils.h:229
void CopyModelPartNodalVar(const TVarType &rVariable, const TVarType &rDestinationVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const unsigned int BuffStep=0)
Copies the nodal value of a variable from an origin model part nodes to the nodes in a destination mo...
Definition: variable_utils.h:163
void CopyModelPartFlaggedNodalNonHistoricalVarToHistoricalVar(const Variable< TDataType > &rVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const Flags &rFlag, const bool CheckValue=true, const unsigned int WriteBufferStep=0)
Definition: variable_utils.h:448
static void SetValue(TEntity &rEntity, const Variable< TValue > &rVariable, std::conditional_t< std::is_arithmetic_v< TValue >, TValue, const TValue & > Value)
Overwrite the value of a variable stored in an entity.
Definition: variable_utils.h:1622
Variable< array_1d< double, 3 > > ArrayVarType
A definition of the array variable.
Definition: variable_utils.h:93
void ApplyFixity(const TVarType &rVariable, const bool IsFixed, NodesContainerType &rNodes, const Flags &rFlag, const bool CheckValue=true)
Fixes/Frees dofs based on a flag.
Definition: variable_utils.h:1121
void ResetFlag(const Flags &rFlag, TContainerType &rContainer)
Flips a flag over a given container.
Definition: variable_utils.h:922
void CopyModelPartNodalVarToNonHistoricalVar(const TVarType &rVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const unsigned int BuffStep=0)
Definition: variable_utils.h:219
void CopyModelPartFlaggedNodalNonHistoricalVarToNonHistoricalVar(const Variable< TDataType > &rOriginVariable, const Variable< TDataType > &rDestinationVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const Flags &rFlag, const bool CheckValue=true)
Definition: variable_utils.h:475
void SetVariable(const TVarType &rVariable, const TDataType &rValue, NodesContainerType &rNodes, const Flags Flag, const bool CheckValue=true)
Sets the nodal value of a scalar variable (considering flag)
Definition: variable_utils.h:701
void CopyModelPartNodalVar(const TVarType &rVariable, const TVarType &rDestinationVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const unsigned int ReadBufferStep, const unsigned int WriteBufferStep)
Copies the nodal value of a variable from an origin model part nodes to the nodes in a destination mo...
Definition: variable_utils.h:126
double SumConditionScalarVariable(const TVarType &rVar, const ModelPart &rModelPart)
Returns the condition value summation of a historical scalar variable.
Definition: variable_utils.h:1291
KRATOS_CLASS_POINTER_DEFINITION(VariableUtils)
We create the Pointer related to VariableUtils.
void CopyModelPartFlaggedConditionVar(const Variable< TDataType > &rOriginVariable, const Variable< TDataType > &rDestinationVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const Flags &rFlag, const bool CheckValue=true)
Definition: variable_utils.h:584
void CopyModelPartElementalVar(const TVarType &rVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart)
Copies the elemental value of a variable from an origin model part elements to the elements in a dest...
Definition: variable_utils.h:646
static void SetHistoricalVariablesToZero(NodesContainerType &rNodes, const TVariableArgs &... rVariableArgs)
Set the Historical Variables To Zero This method sets the provided list of variables to zero in the n...
Definition: variable_utils.h:774
void CopyModelPartFlaggedConditionVar(const Variable< TDataType > &rOriginVariable, const Variable< TDataType > &rDestinationVariable, ModelPart &rModelPart, const Flags &rFlag, const bool CheckValue=true)
Definition: variable_utils.h:613
void CopyModelPartFlaggedElementVar(const Variable< TDataType > &rOriginVariable, const Variable< TDataType > &rDestinationVariable, ModelPart &rModelPart, const Flags &rFlag, const bool CheckValue=true)
Definition: variable_utils.h:560
void SetFlag(const Flags &rFlag, const bool FlagValue, TContainerType &rContainer)
Sets a flag according to a given status over a given container.
Definition: variable_utils.h:900
void CopyModelPartFlaggedNodalHistoricalVarToHistoricalVar(const Variable< TDataType > &rOriginVariable, const Variable< TDataType > &rDestinationVariable, ModelPart &rModelPart, const Flags &rFlag, const bool CheckValue=true, const unsigned int ReadBufferStep=0, const unsigned int WriteBufferStep=0)
Definition: variable_utils.h:283
TDataType SumHistoricalVariable(const TVarType &rVariable, const ModelPart &rModelPart, const unsigned int BuffStep=0)
This method accumulates and return a variable value For a nodal historical variable,...
Definition: variable_utils.h:1253
void CopyModelPartFlaggedNodalNonHistoricalVarToNonHistoricalVar(const Variable< TDataType > &rVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const Flags &rFlag, const bool CheckValue=true)
Definition: variable_utils.h:519
void CopyModelPartFlaggedNodalNonHistoricalVarToHistoricalVar(const Variable< TDataType > &rOriginVariable, const Variable< TDataType > &rDestinationVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const Flags &rFlag, const bool CheckValue=true, const unsigned int WriteBufferStep=0)
Definition: variable_utils.h:398
void CopyModelPartNodalVarToNonHistoricalVar(const TVarType &rVariable, const TVarType &rDestinationVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const unsigned int BuffStep=0)
Definition: variable_utils.h:193
void CopyModelPartFlaggedNodalNonHistoricalVarToHistoricalVar(const Variable< TDataType > &rOriginVariable, const Variable< TDataType > &rDestinationVariable, ModelPart &rModelPart, const Flags &rFlag, const bool CheckValue=true, const unsigned int WriteBufferStep=0)
Definition: variable_utils.h:434
double SumElementScalarVariable(const TVarType &rVar, const ModelPart &rModelPart)
Returns the element value summation of a historical scalar variable.
Definition: variable_utils.h:1332
void FlipFlag(const Flags &rFlag, TContainerType &rContainer)
Flips a flag over a given container.
Definition: variable_utils.h:942
ModelPart::ConditionType ConditionType
The condition type.
Definition: variable_utils.h:72
void CopyModelPartFlaggedNodalHistoricalVarToNonHistoricalVar(const Variable< TDataType > &rOriginVariable, const Variable< TDataType > &rDestinationVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const Flags &rFlag, const bool CheckValue=true, const unsigned int ReadBufferStep=0)
Definition: variable_utils.h:321
int AddDof(VariableData const *pThisDofVariable)
Definition: variables_list.h:282
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
#define KRATOS_ERROR
Definition: exception.h:161
#define KRATOS_ERROR_IF_NOT(conditional)
Definition: exception.h:163
#define KRATOS_ERROR_IF(conditional)
Definition: exception.h:162
#define KRATOS_CHECK_VARIABLE_IN_NODAL_DATA(TheVariable, TheNode)
Definition: checks.h:171
TContainerType & GetContainer(ModelPart::MeshType &rMesh)
void CopyModelPartNodalVarToNonHistoricalVar(VariableUtils &rVariableUtils, const TVarType &rVariable, const ModelPart &rOriginModelPart, ModelPart &rDestinationModelPart, const unsigned int BuffStep=0)
Definition: add_variable_utils_to_python.cpp:104
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
ModelPart::NodesContainerType NodesContainerType
Definition: find_conditions_neighbours_process.h:44
void block_for_each(TIterator itBegin, TIterator itEnd, TFunction &&rFunction)
Execute a functor on all items of a range in parallel.
Definition: parallel_utilities.h:299