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.
auxiliar_model_part_utilities.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: Vicente Mataix Ferrandiz
11 //
12 //
13 
14 #pragma once
15 
16 // System includes
17 
18 // External includes
19 
20 // Project includes
21 #include "includes/model_part.h"
23 
24 namespace Kratos
25 {
28 
32 
34  typedef std::size_t IndexType;
35 
39 
47 class KRATOS_API(KRATOS_CORE) AuxiliarModelPartUtilities
48 {
49 public:
52 
55 
57 
61 
66  mrModelPart(rModelPart)
67  {
68  }
69 
70  virtual ~AuxiliarModelPartUtilities()= default;
71 
75 
79 
83 
87 
91 
97  static void CopySubModelPartStructure(const ModelPart& rModelPartToCopyFromIt, ModelPart& rModelPartToCopyIntoIt);
98 
103  void RecursiveEnsureModelPartOwnsProperties(const bool RemovePreviousProperties = true);
104 
109  void EnsureModelPartOwnsProperties(const bool RemovePreviousProperties = true);
110 
120  void RemoveElementAndBelongings(IndexType ElementId, Flags IdentifierFlag = TO_ERASE, IndexType ThisIndex = 0);
121 
131  void RemoveElementAndBelongings(Element& rThisElement, const Flags IdentifierFlag = TO_ERASE, IndexType ThisIndex = 0);
132 
142  void RemoveElementAndBelongings(Element::Pointer pThisElement, const Flags IdentifierFlag = TO_ERASE, IndexType ThisIndex = 0);
143 
153  void RemoveElementAndBelongingsFromAllLevels(IndexType ElementId, const Flags IdentifierFlag = TO_ERASE, IndexType ThisIndex = 0);
154 
164  void RemoveElementAndBelongingsFromAllLevels(Element& rThisElement, const Flags IdentifierFlag = TO_ERASE, IndexType ThisIndex = 0);
165 
175  void RemoveElementAndBelongingsFromAllLevels(Element::Pointer pThisElement, const Flags IdentifierFlag = TO_ERASE, IndexType ThisIndex = 0);
176 
185  void RemoveElementsAndBelongings(Flags IdentifierFlag = TO_ERASE);
186 
195  void RemoveElementsAndBelongingsFromAllLevels(const Flags IdentifierFlag = TO_ERASE);
196 
206  void RemoveConditionAndBelongings(IndexType ConditionId, Flags IdentifierFlag = TO_ERASE, IndexType ThisIndex = 0);
207 
217  void RemoveConditionAndBelongings(Condition& ThisCondition, const Flags IdentifierFlag = TO_ERASE, IndexType ThisIndex = 0);
218 
226  void RemoveConditionAndBelongings(Condition::Pointer pThisCondition, const Flags IdentifierFlag = TO_ERASE, IndexType ThisIndex = 0);
227 
237  void RemoveConditionAndBelongingsFromAllLevels(IndexType ConditionId, const Flags IdentifierFlag = TO_ERASE, IndexType ThisIndex = 0);
238 
248  void RemoveConditionAndBelongingsFromAllLevels(Condition& rThisCondition, const Flags IdentifierFlag = TO_ERASE, IndexType ThisIndex = 0);
249 
259  void RemoveConditionAndBelongingsFromAllLevels(Condition::Pointer pThisCondition, const Flags IdentifierFlag = TO_ERASE, IndexType ThisIndex = 0);
260 
269  void RemoveConditionsAndBelongings(Flags IdentifierFlag = TO_ERASE);
270 
279  void RemoveConditionsAndBelongingsFromAllLevels(const Flags IdentifierFlag = TO_ERASE);
280 
284  void RemoveOrphanNodesFromSubModelParts();
285 
287  template<class TContainerType>
290  const DataLocation DataLoc,
291  TContainerType& data) const
292  {
293  KRATOS_TRY
294 
295  switch (DataLoc)
296  {
297  case (DataLocation::NodeHistorical):{
298  data.resize(mrModelPart.NumberOfNodes());
299 
300  auto inodebegin = mrModelPart.NodesBegin();
301 
302  IndexPartition<IndexType>(mrModelPart.NumberOfNodes()).for_each([&](IndexType Index){
303  auto inode = inodebegin + Index;
304 
305  data[Index] = inode->FastGetSolutionStepValue(rVariable);
306  });
307 
308  break;
309  }
310  case (DataLocation::NodeNonHistorical):{
311  data.resize(mrModelPart.NumberOfNodes());
312 
313  GetScalarDataFromContainer(mrModelPart.Nodes(), rVariable, data);
314  break;
315  }
316  case (DataLocation::Element):{
317  data.resize(mrModelPart.NumberOfElements());
318 
319  GetScalarDataFromContainer(mrModelPart.Elements(), rVariable, data);
320  break;
321  }
322  case (DataLocation::Condition):{
323  data.resize(mrModelPart.NumberOfConditions());
324 
325  GetScalarDataFromContainer(mrModelPart.Conditions(), rVariable, data);
326  break;
327  }
328  case (DataLocation::ModelPart):{
329  data.resize(1);
330  data[0] = mrModelPart[rVariable];
331  break;
332  }
334  data.resize(1);
335  data[0] = mrModelPart.GetProcessInfo()[rVariable];
336  break;
337  }
338  default:{
339  KRATOS_ERROR << "unknown Datalocation" << std::endl;
340  break;
341  }
342  }
343 
344  KRATOS_CATCH("")
345  }
346 
348  template<class TContainerType, class TVarType>
350  const Variable<TVarType>& rVariable,
351  const DataLocation DataLoc,
352  TContainerType& data) const
353  {
354  KRATOS_TRY
355 
356  switch (DataLoc)
357  {
358  case (DataLocation::NodeHistorical):{
359  unsigned int TSize = mrModelPart.NumberOfNodes() > 0 ? mrModelPart.NodesBegin()->FastGetSolutionStepValue(rVariable).size() : 0;
360 
361  TSize = mrModelPart.GetCommunicator().GetDataCommunicator().MaxAll(TSize);
362  data.resize(mrModelPart.NumberOfNodes()*TSize);
363 
364  auto inodebegin = mrModelPart.NodesBegin();
365 
366  IndexPartition<IndexType>(mrModelPart.NumberOfNodes()).for_each([&](IndexType Index){
367  auto inode = inodebegin + Index;
368 
369  const auto& r_val = inode->FastGetSolutionStepValue(rVariable);
370  for(std::size_t dim = 0 ; dim < TSize ; dim++){
371  data[(Index*TSize) + dim] = r_val[dim];
372  }
373  });
374 
375  break;
376  }
377  case (DataLocation::NodeNonHistorical):{
378  unsigned int TSize = mrModelPart.NumberOfNodes() > 0 ? mrModelPart.NodesBegin()->GetValue(rVariable).size() : 0;
379 
380  TSize = mrModelPart.GetCommunicator().GetDataCommunicator().MaxAll(TSize);
381 
382  data.resize(mrModelPart.NumberOfNodes()*TSize);
383 
384  GetVectorDataFromContainer(mrModelPart.Nodes(), TSize, rVariable, data);
385  break;
386  }
387  case (DataLocation::Element):{
388  unsigned int TSize = mrModelPart.NumberOfElements() > 0 ? mrModelPart.ElementsBegin()->GetValue(rVariable).size() : 0;
389 
390  TSize = mrModelPart.GetCommunicator().GetDataCommunicator().MaxAll(TSize);
391 
392  data.resize(mrModelPart.NumberOfElements()*TSize);
393 
394  GetVectorDataFromContainer(mrModelPart.Elements(), TSize, rVariable, data);
395  break;
396  }
397  case (DataLocation::Condition):{
398  unsigned int TSize = mrModelPart.NumberOfConditions() > 0 ? mrModelPart.ConditionsBegin()->GetValue(rVariable).size() : 0;
399 
400  TSize = mrModelPart.GetCommunicator().GetDataCommunicator().MaxAll(TSize);
401 
402  data.resize(mrModelPart.NumberOfConditions()*TSize);
403 
404  GetVectorDataFromContainer(mrModelPart.Conditions(), TSize, rVariable, data);
405  break;
406  }
407  case (DataLocation::ModelPart):{
408  std::size_t TSize = mrModelPart[rVariable].size();
409  data.resize(TSize);
410 
411  IndexType counter = 0;
412  auto& r_val = mrModelPart[rVariable];
413  for(std::size_t dim = 0 ; dim < TSize ; dim++){
414  data[counter++] = r_val[dim];
415  }
416  break;
417  }
419  const std::size_t TSize = mrModelPart.GetProcessInfo()[rVariable].size();
420  data.resize(TSize);
421 
422  IndexType counter = 0;
423  auto& r_val = mrModelPart.GetProcessInfo()[rVariable];
424  for(std::size_t dim = 0 ; dim < TSize ; dim++){
425  data[counter++] = r_val[dim];
426  }
427  break;
428  }
429  default:{
430  KRATOS_ERROR << "unknown Datalocation" << std::endl;
431  break;
432  }
433  }
434 
435  KRATOS_CATCH("")
436  }
437 
439  template<class TContainerType>
442  const DataLocation DataLoc,
443  const TContainerType& rData)
444  {
445  KRATOS_TRY
446 
447  switch (DataLoc)
448  {
449  case (DataLocation::NodeHistorical):{
450  auto inodebegin = mrModelPart.NodesBegin();
451  IndexPartition<IndexType>(mrModelPart.NumberOfNodes()).for_each([&](IndexType Index){
452  auto inode = inodebegin + Index;
453 
454  auto& r_val = inode->FastGetSolutionStepValue(rVariable);
455  r_val = rData[Index];
456  });
457 
458  break;
459  }
460  case (DataLocation::NodeNonHistorical):{
461  SetScalarDataFromContainer(mrModelPart.Nodes(), rVariable, rData);
462  break;
463  }
464  case (DataLocation::Element):{
465  SetScalarDataFromContainer(mrModelPart.Elements(), rVariable, rData);
466  break;
467  }
468  case (DataLocation::Condition):{
469  SetScalarDataFromContainer(mrModelPart.Conditions(), rVariable, rData);
470  break;
471  }
472  case (DataLocation::ModelPart):{
473  mrModelPart[rVariable]= rData[0];
474  break;
475  }
477  mrModelPart.GetProcessInfo()[rVariable] = rData[0] ;
478  break;
479  }
480  default:{
481  KRATOS_ERROR << "unknown Datalocation" << std::endl;
482  break;
483  }
484  }
485 
486  KRATOS_CATCH("")
487  }
488 
490  template<class TContainerType, class TVarType>
492  const Variable<TVarType>& rVariable,
493  const DataLocation DataLoc,
494  const TContainerType& rData)
495  {
496  KRATOS_TRY
497 
498  switch (DataLoc)
499  {
500  case (DataLocation::NodeHistorical):{
501  unsigned int size = mrModelPart.NumberOfNodes() > 0 ? mrModelPart.NodesBegin()->FastGetSolutionStepValue(rVariable).size() : 0;
502 
503  size = mrModelPart.GetCommunicator().GetDataCommunicator().MaxAll(size);
504 
505  auto inodebegin = mrModelPart.NodesBegin();
506  IndexPartition<IndexType>(mrModelPart.NumberOfNodes()).for_each([&](IndexType Index){
507  auto inode = inodebegin + Index;
508  auto& r_val = inode->FastGetSolutionStepValue(rVariable);
509 
510  KRATOS_DEBUG_ERROR_IF(r_val.size() != size) << "mismatch in size!" << std::endl;
511 
512  for(std::size_t dim = 0 ; dim < size ; dim++){
513  r_val[dim] = rData[(Index*size) + dim];
514  }
515  });
516 
517  break;
518  }
519  case (DataLocation::NodeNonHistorical):{
520  unsigned int size = mrModelPart.NumberOfNodes() > 0 ? mrModelPart.NodesBegin()->GetValue(rVariable).size() : 0;
521 
522  size = mrModelPart.GetCommunicator().GetDataCommunicator().MaxAll(size);
523 
524  SetVectorDataFromContainer(mrModelPart.Nodes(), size, rVariable, rData);
525  break;
526  }
527  case (DataLocation::Element):{
528  unsigned int size = mrModelPart.NumberOfElements() > 0 ? mrModelPart.ElementsBegin()->GetValue(rVariable).size() : 0;
529 
530  size = mrModelPart.GetCommunicator().GetDataCommunicator().MaxAll(size);
531 
532  SetVectorDataFromContainer(mrModelPart.Elements(), size, rVariable, rData);
533  break;
534  }
535  case (DataLocation::Condition):{
536  unsigned int size = mrModelPart.NumberOfConditions() > 0 ? mrModelPart.ConditionsBegin()->GetValue(rVariable).size() : 0;
537 
538  size = mrModelPart.GetCommunicator().GetDataCommunicator().MaxAll(size);
539 
540  SetVectorDataFromContainer(mrModelPart.Conditions(), size, rVariable, rData);
541  break;
542  }
543  case (DataLocation::ModelPart):{
544  const std::size_t size = mrModelPart[rVariable].size();
545 
546  IndexType counter = 0;
547  auto& r_val = mrModelPart[rVariable];
548  for(std::size_t dim = 0 ; dim < size ; dim++){
549  r_val[dim] = rData[counter++];
550  }
551  break;
552  }
554  const std::size_t size = mrModelPart.GetProcessInfo()[rVariable].size();
555 
556  IndexType counter = 0;
557  auto& r_val = mrModelPart.GetProcessInfo()[rVariable];
558  for(std::size_t dim = 0 ; dim < size ; dim++){
559  r_val[dim] = rData[counter++];
560  }
561  break;
562  }
563  default:{
564  KRATOS_ERROR << "unknown Datalocation" << std::endl;
565  break;
566  }
567 
568  }
569 
570  KRATOS_CATCH("")
571  }
572 
581  ModelPart& DeepCopyModelPart(
582  const std::string& rNewModelPartName,
583  Model* pModel = nullptr
584  );
585 
596  template<class TClassContainer, class TReferenceClassContainer>
598  ModelPart& rModelPart,
599  TClassContainer& rEntities,
600  TReferenceClassContainer& rReferenceEntities,
601  std::unordered_map<Geometry<Node>::Pointer,Geometry<Node>::Pointer>& rGeometryPointerDatabase
602  )
603  {
604  KRATOS_TRY
605 
606  auto& r_properties= rModelPart.rProperties();
607  rEntities.SetMaxBufferSize(rReferenceEntities.GetMaxBufferSize());
608  rEntities.SetSortedPartSize(rReferenceEntities.GetSortedPartSize());
609  const auto& r_reference_entities_container = rReferenceEntities.GetContainer();
610  auto& r_entities_container = rEntities.GetContainer();
611  const IndexType number_entities = r_reference_entities_container.size();
612  r_entities_container.resize(number_entities);
613  const auto it_ent_begin = r_reference_entities_container.begin();
614  IndexPartition<std::size_t>(number_entities).for_each([&it_ent_begin,&r_entities_container,&rGeometryPointerDatabase,&r_properties](std::size_t i) {
615  auto it_ent = it_ent_begin + i;
616  auto& p_old_ent = (*it_ent);
617  auto p_new_ent = p_old_ent->Create(p_old_ent->Id(), rGeometryPointerDatabase[p_old_ent->pGetGeometry()], r_properties(p_old_ent->pGetProperties()->Id()));
618  p_new_ent->SetData(p_old_ent->GetData());
619  p_new_ent->Set(Flags(*p_old_ent));
620  r_entities_container[i] = p_new_ent;
621  });
622 
623  KRATOS_CATCH("")
624  }
625 
627  virtual std::string Info() const
628  {
629  return "AuxiliarModelPartUtilities";
630  }
631 
633  virtual void PrintInfo(std::ostream& rOStream) const
634  {
635  rOStream << Info() << std::endl;
636  }
637 
639  virtual void PrintData(std::ostream& rOStream) const
640  {
641  rOStream << Info() << std::endl;
642  }
643 
644 private:
650 
651  ModelPart& mrModelPart;
652 
656 
660 
661  template<typename TDataType, class TContainerType, class TDataContainerType>
662  void GetScalarDataFromContainer(
663  const TContainerType& rContainer,
664  const Variable<TDataType>& rVariable,
665  TDataContainerType& data) const
666  {
667  KRATOS_TRY
668 
669  DataSizeCheck(rContainer.size(), data.size());
670 
671  IndexPartition<std::size_t>(rContainer.size()).for_each([&](std::size_t index){
672  const auto& r_entity = *(rContainer.begin() + index);
673  data[index] = r_entity.GetValue(rVariable);
674  });
675 
676  KRATOS_CATCH("")
677  }
678 
679  template<typename TDataType, class TContainerType, class TDataContainerType>
680  void GetVectorDataFromContainer(
681  const TContainerType& rContainer,
682  const std::size_t VectorSize,
683  const Variable<TDataType>& rVariable,
684  TDataContainerType& data) const
685  {
686  KRATOS_TRY
687 
688  DataSizeCheck(rContainer.size()*VectorSize, data.size());
689 
690  IndexPartition<std::size_t>(rContainer.size()).for_each([&](std::size_t index){
691  const auto& r_entity = *(rContainer.begin() + index);
692  const auto& r_val = r_entity.GetValue(rVariable);
693  for(std::size_t dim = 0 ; dim < VectorSize ; dim++){
694  data[(VectorSize*index) + dim] = r_val[dim];
695  }
696  });
697 
698  KRATOS_CATCH("")
699  }
700 
701  template<typename TDataType, class TContainerType, class TDataContainerType>
702  void SetScalarDataFromContainer(
703  TContainerType& rContainer,
704  const Variable<TDataType>& rVariable,
705  const TDataContainerType& rData) const
706  {
707  KRATOS_TRY
708 
709  DataSizeCheck(rContainer.size(), rData.size());
710 
711  IndexPartition<std::size_t>(rContainer.size()).for_each([&](std::size_t index){
712  auto& r_entity = *(rContainer.begin() + index);
713  r_entity.SetValue(rVariable, rData[index]);
714  });
715 
716  KRATOS_CATCH("")
717  }
718 
719  template<typename TDataType, class TContainerType, class TDataContainerType>
720  void SetVectorDataFromContainer(
721  TContainerType& rContainer,
722  const std::size_t VectorSize,
723  const Variable<TDataType>& rVariable,
724  const TDataContainerType& rData) const
725  {
726  KRATOS_TRY
727 
728  DataSizeCheck(rContainer.size()*VectorSize, rData.size());
729 
730  IndexPartition<std::size_t>(rContainer.size()).for_each([&](std::size_t index){
731  auto& r_entity = *(rContainer.begin() + index);
732  TDataType aux;
733  KRATOS_DEBUG_ERROR_IF(aux.size() != VectorSize) << "mismatch in size!" << std::endl;
734  for(std::size_t dim = 0 ; dim < VectorSize ; dim++){
735  aux[dim] = rData[(VectorSize*index) + dim];
736  }
737  r_entity.SetValue(rVariable, aux);
738  });
739 
740  KRATOS_CATCH("")
741  }
742 
743  void DataSizeCheck(
744  const std::size_t ContainerSize,
745  const std::size_t DataSize) const
746  {
747  KRATOS_ERROR_IF(ContainerSize != DataSize) << "Mismatch in size! Container size: " << ContainerSize << " | Data size: " << DataSize << std::endl;
748  }
749 
756  void DeepCopySubModelPart(
757  const ModelPart& rOldModelPart,
758  ModelPart& rNewModelPart
759  );
760 
765 
767 
771 
775 };// class AuxiliarModelPartUtilities
776 
780 
781 
785 
787 
788 } // namespace Kratos.
std::string Info() const override
Turn back information as a string.
Definition: periodic_interface_process.hpp:93
This utility includes auxiliar methods not included in the model part to avoid increase more than nec...
Definition: auxiliar_model_part_utilities.h:48
void DeepCopyEntities(ModelPart &rModelPart, TClassContainer &rEntities, TReferenceClassContainer &rReferenceEntities, std::unordered_map< Geometry< Node >::Pointer, Geometry< Node >::Pointer > &rGeometryPointerDatabase)
This method deep copies a entities.
Definition: auxiliar_model_part_utilities.h:597
void GetVectorData(const Variable< TVarType > &rVariable, const DataLocation DataLoc, TContainerType &data) const
To Export a Vector data (std::vector/array/..)
Definition: auxiliar_model_part_utilities.h:349
virtual void PrintData(std::ostream &rOStream) const
Print object's data.
Definition: auxiliar_model_part_utilities.h:639
void SetScalarData(const Variable< typename TContainerType::value_type > &rVariable, const DataLocation DataLoc, const TContainerType &rData)
To Import a Scalar data (Double/int/...)
Definition: auxiliar_model_part_utilities.h:440
void SetVectorData(const Variable< TVarType > &rVariable, const DataLocation DataLoc, const TContainerType &rData)
To Import a Vector data (std::vector/array/..)
Definition: auxiliar_model_part_utilities.h:491
AuxiliarModelPartUtilities(ModelPart &rModelPart)
Definition: auxiliar_model_part_utilities.h:65
KRATOS_CLASS_POINTER_DEFINITION(AuxiliarModelPartUtilities)
Counted pointer of AuxiliarModelPartUtilities.
virtual ~AuxiliarModelPartUtilities()=default
virtual void PrintInfo(std::ostream &rOStream) const
Print information about this object.
Definition: auxiliar_model_part_utilities.h:633
void GetScalarData(const Variable< typename TContainerType::value_type > &rVariable, const DataLocation DataLoc, TContainerType &data) const
To Export a Scalar data (Double/int/...)
Definition: auxiliar_model_part_utilities.h:288
virtual std::string Info() const
Turn back information as a string.
Definition: auxiliar_model_part_utilities.h:627
Base class for all Conditions.
Definition: condition.h:59
Base class for all Elements.
Definition: element.h:60
Definition: flags.h:58
Geometry base class.
Definition: geometry.h:71
This class is useful for index iteration over containers.
Definition: parallel_utilities.h:451
void for_each(TUnaryFunction &&f)
Definition: parallel_utilities.h:514
This class aims to manage different model parts across multi-physics simulations.
Definition: model.h:60
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
PropertiesContainerType & rProperties(IndexType ThisIndex=0)
Definition: model_part.h:1003
void SetMaxBufferSize(const size_type NewSize)
Set the maximum size of buffer used in the container.
Definition: pointer_vector_set.h:802
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
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(conditional)
Definition: exception.h:162
#define KRATOS_DEBUG_ERROR_IF(conditional)
Definition: exception.h:171
Kratos::ModelPart ModelPart
Definition: kratos_wrapper.h:31
DataLocation
Enum for location of data.
Definition: global_variables.h:48
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
ProcessInfo
Definition: edgebased_PureConvection.py:116
def Index()
Definition: hdf5_io_tools.py:38
data
Definition: mesh_to_mdpa_converter.py:59
int counter
Definition: script_THERMAL_CORRECT.py:218
int dim
Definition: sensitivityMatrix.py:25
integer i
Definition: TensorModule.f:17