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.
assign_parent_element_conditions_process.h
Go to the documentation of this file.
1 // KRATOS ______ __ __ _____ __ __ __
2 // / ____/___ ____ / /_____ ______/ /_/ ___// /________ _______/ /___ ___________ _/ /
3 // / / / __ \/ __ \/ __/ __ `/ ___/ __/\__ \/ __/ ___/ / / / ___/ __/ / / / ___/ __ `/ /
4 // / /___/ /_/ / / / / /_/ /_/ / /__/ /_ ___/ / /_/ / / /_/ / /__/ /_/ /_/ / / / /_/ / /
5 // \____/\____/_/ /_/\__/\__,_/\___/\__//____/\__/_/ \__,_/\___/\__/\__,_/_/ \__,_/_/ MECHANICS
6 //
7 // License: BSD License
8 // license: ContactStructuralMechanicsApplication/license.txt
9 //
10 // Main authors: Vicente Mataix Ferrandiz
11 //
12 
13 #pragma once
14 
15 // System includes
16 
17 // External includes
18 
19 // Project includes
20 #include "processes/process.h"
21 #include "includes/model_part.h"
22 #include "includes/key_hash.h"
23 
24 namespace Kratos
25 {
28 
32 
36 
40 
47 class KRATOS_API(CONTACT_STRUCTURAL_MECHANICS_APPLICATION) AssignParentElementConditionsProcess
48  : public Process
49 {
50 public:
53 
56 
58  using VectorIndexType = std::vector<std::size_t>;
59 
62 
65 
67  using HashMapVectorIntType = std::unordered_map<VectorIndexType, std::size_t, VectorIndexHasherType, VectorIndexComparorType>;
68 
71 
74 
77 
81 
88  ModelPart& rConditionsModelPart,
89  ModelPart& rElementsModelPart
90  )
91  : mrConditionsModelPart(rConditionsModelPart),
92  mrElementsModelPart(rElementsModelPart)
93  {
94  KRATOS_TRY;
95 
96  KRATOS_CATCH("");
97  }
98 
105  Model& rModel,
106  Parameters ThisParameters
107  ) : mrConditionsModelPart(rModel.GetModelPart(ThisParameters["conditions_model_part_name"].GetString())),
108  mrElementsModelPart(rModel.GetModelPart(ThisParameters["elements_model_part_name"].GetString()))
109  {
110  KRATOS_TRY;
111 
112  Parameters default_parameters = GetDefaultParameters();
113  ThisParameters.RecursivelyValidateAndAssignDefaults(default_parameters);
114 
115  mEchoLevel = ThisParameters["echo_level"].GetInt();
116 
117  KRATOS_CATCH("");
118  }
119 
122 
125 
129 
133 
137 
141 
145 
149  void operator()()
150  {
151  Execute();
152  }
153 
157 
161  void Execute() override;
162 
166  void ExecuteInitialize() override;
167 
171  void ExecuteInitializeSolutionStep() override;
172 
176  const Parameters GetDefaultParameters() const override;
177 
181 
185 
189 
191  std::string Info() const override
192  {
193  return "AssignParentElementConditionsProcess";
194  }
195 
197  void PrintInfo(std::ostream& rOStream) const override
198  {
199  rOStream << "AssignParentElementConditionsProcess";
200  }
201 
203  void PrintData(std::ostream& rOStream) const override
204  {
205  }
206 
208 private:
211 
215 
216  ModelPart& mrConditionsModelPart;
217  ModelPart& mrElementsModelPart;
218  int mEchoLevel;
219  HashMapVectorIntType mMapFaceIds;
220 
224 
228 
232 
236 
240 
243 
245  //AssignParentElementConditionsProcess(AssignParentElementConditionsProcess const& rOther);
246 
248 
249 }; // Class AssignParentElementConditionsProcess
250 
252 
255 
259 
261 inline std::istream& operator >> (std::istream& rIStream,
263 
265 inline std::ostream& operator << (std::ostream& rOStream,
267 {
268  rThis.PrintInfo(rOStream);
269  rOStream << std::endl;
270  rThis.PrintData(rOStream);
271 
272  return rOStream;
273 }
274 
275 } // namespace Kratos.
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
void ExecuteInitialize() override
Definition: periodic_interface_process.hpp:37
This process initializes the variables related with the ALM.
Definition: assign_parent_element_conditions_process.h:49
AssignParentElementConditionsProcess(ModelPart &rConditionsModelPart, ModelPart &rElementsModelPart)
Constructor that takes the conditions and elements model parts.
Definition: assign_parent_element_conditions_process.h:87
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: assign_parent_element_conditions_process.h:203
std::unordered_map< VectorIndexType, std::size_t, VectorIndexHasherType, VectorIndexComparorType > HashMapVectorIntType
Define the map considered for face ids.
Definition: assign_parent_element_conditions_process.h:67
void operator()()
Executes the function.
Definition: assign_parent_element_conditions_process.h:149
~AssignParentElementConditionsProcess() override=default
Destructor.
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: assign_parent_element_conditions_process.h:197
ModelPart::NodesContainerType NodesArrayType
Nodes array type definition.
Definition: assign_parent_element_conditions_process.h:73
AssignParentElementConditionsProcess(Model &rModel, Parameters ThisParameters)
Constructs an AssignParentElementConditionsProcess object.
Definition: assign_parent_element_conditions_process.h:104
KRATOS_CLASS_POINTER_DEFINITION(AssignParentElementConditionsProcess)
Pointer definition of AssignParentElementConditionsProcess.
std::string Info() const override
Turn back information as a string.
Definition: assign_parent_element_conditions_process.h:191
std::vector< std::size_t > VectorIndexType
Definition of the vector indexes considered.
Definition: assign_parent_element_conditions_process.h:58
ModelPart::ConditionsContainerType ConditionsArrayType
Conditions array type definition.
Definition: assign_parent_element_conditions_process.h:76
AssignParentElementConditionsProcess(AssignParentElementConditionsProcess const &rOther)=delete
Copy constructor.
Geometry base class.
Definition: geometry.h:71
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
MeshType::ConditionsContainerType ConditionsContainerType
Condintions container. A vector set of Conditions with their Id's as key.
Definition: model_part.h:183
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
int GetInt() const
This method returns the integer contained in the current Parameter.
Definition: kratos_parameters.cpp:666
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
The base class for all processes in Kratos.
Definition: process.h:49
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
This is a key comparer between two vectors of indexes.
Definition: key_hash.h:351
This is a hasher between two vectors of indexes.
Definition: key_hash.h:333