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.
settle_contact_model_structure_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_SETTLE_CONTACT_MODEL_STRUCTURE_PROCESS_H_INCLUDED )
11 #define KRATOS_SETTLE_CONTACT_MODEL_STRUCTURE_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 
63  {
64  public:
67 
70 
74 
77  Flags Options,
78  int EchoLevel = 0)
79  : SettleModelStructureProcess(rMainModelPart, Options, EchoLevel)
80  {
81  }
82 
85  {
86  }
87 
88 
92 
93  void operator()()
94  {
95  Execute();
96  }
97 
98 
102 
103  void Execute() override
104  {
105 
106  };
107 
110  void ExecuteInitialize() override
111  {
112 
113  KRATOS_TRY
114 
115  this->ClearContactConditions();
116 
117  this->ClearContactFlags();
118 
119  //CONDITIONS MASTER_ELEMENTS and MASTER_NODES SEARCH
120  if( mrMainModelPart.GetProcessInfo()[IS_RESTARTED] == true ){
121  std::cout<<" Search conditions masters "<<std::endl;
123  BuildBoundaryProcess.SearchConditionMasters();
124  }
125 
126  //Update Boundary Normals before Contact Search
127  //(needed when meshing of the domains is not performed:
128  // normal directions change with mesh movement)
129  BoundaryNormalsCalculationUtilities BoundaryComputation;
131 
132  //mrMainModelPart.Conditions().Sort();
133  //mrMainModelPart.Conditions().Unique();
134 
135 
136  KRATOS_CATCH(" ")
137  }
138 
139 
142  void ExecuteFinalize() override
143  {
144  KRATOS_TRY
145 
146  // Add contact conditions
147  this->AddContactConditions();
148 
149  // Restore contact nodal flags
150  this->RestoreContactFlags();
151 
152  // Renumerate conditions
153  this->RenumerateConditions();
154 
155  // Clear Contact Forces
156  this->ClearContactForces();
157 
158  // Clear Contact Normals
159  this->ClearContactNormals();
160 
161 
162  KRATOS_CATCH(" ")
163  }
164 
168 
169 
173 
174 
178 
180  std::string Info() const override
181  {
182  return "SettleContactModelStructureProcess";
183  }
184 
186  void PrintInfo(std::ostream& rOStream) const override
187  {
188  rOStream << "SettleContactModelStructureProcess";
189  }
190 
192  void PrintData(std::ostream& rOStream) const override
193  {
194  }
195 
196 
200 
201 
203 
204  protected:
207 
208 
212 
213 
217 
218 
222 
224  {
225 
226  KRATOS_TRY
227 
228  unsigned int consecutive_index = 0;
229  //Renumerate conditions to add in the end of the Elements array (for writing purposes in the ID)
230  //mrMainModelPart.Elements().Sort();
231  int LastElementId = mrMainModelPart.Elements().back().Id();
232  int LastConditionId = mrMainModelPart.Conditions().back().Id();
233  if(LastElementId>LastConditionId){
234  consecutive_index = LastElementId+1;
235  }
236  else{
237  consecutive_index = LastConditionId+1;
238  }
239 
240  for(ModelPart::ConditionsContainerType::iterator it = mrMainModelPart.ConditionsBegin(); it!=mrMainModelPart.ConditionsEnd(); ++it){
241  if(it->Is(CONTACT)){
242  it->SetId(consecutive_index);
243  consecutive_index++;
244  }
245  }
246 
247  KRATOS_CATCH(" ")
248 
249  }
250 
251 
255 
256 
260 
261 
265 
266 
268 
269  private:
272 
276 
280 
284 
285 
286 
287  //**************************************************************************
288  //**************************************************************************
289 
290  void ClearContactForces()
291  {
292  KRATOS_TRY
293 
294 
296 
297  // create contact condition for rigid and deformable bodies
298  for(ModelPart::NodesContainerType::ptr_iterator nd = rNodes.ptr_begin(); nd != rNodes.ptr_end(); ++nd)
299  {
300  if( (*nd)->Is(BOUNDARY) && (*nd)->IsNot(CONTACT) ){
301  array_1d<double, 3> & ContactForce = (*nd)->FastGetSolutionStepValue(CONTACT_FORCE);
302  ContactForce.clear();
303  }
304 
305  }
306 
307  KRATOS_CATCH( "" )
308 
309  }
310 
311  //**************************************************************************
312  //**************************************************************************
313 
314  void ClearContactNormals()
315  {
316  KRATOS_TRY
317 
318 
320 
321  // create contact condition for rigid and deformable bodies
322  for(ModelPart::NodesContainerType::ptr_iterator nd = rNodes.ptr_begin(); nd != rNodes.ptr_end(); ++nd)
323  {
324  if( (*nd)->Is(BOUNDARY) && (*nd)->IsNot(CONTACT) ){
325  array_1d<double, 3> & ContactNormal = (*nd)->FastGetSolutionStepValue(CONTACT_NORMAL);
326  ContactNormal.clear();
327  }
328 
329  }
330 
331  KRATOS_CATCH( "" )
332 
333  }
334 
335  //**************************************************************************
336  //**************************************************************************
337 
338  void ClearContactFlags ( )
339  {
340  KRATOS_TRY
341 
342  for(ModelPart::NodesContainerType::iterator i_node = mrMainModelPart.NodesBegin(); i_node!= mrMainModelPart.NodesEnd(); ++i_node)
343  {
344  if( i_node->Is(CONTACT) ){
345  i_node->Set(CONTACT,false);
346  }
347  }
348 
349  KRATOS_CATCH( "" )
350  }
351 
352 
353  //**************************************************************************
354  //**************************************************************************
355 
356  void RestoreContactFlags ( )
357  {
358 
359  KRATOS_TRY
360 
361  for(ModelPart::ConditionsContainerType::iterator i_cond = mrMainModelPart.ConditionsBegin(); i_cond!= mrMainModelPart.ConditionsEnd(); ++i_cond)
362  {
363  if( i_cond->Is(CONTACT) ){
364  for(unsigned int i=0; i<i_cond->GetGeometry().size(); ++i)
365  {
366  i_cond->GetGeometry()[i].Set(CONTACT,true);
367  }
368  }
369  }
370 
371  KRATOS_CATCH( "" )
372  }
373 
374  //**************************************************************************
375  //**************************************************************************
376 
377  void AddContactConditions()
378  {
379 
380  KRATOS_TRY
381 
382  //Add contact conditions from contact domain
383  std::string ModelPartName;
385  {
386  if(i_mp->Is(CONTACT))
387  ModelPartName = i_mp->Name();
388  }
389 
390  ModelPart& ContactModelPart = mrMainModelPart.GetSubModelPart(ModelPartName);
391 
392  //Add contact conditions to computing domain
394  {
395  if(i_mp->Is(SOLID) && i_mp->Is(ACTIVE))
396  ModelPartName = i_mp->Name();
397  }
398 
399 
400  AddContactConditions(ContactModelPart, mrMainModelPart.GetSubModelPart(ModelPartName));
401 
402  //Add contact conditions to main domain (if added in computing domaing with AddCondition, are automatically added to the main domain, else use push_back)
403  AddContactConditions(ContactModelPart, mrMainModelPart);
404 
405  KRATOS_CATCH( "" )
406 
407  }
408 
409  //**************************************************************************
410  //**************************************************************************
411 
412  void AddContactConditions(ModelPart& rOriginModelPart, ModelPart& rDestinationModelPart)
413  {
414 
415  KRATOS_TRY
416 
417  //*******************************************************************
418  //adding contact conditions
419  //
420 
421  if( mEchoLevel >= 1 ){
422  std::cout<<" ["<<rDestinationModelPart.Name()<<" :: CONDITIONS [OLD:"<<rDestinationModelPart.NumberOfConditions();
423  }
424 
425  for(ModelPart::ConditionsContainerType::iterator ic = rOriginModelPart.ConditionsBegin(); ic!= rOriginModelPart.ConditionsEnd(); ++ic)
426  {
427 
428  if(ic->Is(CONTACT))
429  rDestinationModelPart.AddCondition(*(ic.base()));
430 
431  }
432 
433  if( mEchoLevel >= 1 ){
434  std::cout<<" / NEW:"<<rDestinationModelPart.NumberOfConditions()<<"] "<<std::endl;
435  }
436 
437  KRATOS_CATCH( "" )
438 
439  }
440 
441 
442 
443 
444  //**************************************************************************
445  //**************************************************************************
446 
447  void ClearContactConditions()
448  {
449 
450  KRATOS_TRY
451 
452  //Clear contact conditions from computing domain
453  std::string ModelPartName;
455  {
456  if(i_mp->Is(SOLID) && i_mp->Is(ACTIVE))
457  ModelPartName = i_mp->Name();
458  }
459 
460  ClearContactConditions(mrMainModelPart.GetSubModelPart(ModelPartName));
461 
462  //Clear contact conditions from the main domain
463  ClearContactConditions(mrMainModelPart);
464 
465  KRATOS_CATCH( "" )
466 
467  }
468 
469 
470  //**************************************************************************
471  //**************************************************************************
472 
473  void ClearContactConditions(ModelPart& rModelPart)
474  {
475 
476  KRATOS_TRY
477 
478  //*******************************************************************
479  //clearing contact conditions
480  //
481 
482  if( mEchoLevel >= 1 ){
483  std::cout<<" ["<<rModelPart.Name()<<" :: CONDITIONS [OLD:"<<rModelPart.NumberOfConditions();
484  }
485 
486  ModelPart::ConditionsContainerType PreservedConditions;
487 
488  for(ModelPart::ConditionsContainerType::iterator ic = rModelPart.ConditionsBegin(); ic!= rModelPart.ConditionsEnd(); ++ic)
489  {
490 
491  if(ic->IsNot(CONTACT) && ic->GetGeometry().size() > 1){
492  PreservedConditions.push_back(*(ic.base()));
493  }
494  }
495 
496  rModelPart.Conditions().swap(PreservedConditions);
497 
498  //rModelPart.Conditions().Sort();
499  //rModelPart.Conditions().Unique();
500 
501  if( mEchoLevel >= 1 ){
502  std::cout<<" / NEW:"<<rModelPart.NumberOfConditions()<<"] "<<std::endl;
503  }
504 
505  KRATOS_CATCH( "" )
506 
507  }
508 
512 
513 
517 
518 
522 
525 
527  //SettleContactModelStructureProcess(SettleContactModelStructureProcess const& rOther);
528 
529 
531 
532  }; // Class SettleContactModelStructureProcess
533 
535 
538 
539 
543 
544 
546  inline std::istream& operator >> (std::istream& rIStream,
548 
550  inline std::ostream& operator << (std::ostream& rOStream,
552  {
553  rThis.PrintInfo(rOStream);
554  rOStream << std::endl;
555  rThis.PrintData(rOStream);
556 
557  return rOStream;
558  }
560 
561 
562 } // namespace Kratos.
563 
564 #endif // KRATOS_SETTLE_CONTACT_MODEL_STRUCTURE_PROCESS_H_INCLUDED defined
Definition: boundary_normals_calculation_utilities.hpp:48
void CalculateWeightedBoundaryNormals(ModelPart &rModelPart, int EchoLevel=0)
Calculates the area normal (unitary vector oriented as the normal) and weight the normal to shrink.
Definition: boundary_normals_calculation_utilities.hpp:123
Short class definition.
Definition: build_model_part_boundary_process.hpp:79
bool SearchConditionMasters()
Definition: build_model_part_boundary_process.hpp:240
Definition: flags.h:58
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType::ConditionsContainerType ConditionsContainerType
Condintions container. A vector set of Conditions with their Id's as key.
Definition: model_part.h:183
SubModelPartIterator SubModelPartsEnd()
Definition: model_part.h:1708
ConditionIterator ConditionsBegin(IndexType ThisIndex=0)
Definition: model_part.h:1361
SubModelPartIterator SubModelPartsBegin()
Definition: model_part.h:1698
MeshType::ElementsContainerType ElementsContainerType
Element container. A vector set of Elements with their Id's as key.
Definition: model_part.h:168
std::string & Name()
Definition: model_part.h:1811
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
SubModelPartsContainerType::iterator SubModelPartIterator
Iterator over the sub model parts of this model part.
Definition: model_part.h:258
ConditionsContainerType & Conditions(IndexType ThisIndex=0)
Definition: model_part.h:1381
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
ElementsContainerType & Elements(IndexType ThisIndex=0)
Definition: model_part.h:1189
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
ModelPart & GetSubModelPart(std::string const &SubModelPartName)
Definition: model_part.cpp:2029
NodeIterator NodesEnd(IndexType ThisIndex=0)
Definition: model_part.h:497
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
ConditionIterator ConditionsEnd(IndexType ThisIndex=0)
Definition: model_part.h:1371
void push_back(TPointerType x)
Adds a pointer to the end of the set.
Definition: pointer_vector_set.h:544
Short class definition.
Definition: settle_contact_model_structure_process.hpp:63
void ExecuteInitialize() override
Definition: settle_contact_model_structure_process.hpp:110
virtual ~SettleContactModelStructureProcess()
Destructor.
Definition: settle_contact_model_structure_process.hpp:84
void ExecuteFinalize() override
Definition: settle_contact_model_structure_process.hpp:142
SettleContactModelStructureProcess(ModelPart &rMainModelPart, Flags Options, int EchoLevel=0)
Default constructor.
Definition: settle_contact_model_structure_process.hpp:76
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: settle_contact_model_structure_process.hpp:192
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: settle_contact_model_structure_process.hpp:103
KRATOS_CLASS_POINTER_DEFINITION(SettleContactModelStructureProcess)
Pointer definition of SettleContactModelStructureProcess.
void RenumerateConditions()
Definition: settle_contact_model_structure_process.hpp:223
void operator()()
Definition: settle_contact_model_structure_process.hpp:93
std::string Info() const override
Turn back information as a string.
Definition: settle_contact_model_structure_process.hpp:180
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: settle_contact_model_structure_process.hpp:186
Short class definition.
Definition: settle_model_structure_process.hpp:71
int mEchoLevel
Definition: settle_model_structure_process.hpp:193
ModelPart & mrMainModelPart
Definition: settle_model_structure_process.hpp:189
BOOST_UBLAS_INLINE void clear()
Definition: array_1d.h:325
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
Kratos::ModelPart ModelPart
Definition: kratos_wrapper.h:31
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
integer i
Definition: TensorModule.f:17