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.
multiscale_refining_process.h
Go to the documentation of this file.
1 // KRATOS __ __ _____ ____ _ _ ___ _ _ ____
2 // | \/ | ____/ ___|| | | |_ _| \ | |/ ___|
3 // | |\/| | _| \___ \| |_| || || \| | | _
4 // | | | | |___ ___) | _ || || |\ | |_| |
5 // |_| |_|_____|____/|_| |_|___|_| \_|\____| APPLICATION
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Miguel Maso Sotomayor
11 //
12 
13 #if !defined( KRATOS_MULTISCALE_REFINING_PROCESS_H_INCLUDED )
14 #define KRATOS_MULTISCALE_REFINING_PROCESS_H_INCLUDED
15 
16 
17 // System includes
18 #include <string>
19 #include <iostream>
20 #include <algorithm>
21 
22 // External includes
23 
24 // Project includes
25 #include "processes/process.h"
26 #include "includes/model_part.h"
31 
32 
33 namespace Kratos {
34 
37 
41 
45 
49 
53 
85 class KRATOS_API(MESHING_APPLICATION) MultiscaleRefiningProcess : public Process
86 {
87 public:
88 
91 
95  typedef Node NodeType;
96 
100  typedef std::size_t IndexType;
101 
105  typedef std::vector<IndexType> IndexVectorType;
106 
110  typedef std::vector<std::string> StringVectorType;
111 
116 
120  typedef std::unordered_map<IndexType, IndexType> IndexIndexMapType;
121  typedef std::unordered_map<IndexType, std::vector<std::string>> IndexStringMapType;
122  typedef std::unordered_map<IndexType, std::vector<IndexType>> IndexVectorMapType;
123 
128 
132 
134  ModelPart& rThisCoarseModelPart,
135  ModelPart& rThisRefinedModelPart,
136  ModelPart& rThisVisualizationModelPart,
137  Parameters ThisParameters = Parameters(R"({})")
138  );
139 
141  ~MultiscaleRefiningProcess() override = default;
142 
146 
147  void operator()()
148  {
149  Execute();
150  }
151 
155 
156  void Execute() override {}
157 
161  const Parameters GetDefaultParameters() const override;
162 
166  int Check() override;
167 
172  void ExecuteRefinement();
173 
178  void ExecuteCoarsening();
179 
186  static void InitializeVisualizationModelPart(ModelPart& rReferenceModelPart, ModelPart& rNewModelPart);
187 
195  static void InitializeRefinedModelPart(ModelPart& rReferenceModelPart, ModelPart& rNewModelPart);
196 
201  template<class TVarType>
202  void TransferLastStepToCoarseModelPart(const TVarType& rVariable)
203  {
204  ModelPart::NodeIterator refined_begin = mrRefinedModelPart.NodesBegin();
205  for (int i = 0; i < static_cast<int>(mrRefinedModelPart.Nodes().size()); i++)
206  {
207  auto refined_node = refined_begin + i;
208  if (refined_node->GetValue(FATHER_NODES).size() == 1)
209  {
210  auto& value = refined_node->GetValue(FATHER_NODES)[0].FastGetSolutionStepValue(rVariable);
211  value = refined_node->FastGetSolutionStepValue(rVariable); // we only copy the current step
212  }
213  }
214  }
215 
225  template<class TVarType>
226  void TransferSubstepToRefinedInterface(const TVarType& rVariable, const double& rSubstepFraction)
227  {
228  ModelPart::NodeIterator refined_begin = mRefinedInterfaceContainer.begin();
229  for (int i = 0; i < static_cast<int>(mRefinedInterfaceContainer.size()); i++)
230  {
231  auto refined_node = refined_begin + i;
232  GlobalPointersVector<NodeType>& father_nodes = refined_node->GetValue(FATHER_NODES);
233  IndexType number_of_father_nodes = father_nodes.size();
234  std::vector<double> weights = refined_node->GetValue(FATHER_NODES_WEIGHTS);
235 
236  // Transfer the data
237  auto& value = refined_node->FastGetSolutionStepValue(rVariable);
238  value = weights[0] * rSubstepFraction * father_nodes[0].FastGetSolutionStepValue(rVariable);
239  value += weights[0] * (1-rSubstepFraction) * father_nodes[0].FastGetSolutionStepValue(rVariable, 1);
240 
241  if (number_of_father_nodes > 1)
242  {
243  for (IndexType j = 1; j < number_of_father_nodes; j++)
244  {
245  value += weights[j] * rSubstepFraction * father_nodes[j].FastGetSolutionStepValue(rVariable);
246  value += weights[j] * (1-rSubstepFraction) * father_nodes[j].FastGetSolutionStepValue(rVariable, 1);
247  }
248  }
249  }
250  }
251 
259  template<class TVarType>
260  void FixRefinedInterface(const TVarType& rVariable, bool IsFixed)
261  {
262  VariableUtils().ApplyFixity(rVariable, IsFixed, mRefinedInterfaceContainer);
263  }
264 
268 
272 
274  {
275  return mrCoarseModelPart;
276  }
277 
279  {
280  return mrRefinedModelPart;
281  }
282 
284  {
285  return mrVisualizationModelPart;
286  }
287 
291 
293  std::string Info() const override
294  {
295  return "MultiscaleRefiningProcess";
296  }
297 
299  void PrintInfo(std::ostream& rOStream) const override
300  {
301  rOStream << Info();
302  }
303 
305  void PrintData(std::ostream& rOStream) const override
306  {
307  }
308 
312 
314 
315  protected:
316 
319 
323 
327 
331 
335 
339 
343 
345 
346  private:
347 
350 
354 
355  ModelPart& mrCoarseModelPart;
356  ModelPart& mrRefinedModelPart;
357  ModelPart& mrVisualizationModelPart;
358  Parameters mParameters;
359 
360  unsigned int mEchoLevel;
361  int mDivisionsAtSubscale;
362  IndexType mStepDataSize;
363 
364  UniformRefinementUtility mUniformRefinement;
365 
366  NodesArrayType mRefinedInterfaceContainer;
367 
368  std::string mRefinedInterfaceName;
369  std::string mInterfaceConditionName;
370 
371  IndexStringMapType mCollections;
372 
376 
380  void InterpolateLevelBoundaryValuesAtSubStep(const int& rSubStep, const int& rSubSteps);
381 
385  void UpdateSubLevel();
386 
390  void TransferDataToCoarseLevel();
391 
395  void InitializeCoarseModelPartInterface();
396 
400  void InitializeRefinedModelPartInterface();
401 
408  static void InitializeNewModelPart(ModelPart& rReferenceModelPart, ModelPart& rNewModelPart);
409 
416  static void AddAllPropertiesToModelPart(ModelPart& rOriginModelPart, ModelPart& rDestinationModelPart);
417 
424  static void AddAllTablesToModelPart(ModelPart& rOriginModelPart, ModelPart& rDestinationModelPart);
425 
432  void MarkElementsFromNodalFlag();
433 
440  void MarkConditionsFromNodalFlag();
441 
448  void CloneNodesToRefine(IndexType& rNodeId);
449 
454  void CreateElementsToRefine(IndexType& rElemId, IndexIndexMapType& rElemTag);
455 
460  void CreateConditionsToRefine(IndexType& rCondId, IndexIndexMapType& rCondTag);
461 
466  void UpdateVisualizationAfterRefinement();
467 
471  void UpdateVisualizationAfterCoarsening();
472 
481  void IdentifyParentNodesToCoarsen();
482 
490  void IdentifyElementsToErase();
491 
499  void IdentifyConditionsToErase();
500 
509  void IdentifyRefinedNodesToErase();
510 
517  void FinalizeRefinement();
518 
524  void FinalizeCoarsening();
525 
531  void IdentifyCurrentInterface();
532 
537  void UpdateRefinedInterface();
538 
543  void GetLastId(IndexType& rNodesId, IndexType& rElemsId, IndexType& rCondsId);
544 
548 
552 
556 
560 
562  // MultiscaleRefiningProcess& operator=(MultiscaleRefiningProcess const& rOther);
563 
565  // MultiscaleRefiningProcess(MultiscaleRefiningProcess const& rOther);
566 
568 
569 }; // Class MultiscaleRefiningProcess
570 
572 
575 
579 
581 inline std::istream& operator >> (std::istream& rIStream,
583 
585 inline std::ostream& operator << (std::ostream& rOStream,
586  const MultiscaleRefiningProcess& rThis)
587 {
588  rThis.PrintInfo(rOStream);
589  rOStream << std::endl;
590  rThis.PrintData(rOStream);
591 
592  return rOStream;
593 }
595 
596 } // namespace Kratos.
597 
598 #endif // KRATOS_MULTISCALE_REFINING_PROCESS_H_INCLUDED
std::string Info() const override
Turn back information as a string.
Definition: periodic_interface_process.hpp:93
std::size_t IndexType
Definition: flags.h:74
This class is a vector which stores global pointers.
Definition: global_pointers_vector.h:61
size_type size() const
Definition: global_pointers_vector.h:307
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType::NodeIterator NodeIterator
Definition: model_part.h:134
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 a non conforming refinement to perform multi scale analysis @detail This process ...
Definition: multiscale_refining_process.h:86
void TransferLastStepToCoarseModelPart(const TVarType &rVariable)
Copies all the last nodal step data from the refined model part to the coarse one.
Definition: multiscale_refining_process.h:202
~MultiscaleRefiningProcess() override=default
Destructor.
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: multiscale_refining_process.h:156
void operator()()
Definition: multiscale_refining_process.h:147
ModelPart & GetVisualizationModelPart()
Definition: multiscale_refining_process.h:283
std::unordered_map< IndexType, std::vector< std::string > > IndexStringMapType
Definition: multiscale_refining_process.h:121
std::unordered_map< IndexType, IndexType > IndexIndexMapType
Definition: multiscale_refining_process.h:120
std::vector< IndexType > IndexVectorType
Definition: multiscale_refining_process.h:105
void TransferSubstepToRefinedInterface(const TVarType &rVariable, const double &rSubstepFraction)
Copies the nodal step data with a linear interpolation between the last nodal steps of the given vari...
Definition: multiscale_refining_process.h:226
KRATOS_CLASS_POINTER_DEFINITION(MultiscaleRefiningProcess)
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: multiscale_refining_process.h:305
std::size_t IndexType
Definition: multiscale_refining_process.h:100
std::vector< std::string > StringVectorType
Definition: multiscale_refining_process.h:110
ModelPart & GetCoarseModelPart()
Definition: multiscale_refining_process.h:273
std::string Info() const override
Turn back information as a string.
Definition: multiscale_refining_process.h:293
std::unordered_map< IndexType, std::vector< IndexType > > IndexVectorMapType
Definition: multiscale_refining_process.h:122
Node NodeType
Definition: multiscale_refining_process.h:95
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: multiscale_refining_process.h:299
ModelPart & GetRefinedModelPart()
Definition: multiscale_refining_process.h:278
void FixRefinedInterface(const TVarType &rVariable, bool IsFixed)
Applies fixity forthe given variable at the nodes which define the interface.
Definition: multiscale_refining_process.h:260
ModelPart::NodesContainerType NodesArrayType
Definition: multiscale_refining_process.h:115
This class defines the node.
Definition: node.h:65
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
The base class for all processes in Kratos.
Definition: process.h:49
This class splits all the elements until NUMBER_OF_DIVISIONS reaches the desired value @detail A node...
Definition: uniform_refinement_utility.h:83
This class implements a set of auxiliar, already parallelized, methods to perform some common tasks r...
Definition: variable_utils.h:63
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
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
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
int j
Definition: quadrature.py:648
integer i
Definition: TensorModule.f:17