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.
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Kratos::TrilinosRefineMesh Class Reference

#include <trilinos_refine_mesh.h>

Collaboration diagram for Kratos::TrilinosRefineMesh:

Public Types

typedef ModelPart::NodesContainerType NodesArrayType
 
typedef ModelPart::ElementsContainerType ElementsArrayType
 
typedef ModelPart::ConditionsContainerType ConditionsArrayType
 
typedef vector< MatrixMatrix_Order_Tensor
 
typedef vector< VectorVector_Order_Tensor
 
typedef vector< Vector_Order_TensorNode_Vector_Order_Tensor
 
typedef Node PointType
 
typedef Node ::Pointer PointPointerType
 
typedef std::vector< PointType::Pointer > PointVector
 
typedef PointVector::iterator PointIterator
 

Public Member Functions

 TrilinosRefineMesh (ModelPart &model_part, Epetra_MpiComm &Comm)
 
 ~TrilinosRefineMesh ()
 
void Local_Refine_Mesh (bool refine_on_reference, bool interpolate_internal_variables, int domain_size)
 
void PrintDebugInfo ()
 
void Clear ()
 This function frees all of the memory used. More...
 
void CSR_Row_Matrix (ModelPart &this_model_part, Kratos::shared_ptr< Epetra_FECrsMatrix > &p_edge_ids)
 
void Search_Edge_To_Be_Refined (ModelPart &this_model_part, Kratos::shared_ptr< Epetra_FECrsMatrix > &p_edge_ids, Kratos::shared_ptr< Epetra_FECrsMatrix > &p_partition_ids)
 
void Create_List_Of_New_Nodes (ModelPart &this_model_part, Kratos::shared_ptr< Epetra_FECrsMatrix > &p_edge_ids, Kratos::shared_ptr< Epetra_FECrsMatrix > &p_partition_ids, vector< int > &List_New_Nodes, vector< int > &partition_new_nodes, vector< array_1d< int, 2 > > &father_node_ids)
 
void Calculate_Coordinate_And_Insert_New_Nodes (ModelPart &this_model_part, const vector< array_1d< int, 2 > > &father_node_ids, const vector< int > &List_New_Nodes, const vector< int > &partition_new_nodes)
 
void Erase_Old_Element_And_Create_New_Triangle_Element (ModelPart &this_model_part, const Kratos::shared_ptr< Epetra_FECrsMatrix > p_edge_ids, PointerVector< Element > &New_Elements, bool interpolate_internal_variables)
 
void Calculate_Triangle_Edges (Element::GeometryType &geom, const Kratos::shared_ptr< Epetra_FECrsMatrix > p_edge_ids, int *edge_ids, int aux[6])
 
void Renumbering_Elements (ModelPart &this_model_part, PointerVector< Element > &New_Elements)
 
void Renumbering_Conditions (ModelPart &this_model_part, PointerVector< Condition > &New_Conditions)
 

Protected Member Functions

void InterpolateInternalVariables (const int &nel, const Element::Pointer father_elem, Element::Pointer child_elem, const ProcessInfo &rCurrentProcessInfo)
 
double GetValueFromRow (int row, int j, int row_size, int *indices, double *values)
 
double GetUpperTriangularMatrixValue (const Kratos::shared_ptr< Epetra_FECrsMatrix > &p_edge_ids, int index_0, int index_1, int &MaxNumEntries, int &NumEntries, int *Indices, double *values)
 
void Erase_Old_Element_And_Create_New_Tetra_Element (ModelPart &this_model_part, const Kratos::shared_ptr< Epetra_FECrsMatrix > p_edge_ids, PointerVector< Element > &New_Elements, bool interpolate_internal_variables)
 
unsigned int CreateCenterNode (Geometry< Node > &geom, ModelPart &model_part, ModelPart::NodesContainerType &center_nodes)
 
void AssignIdToCenterNode (ModelPart::NodesContainerType &center_nodes)
 
void Erase_Old_Condition_And_Create_New_Triangle_Conditions (ModelPart &this_model_part, const Kratos::shared_ptr< Epetra_FECrsMatrix > p_edge_ids, PointerVector< Condition > &New_Conditions, bool interpolate_internal_variables)
 
Node::Pointer AuxCreateNewNode (ModelPart &r_model_part, int Id, double x, double y, double z)
 

Protected Attributes

ModelPartmr_model_part
 
Epetra_MpiComm & mrComm
 
Kratos::shared_ptr< Epetra_Map > mp_overlapping_map
 
Kratos::shared_ptr< Epetra_Map > mp_non_overlapping_map
 
int mtotal_number_of_existing_nodes
 
Kratos::shared_ptr< Epetra_FECrsGraph > mp_non_overlapping_graph
 
Kratos::shared_ptr< Epetra_CrsGraph > mp_overlapping_graph
 

Detailed Description

