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.
static_condensation_utility.h
Go to the documentation of this file.
1 // KRATOS ___| | | |
2 // \___ \ __| __| | | __| __| | | __| _` | |
3 // | | | | | ( | | | | ( | |
4 // _____/ \__|_| \__,_|\___|\__|\__,_|_| \__,_|_| MECHANICS
5 //
6 // License: BSD License
7 // license: StructuralMechanicsApplication/license.txt
8 //
9 // Main authors: Klauss B Sautter
10 // Philipp Bucher
11 //
12 
13 #pragma once
14 
15 
16 // System includes
17 
18 
19 // External includes
20 
21 
22 // Project includes
23 #include "includes/define.h"
24 #include "utilities/math_utils.h"
25 
26 
27 namespace Kratos
28 {
29 
41  {
43  typedef std::size_t SizeType;
44  typedef Matrix MatrixType;
45 
54  ElementType& rTheElement,
55  MatrixType& rLeftHandSideMatrix,
56  const std::vector<int> & rDofList);
57 
65  std::vector<MatrixType> CalculateSchurComplements(
66  ElementType& rTheElement,
67  const MatrixType& rLeftHandSideMatrix,
68  const std::vector<int> & rDofList);
69 
76  std::vector<int> CreateRemainingDofList(
77  ElementType& rTheElement,
78  const std::vector<int> & rDofList);
79 
91  MatrixType& Submatrix,
92  const MatrixType& rLeftHandSideMatrix,
93  const std::vector<int>& rVecA,
94  const std::vector<int>& rVecB,
95  const SizeType& rSizeA,
96  const SizeType& rSizeB); //maybe inline
97 
108  ElementType& rTheElement,
109  Vector& rLocalizedDofVector,
110  Vector& rValues,
111  const std::vector<int>& rDofList,
112  const MatrixType& rLeftHandSideMatrix);
113 
119  SizeType GetNumDofsElement(const ElementType& rTheElement);
120 
121  } // namespace StaticCondensationUtility
122 
123 } // namespace Kratos.
124 
125 
Base class for all Elements.
Definition: element.h:60
void CondenseLeftHandSide(ElementType &rTheElement, MatrixType &rLeftHandSideMatrix, const std::vector< int > &rDofList)
This function is the main operation of this utility. It sorts the reference matrix w....
Definition: static_condensation_utility.cpp:25
std::vector< MatrixType > CalculateSchurComplements(ElementType &rTheElement, const MatrixType &rLeftHandSideMatrix, const std::vector< int > &rDofList)
This function calculates the 4 schur-complements linking the dofs to be condensed to the dofs to rema...
Definition: static_condensation_utility.cpp:66
std::size_t SizeType
Definition: static_condensation_utility.h:43
SizeType GetNumDofsElement(const ElementType &rTheElement)
This function returns the number of dofs of the respective element by using the following input:
Definition: static_condensation_utility.cpp:223
std::vector< int > CreateRemainingDofList(ElementType &rTheElement, const std::vector< int > &rDofList)
This function creates a list containing all dofs to remain by using the following inputs:
Definition: static_condensation_utility.cpp:103
Matrix MatrixType
Definition: static_condensation_utility.h:44
void FillSchurComplements(MatrixType &Submatrix, const MatrixType &rLeftHandSideMatrix, const std::vector< int > &rVecA, const std::vector< int > &rVecB, const SizeType &rSizeA, const SizeType &rSizeB)
This function creates the single schur-complements, called by CalculateSchurComplements,...
Definition: static_condensation_utility.cpp:128
Element ElementType
Definition: static_condensation_utility.h:42
void ConvertingCondensation(ElementType &rTheElement, Vector &rLocalizedDofVector, Vector &rValues, const std::vector< int > &rDofList, const MatrixType &rLeftHandSideMatrix)
This function re-calculates the condensed degree of freedom in relation to the remaining dofs by usin...
Definition: static_condensation_utility.cpp:152
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
This utilitiy condenses given degrees of freedom from any element stiffness matrix to model e....