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.
clear_contact_conditions_mesher_process.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosContactMechanicsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: July 2016 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_CLEAR_CONTACT_CONDITIONS_MESHER_PROCESS_H_INCLUDED )
11 #define KRATOS_CLEAR_CONTACT_CONDITIONS_MESHER_PROCESS_H_INCLUDED
12 
13 
14 // System includes
15 
16 // External includes
17 
18 // Project includes
19 
20 #include "includes/model_part.h"
23 
25 //Data:
26 //StepData:
27 //Flags: (checked)
28 // (set)
29 // (modified)
30 // (reset)
31 
32 
33 namespace Kratos
34 {
35 
38 
45 
49 
53 
57 
59 
62  : public MesherProcess
63  {
64  public:
67 
70 
74 
77  int EchoLevel = 0)
78  : mrModelPart(rModelPart)
79  {
80  mEchoLevel = EchoLevel;
81  }
82 
85  {
86  }
87 
88 
92 
93  void operator()()
94  {
95  Execute();
96  }
97 
98 
102 
103  void Execute() override
104  {
105  KRATOS_TRY
106 
107  this->ClearContactConditions();
108 
109  KRATOS_CATCH(" ")
110  };
111 
112 
116 
117 
121 
122 
126 
128  std::string Info() const override
129  {
130  return "ClearContactConditionsMesherProcess";
131  }
132 
134  void PrintInfo(std::ostream& rOStream) const override
135  {
136  rOStream << "ClearContactConditionsMesherProcess";
137  }
138 
140  void PrintData(std::ostream& rOStream) const override
141  {
142  }
143 
144 
148 
149 
151 
152  protected:
155 
156 
160 
161 
165 
166 
170 
171 
175 
176 
180 
181 
185 
186 
188 
189  private:
192 
196 
197  ModelPart& mrModelPart;
198 
199  int mEchoLevel;
200 
204 
208 
209  //**************************************************************************
210  //**************************************************************************
211 
212  void ClearContactConditions()
213  {
214 
215  KRATOS_TRY
216 
217  //Clear contact conditions from contact domain
218  if( mrModelPart.IsNot(CONTACT) )
219  std::cout<<" ModelPart Supplied do not corresponds to the Contact Domain: ("<<mrModelPart.Name()<<")"<<std::endl;
220 
221 
222  ClearContactConditions(mrModelPart);
223 
224 
225  KRATOS_CATCH( "" )
226 
227  }
228 
229  //**************************************************************************
230  //**************************************************************************
231 
232  void ClearContactConditions(ModelPart& rModelPart)
233  {
234 
235  KRATOS_TRY
236 
237  //*******************************************************************
238  //clearing contact conditions
239  //
240 
241  if( mEchoLevel >= 1 ){
242  std::cout<<" ["<<rModelPart.Name()<<" :: CONDITIONS [OLD:"<<rModelPart.NumberOfConditions();
243  }
244 
245  ModelPart::ConditionsContainerType PreservedConditions;
246 
247  for(ModelPart::ConditionsContainerType::iterator ic = rModelPart.ConditionsBegin(); ic!= rModelPart.ConditionsEnd(); ++ic)
248  {
249  if(ic->IsNot(CONTACT)){
250  PreservedConditions.push_back(*(ic.base()));
251  }
252  }
253 
254  rModelPart.Conditions().swap(PreservedConditions);
255 
256  if( mEchoLevel >= 1 ){
257  std::cout<<" / PRE:"<<rModelPart.NumberOfConditions();
258  }
259 
260  //rModelPart.Conditions().Sort();
261  //rModelPart.Conditions().Unique();
262 
263  if( mEchoLevel >= 1 ){
264  std::cout<<" / NEW:"<<rModelPart.NumberOfConditions()<<"] "<<std::endl;
265  }
266 
267  KRATOS_CATCH( "" )
268 
269  }
270 
271 
275 
276 
280 
281 
285 
288 
290  //ClearContactConditionsMesherProcess(ClearContactConditionsMesherProcess const& rOther);
291 
292 
294 
295  }; // Class ClearContactConditionsMesherProcess
296 
298 
301 
302 
306 
307 
309  inline std::istream& operator >> (std::istream& rIStream,
311 
313  inline std::ostream& operator << (std::ostream& rOStream,
315  {
316  rThis.PrintInfo(rOStream);
317  rOStream << std::endl;
318  rThis.PrintData(rOStream);
319 
320  return rOStream;
321  }
323 
324 
325 } // namespace Kratos.
326 
327 #endif // CONTACT_CLEAR_CONTACT_CONDITIONS_MESHER_PROCESS_H_INCLUDED defined
Short class definition.
Definition: clear_contact_conditions_mesher_process.hpp:63
void operator()()
Definition: clear_contact_conditions_mesher_process.hpp:93
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: clear_contact_conditions_mesher_process.hpp:103
KRATOS_CLASS_POINTER_DEFINITION(ClearContactConditionsMesherProcess)
Pointer definition of ClearContactConditionsMesherProcess.
ClearContactConditionsMesherProcess(ModelPart &rModelPart, int EchoLevel=0)
Default constructor.
Definition: clear_contact_conditions_mesher_process.hpp:76
virtual ~ClearContactConditionsMesherProcess()
Destructor.
Definition: clear_contact_conditions_mesher_process.hpp:84
std::string Info() const override
Turn back information as a string.
Definition: clear_contact_conditions_mesher_process.hpp:128
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: clear_contact_conditions_mesher_process.hpp:134
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: clear_contact_conditions_mesher_process.hpp:140
The base class for processes passed to the solution scheme.
Definition: mesher_process.hpp:37
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType::ConditionsContainerType ConditionsContainerType
Condintions container. A vector set of Conditions with their Id's as key.
Definition: model_part.h:183
MeshType::ElementsContainerType ElementsContainerType
Element container. A vector set of Elements with their Id's as key.
Definition: model_part.h:168
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
void push_back(TPointerType x)
Adds a pointer to the end of the set.
Definition: pointer_vector_set.h:544
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
static int EchoLevel
Definition: co_sim_EMPIRE_API.h:42
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
ModelPart::NodesContainerType NodesContainerType
Definition: find_conditions_neighbours_process.h:44
ModelPart::ConditionsContainerType ConditionsContainerType
Definition: find_conditions_neighbours_process.h:45
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
ModelPart::ElementsContainerType ElementsContainerType
Definition: clear_contact_conditions_mesher_process.hpp:43
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432