This Function is designed to refine a mesh of triangles or tetrahedra in MPI This is achieved by using the trilinos epetra facilities. Please note that Trilinos has to be patched to work correctly, meaning that versions > 10.6 are needed to compile this

Member Typedef Documentation

◆ ConditionsArrayType

◆ ElementsArrayType

◆ Matrix_Order_Tensor

◆ Node_Vector_Order_Tensor

◆ NodesArrayType

◆ PointIterator

typedef PointVector::iterator Kratos::TrilinosRefineMesh::PointIterator

◆ PointPointerType

◆ PointType

◆ PointVector

typedef std::vector<PointType::Pointer> Kratos::TrilinosRefineMesh::PointVector

◆ Vector_Order_Tensor

Constructor & Destructor Documentation

◆ TrilinosRefineMesh()

Kratos::TrilinosRefineMesh::TrilinosRefineMesh ( ModelPart model_part,
Epetra_MpiComm &  Comm 
)
inline

constructor:

Parameters
ModelPart&the model part to be refined
Epetra_MpiCommthe Epetra Communicator to be used

◆ ~TrilinosRefineMesh()

Kratos::TrilinosRefineMesh::~TrilinosRefineMesh ( )
inline

Member Function Documentation

◆ AssignIdToCenterNode()

void Kratos::TrilinosRefineMesh::AssignIdToCenterNode ( ModelPart::NodesContainerType center_nodes)
inlineprotected

◆ AuxCreateNewNode()

Node::Pointer Kratos::TrilinosRefineMesh::AuxCreateNewNode ( ModelPart r_model_part,
int  Id,
double  x,
double  y,
double  z 
)
inlineprotected

◆ Calculate_Coordinate_And_Insert_New_Nodes()

void Kratos::TrilinosRefineMesh::Calculate_Coordinate_And_Insert_New_Nodes ( ModelPart this_model_part,
const vector< array_1d< int, 2 > > &  father_node_ids,
const vector< int > &  List_New_Nodes,
const vector< int > &  partition_new_nodes 
)
inline

calculating the coordinate of the news nodes

inserting the news node in the model part

  • intepolating the data

*copying this data in the position of the vector we are interested in

WARNING = only for reactions;

◆ Calculate_Triangle_Edges()

void Kratos::TrilinosRefineMesh::Calculate_Triangle_Edges ( Element::GeometryType geom,
const Kratos::shared_ptr< Epetra_FECrsMatrix >  p_edge_ids,
int edge_ids,
int  aux[6] 
)
inline


◆ Clear()

void Kratos::TrilinosRefineMesh::Clear ( )
inline

This function frees all of the memory used.

◆ Create_List_Of_New_Nodes()

void Kratos::TrilinosRefineMesh::Create_List_Of_New_Nodes ( ModelPart this_model_part,
Kratos::shared_ptr< Epetra_FECrsMatrix > &  p_edge_ids,
Kratos::shared_ptr< Epetra_FECrsMatrix > &  p_partition_ids,
vector< int > &  List_New_Nodes,
vector< int > &  partition_new_nodes,
vector< array_1d< int, 2 > > &  father_node_ids 
)
inline
  • New Id de los Nodos

◆ CreateCenterNode()

unsigned int Kratos::TrilinosRefineMesh::CreateCenterNode ( Geometry< Node > &  geom,
ModelPart model_part,
ModelPart::NodesContainerType center_nodes 
)
inlineprotected
  • intepolating the data

*copying this data in the position of the vector we are interested in

◆ CSR_Row_Matrix()

void Kratos::TrilinosRefineMesh::CSR_Row_Matrix ( ModelPart this_model_part,
Kratos::shared_ptr< Epetra_FECrsMatrix > &  p_edge_ids 
)
inline

◆ Erase_Old_Condition_And_Create_New_Triangle_Conditions()

void Kratos::TrilinosRefineMesh::Erase_Old_Condition_And_Create_New_Triangle_Conditions ( ModelPart this_model_part,
const Kratos::shared_ptr< Epetra_FECrsMatrix >  p_edge_ids,
PointerVector< Condition > &  New_Conditions,
bool  interpolate_internal_variables 
)
inlineprotected


as the name says...

  • crea las nuevas conectividades
  • crea los nuevos Conditionos

setting the internal variables in the child elem

  • all of the Conditions to be erased are at the end

*now remove all of the "old" Conditions

  • adding news Conditions to the model part

◆ Erase_Old_Element_And_Create_New_Tetra_Element()

void Kratos::TrilinosRefineMesh::Erase_Old_Element_And_Create_New_Tetra_Element ( ModelPart this_model_part,
const Kratos::shared_ptr< Epetra_FECrsMatrix >  p_edge_ids,
PointerVector< Element > &  New_Elements,
bool  interpolate_internal_variables 
)
inlineprotected

setting the internal variables in the child elem

  • all of the elements to be erased are at the end

*now remove all of the "old" elements

  • adding news elements to the model part

◆ Erase_Old_Element_And_Create_New_Triangle_Element()

