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.
eigensolver_scheme.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosSolidMechanicsApplication $
3 // Created by: $Author: michael.andre@tum.de $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: September 2016 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 
11 #if !defined(KRATOS_EIGENSOLVER_SCHEME_H_INCLUDED)
12 #define KRATOS_EIGENSOLVER_SCHEME_H_INCLUDED
13 
14 
15 // System includes
16 
17 // External includes
18 
19 // Project includes
20 #include "includes/define.h"
21 #include "includes/element.h"
22 #include "includes/condition.h"
23 #include "includes/process_info.h"
26 
27 // Application includes
29 
30 namespace Kratos
31 {
32 
35 
39 
41 
44 
48 
52 
54 template<class TSparseSpace,
55  class TDenseSpace
56  >
57 class EigensolverScheme : public SolutionScheme<TSparseSpace,TDenseSpace>
58 {
59  public:
62 
64 
67 
70 
74 
77  :BaseType()
78  {
79  }
80 
83  :BaseType(rOptions)
84  {
85  }
86 
88  ~EigensolverScheme() override {}
89 
93 
97 
99  Element::Pointer pCurrentElement,
100  LocalSystemMatrixType& rLHS_Contribution,
101  LocalSystemVectorType& rRHS_Contribution,
102  Element::EquationIdVectorType& rEquationId,
103  ProcessInfo& rCurrentProcessInfo
104  ) override
105  {
106  KRATOS_TRY
107 
108  if (rCurrentProcessInfo[BUILD_LEVEL] == 1)
109  { // mass matrix
110  pCurrentElement->CalculateMassMatrix(rLHS_Contribution,rCurrentProcessInfo);
111  std::size_t LocalSize = rLHS_Contribution.size1();
112  if (rRHS_Contribution.size() != LocalSize)
113  rRHS_Contribution.resize(LocalSize,false);
114  noalias(rRHS_Contribution) = ZeroVector(LocalSize);
115  }
116  else if (rCurrentProcessInfo[BUILD_LEVEL] == 2) // stiffness matrix
117  {
118  pCurrentElement->CalculateLocalSystem(rLHS_Contribution,rRHS_Contribution,rCurrentProcessInfo);
119  }
120  else
121  {
122  KRATOS_ERROR <<"Invalid BUILD_LEVEL" << std::endl;
123  }
124 
125  pCurrentElement->EquationIdVector(rEquationId,rCurrentProcessInfo);
126 
127  KRATOS_CATCH("")
128  }
129 
131  Element::Pointer pCurrentElement,
132  LocalSystemMatrixType& rLHS_Contribution,
133  Element::EquationIdVectorType& rEquationId,
134  ProcessInfo& rCurrentProcessInfo) override
135  {
136  KRATOS_TRY
137 
138  LocalSystemVectorType RHS_Contribution;
139  RHS_Contribution.resize(rLHS_Contribution.size1(), false);
141  pCurrentElement,
142  rLHS_Contribution,
143  RHS_Contribution,
144  rEquationId,
145  rCurrentProcessInfo);
146 
147  KRATOS_CATCH("")
148  }
149 
151  Condition::Pointer pCurrentCondition,
152  LocalSystemMatrixType& rLHS_Contribution,
153  LocalSystemVectorType& rRHS_Contribution,
154  Condition::EquationIdVectorType& rEquationId,
155  ProcessInfo& rCurrentProcessInfo) override
156  {
157  KRATOS_TRY
158 
159  if (rCurrentProcessInfo[BUILD_LEVEL] == 1)
160  { // mass matrix
161  pCurrentCondition->CalculateMassMatrix(rLHS_Contribution,rCurrentProcessInfo);
162  std::size_t LocalSize = rLHS_Contribution.size1();
163  if (rRHS_Contribution.size() != LocalSize)
164  {
165  rRHS_Contribution.resize(LocalSize,false);
166  }
167  noalias(rRHS_Contribution) = ZeroVector(LocalSize);
168  }
169  else if (rCurrentProcessInfo[BUILD_LEVEL] == 2) // stiffness matrix
170  {
171  pCurrentCondition->CalculateLocalSystem(rLHS_Contribution,rRHS_Contribution,rCurrentProcessInfo);
172  }
173  else
174  {
175  KRATOS_ERROR <<"Invalid BUILD_LEVEL" << std::endl;
176  }
177 
178  pCurrentCondition->EquationIdVector(rEquationId,rCurrentProcessInfo);
179 
180  KRATOS_CATCH("")
181  }
182 
184  Condition::Pointer pCurrentCondition,
185  LocalSystemMatrixType& rLHS_Contribution,
186  Condition::EquationIdVectorType& rEquationId,
187  ProcessInfo& rCurrentProcessInfo) override
188  {
189  KRATOS_TRY
190 
191  LocalSystemVectorType RHS_Contribution;
192  RHS_Contribution.resize(rLHS_Contribution.size1(), false);
194  pCurrentCondition,
195  rLHS_Contribution,
196  RHS_Contribution,
197  rEquationId,
198  rCurrentProcessInfo);
199 
200  KRATOS_CATCH("")
201  }
202 
206 
210 
214 
216 
217  protected:
220 
224 
228 
232 
236 
240 
244 
246 
247  private:
250 
254 
258 
262 
266 
270 
274 
276 
277 }; // Class Eigen Solver Scheme
278 
280 
283 
285 
286 } // namespace Kratos.
287 
288 #endif // KRATOS_EIGENSOLVER_SCHEME_H_INCLUDED defined
289 
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
An adapter scheme for obtaining mass and stiffness matrices for dynamic eigenvalue problems.
Definition: eigensolver_scheme.hpp:58
void Condition_CalculateSystemContributions(Condition::Pointer pCurrentCondition, LocalSystemMatrixType &rLHS_Contribution, LocalSystemVectorType &rRHS_Contribution, Condition::EquationIdVectorType &rEquationId, ProcessInfo &rCurrentProcessInfo) override
Definition: eigensolver_scheme.hpp:150
EigensolverScheme()
Default Constructor.
Definition: eigensolver_scheme.hpp:76
void Calculate_LHS_Contribution(Element::Pointer pCurrentElement, LocalSystemMatrixType &rLHS_Contribution, Element::EquationIdVectorType &rEquationId, ProcessInfo &rCurrentProcessInfo) override
Definition: eigensolver_scheme.hpp:130
void Condition_Calculate_LHS_Contribution(Condition::Pointer pCurrentCondition, LocalSystemMatrixType &rLHS_Contribution, Condition::EquationIdVectorType &rEquationId, ProcessInfo &rCurrentProcessInfo) override
Definition: eigensolver_scheme.hpp:183
BaseType::LocalSystemVectorType LocalSystemVectorType
Definition: eigensolver_scheme.hpp:68
BaseType::LocalSystemMatrixType LocalSystemMatrixType
Definition: eigensolver_scheme.hpp:69
SolutionScheme< TSparseSpace, TDenseSpace > BaseType
Definition: eigensolver_scheme.hpp:65
EigensolverScheme(Flags &rOptions)
Constructor.
Definition: eigensolver_scheme.hpp:82
KRATOS_CLASS_POINTER_DEFINITION(EigensolverScheme)
void CalculateSystemContributions(Element::Pointer pCurrentElement, LocalSystemMatrixType &rLHS_Contribution, LocalSystemVectorType &rRHS_Contribution, Element::EquationIdVectorType &rEquationId, ProcessInfo &rCurrentProcessInfo) override
Definition: eigensolver_scheme.hpp:98
~EigensolverScheme() override
Destructor.
Definition: eigensolver_scheme.hpp:88
BaseType::SolutionSchemePointerType BasePointerType
Definition: eigensolver_scheme.hpp:66
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
Definition: flags.h:58
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Solution scheme base class.
Definition: solution_scheme.hpp:54
SolutionSchemeType::Pointer SolutionSchemePointerType
Definition: solution_scheme.hpp:60
TDenseSpace::MatrixType LocalSystemMatrixType
Definition: solution_scheme.hpp:66
TDenseSpace::VectorType LocalSystemVectorType
Definition: solution_scheme.hpp:67
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR
Definition: exception.h:161
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484