void Kratos::TrilinosRefineMesh::Erase_Old_Element_And_Create_New_Triangle_Element ( ModelPart this_model_part,
const Kratos::shared_ptr< Epetra_FECrsMatrix >  p_edge_ids,
PointerVector< Element > &  New_Elements,
bool  interpolate_internal_variables 
)
inline


as the name says...

  • crea las nuevas conectividades
  • crea los nuevos elementos

setting the internal variables in the child elem

  • adding news elements to the model part
  • all of the elements to be erased are at the end

*now remove all of the "old" elements

◆ GetUpperTriangularMatrixValue()

double Kratos::TrilinosRefineMesh::GetUpperTriangularMatrixValue ( const Kratos::shared_ptr< Epetra_FECrsMatrix > &  p_edge_ids,
int  index_0,
int  index_1,
int MaxNumEntries,
int NumEntries,
int Indices,
double values 
)
inlineprotected

◆ GetValueFromRow()

double Kratos::TrilinosRefineMesh::GetValueFromRow ( int  row,
int  j,
int  row_size,
int indices,
double values 
)
inlineprotected

◆ InterpolateInternalVariables()

void Kratos::TrilinosRefineMesh::InterpolateInternalVariables ( const int nel,
const Element::Pointer  father_elem,
Element::Pointer  child_elem,
const ProcessInfo rCurrentProcessInfo 
)
inlineprotected

this function transfers the Constitutive Law internal variables from the father to the child. note that this is done through the vector Variable INTERNAL_VARIABLES which should also contain the geometric data needed for this.

◆ Local_Refine_Mesh()

void Kratos::TrilinosRefineMesh::Local_Refine_Mesh ( bool  refine_on_reference,
bool  interpolate_internal_variables,
int  domain_size 
)
inline

Local_Refine_Mesh is the umbrella function that performs the refinement. The elements to be refined are identified by the flag SPLIT_ELEMENT=true current level of refinement is stored element by element in the flag REFINEMENT_LEVEL of type integer NOTE: the refinement is performed with the ONLY AIM of delivering a conformant mesh, that is, no step is performed to ensure the quality of the resultign mesh. The user is expected to select the elements to be refined so to guarantee that some quality level is retained. No effort is made to keep the resulting mesh balanced, consequently a load balancing step should be performed after the refinement step.

Parameters
refine_on_reference--> it controls if the interpolation of data (both internal and nodal) should be performed on the initial domain (if set it to true) OR on the deformed one. this is very important for example for refining correctly total lagrangian elements. If the parameter is set to two, the model_part is expected to contain the DISPLACEMENT variable
interpolate_internal_variables--> this flag controls if the internal variables (Constitutive Law vars) should be interpolated or not
  • the news nodes
  • the news nodes
  • edges where are the news nodes
  • the coordinate of the new nodes

◆ PrintDebugInfo()

void Kratos::TrilinosRefineMesh::PrintDebugInfo ( )
inline

◆ Renumbering_Conditions()

void Kratos::TrilinosRefineMesh::Renumbering_Conditions ( ModelPart this_model_part,
PointerVector< Condition > &  New_Conditions 
)
inline


function to ensure that the Conditions and conditions are numbered consecutively

◆ Renumbering_Elements()

void Kratos::TrilinosRefineMesh::Renumbering_Elements ( ModelPart this_model_part,
PointerVector< Element > &  New_Elements 
)
inline


function to ensure that the elements and conditions are numbered consecutively

◆ Search_Edge_To_Be_Refined()

void Kratos::TrilinosRefineMesh::Search_Edge_To_Be_Refined ( ModelPart this_model_part,
Kratos::shared_ptr< Epetra_FECrsMatrix > &  p_edge_ids,
Kratos::shared_ptr< Epetra_FECrsMatrix > &  p_partition_ids 
)
inline

Member Data Documentation

◆ mp_non_overlapping_graph

Kratos::shared_ptr<Epetra_FECrsGraph> Kratos::TrilinosRefineMesh::mp_non_overlapping_graph
protected

◆ mp_non_overlapping_map

Kratos::shared_ptr<Epetra_Map> Kratos::TrilinosRefineMesh::mp_non_overlapping_map
protected

◆ mp_overlapping_graph

Kratos::shared_ptr<Epetra_CrsGraph> Kratos::TrilinosRefineMesh::mp_overlapping_graph
protected

◆ mp_overlapping_map

Kratos::shared_ptr<Epetra_Map> Kratos::TrilinosRefineMesh::mp_overlapping_map
protected

◆ mr_model_part

ModelPart& Kratos::TrilinosRefineMesh::mr_model_part
protected

◆ mrComm

Epetra_MpiComm& Kratos::TrilinosRefineMesh::mrComm
protected

◆ mtotal_number_of_existing_nodes

int Kratos::TrilinosRefineMesh::mtotal_number_of_existing_nodes
protected

The documentation for this class was generated from the following file: