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_dynamic_scheme.hpp
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 // Project Name: $StructuralMechanicsApplication $
10 // Last modified by: $Author: michael.andre@tum.de $
11 // Date: $Date: September 2016 $
12 // Revision: $Revision: 0.0 $
13 
14 #pragma once
15 
16 
17 // System includes
18 
19 // External includes
20 
21 // Project includes
22 #include "includes/define.h"
23 #include "includes/element.h"
24 #include "includes/condition.h"
25 #include "includes/process_info.h"
28 
29 // Application includes
31 
32 namespace Kratos
33 {
34 
37 
41 
43 
46 
50 
54 
56 template<class TSparseSpace,
57  class TDenseSpace
58  >
59 class EigensolverDynamicScheme : public Scheme<TSparseSpace,TDenseSpace>
60 {
61 public:
64 
66 
68 
70 
72 
76 
78  EigensolverDynamicScheme() : Scheme<TSparseSpace,TDenseSpace>() {}
79 
82 
86 
90 
92  Element& rCurrentElement,
93  LocalSystemMatrixType& LHS_Contribution,
94  LocalSystemVectorType& RHS_Contribution,
96  const ProcessInfo& CurrentProcessInfo
97  ) override
98  {
100 
101  if (CurrentProcessInfo[BUILD_LEVEL] == 1)
102  { // mass matrix
103  rCurrentElement.CalculateMassMatrix(LHS_Contribution,CurrentProcessInfo);
104  std::size_t LocalSize = LHS_Contribution.size1();
105  if (RHS_Contribution.size() != LocalSize)
106  RHS_Contribution.resize(LocalSize,false);
107  noalias(RHS_Contribution) = ZeroVector(LocalSize);
108  }
109  else if (CurrentProcessInfo[BUILD_LEVEL] == 2) // stiffness matrix
110  {
111  rCurrentElement.CalculateLocalSystem(LHS_Contribution,RHS_Contribution,CurrentProcessInfo);
112  }
113  else
114  {
115  KRATOS_ERROR <<"Invalid BUILD_LEVEL" << std::endl;
116  }
117 
118  rCurrentElement.EquationIdVector(EquationId,CurrentProcessInfo);
119 
120  KRATOS_CATCH("")
121  }
122 
124  Element& rCurrentElement,
125  LocalSystemMatrixType& LHS_Contribution,
127  const ProcessInfo& CurrentProcessInfo) override
128  {
129  KRATOS_TRY
130 
131  LocalSystemVectorType RHS_Contribution;
132  RHS_Contribution.resize(LHS_Contribution.size1(), false);
134  rCurrentElement,
135  LHS_Contribution,
136  RHS_Contribution,
137  EquationId,
138  CurrentProcessInfo);
139 
140  KRATOS_CATCH("")
141  }
142 
144  Element& rCurrentElement,
145  LocalSystemVectorType& RHS_Contribution,
147  const ProcessInfo& CurrentProcessInfo) override
148  {
149  KRATOS_TRY
150 
151  rCurrentElement.CalculateRightHandSide(RHS_Contribution,CurrentProcessInfo);
152 
153  rCurrentElement.EquationIdVector(EquationId,CurrentProcessInfo);
154 
155  KRATOS_CATCH("")
156  }
157 
159  Condition& rCurrentCondition,
160  LocalSystemMatrixType& LHS_Contribution,
161  LocalSystemVectorType& RHS_Contribution,
163  const ProcessInfo& CurrentProcessInfo) override
164  {
165  KRATOS_TRY
166 
167  if (CurrentProcessInfo[BUILD_LEVEL] == 1)
168  { // mass matrix
169  rCurrentCondition.CalculateMassMatrix(LHS_Contribution,CurrentProcessInfo);
170  std::size_t LocalSize = LHS_Contribution.size1();
171  if (RHS_Contribution.size() != LocalSize)
172  {
173  RHS_Contribution.resize(LocalSize,false);
174  }
175  noalias(RHS_Contribution) = ZeroVector(LocalSize);
176  }
177  else if (CurrentProcessInfo[BUILD_LEVEL] == 2) // stiffness matrix
178  {
179  rCurrentCondition.CalculateLocalSystem(LHS_Contribution,RHS_Contribution,CurrentProcessInfo);
180  }
181  else
182  {
183  KRATOS_ERROR <<"Invalid BUILD_LEVEL" << std::endl;
184  }
185 
186  rCurrentCondition.EquationIdVector(EquationId,CurrentProcessInfo);
187 
188  KRATOS_CATCH("")
189  }
190 
192  Condition& rCurrentCondition,
193  LocalSystemMatrixType& LHS_Contribution,
195  const ProcessInfo& CurrentProcessInfo) override
196  {
197  KRATOS_TRY
198 
199  LocalSystemVectorType RHS_Contribution;
200  RHS_Contribution.resize(LHS_Contribution.size1(), false);
202  rCurrentCondition,
203  LHS_Contribution,
204  RHS_Contribution,
205  EquationId,
206  CurrentProcessInfo);
207 
208  KRATOS_CATCH("")
209  }
210 
212  Condition& rCurrentCondition,
213  LocalSystemVectorType& RHS_Contribution,
215  const ProcessInfo& CurrentProcessInfo) override
216  {
217  KRATOS_TRY
218 
219  rCurrentCondition.CalculateRightHandSide(RHS_Contribution,CurrentProcessInfo);
220 
221  rCurrentCondition.EquationIdVector(EquationId,CurrentProcessInfo);
222 
223  KRATOS_CATCH("")
224  }
225 
227 
228 }; /* Class Scheme */
229 
231 
234 
236 
237 } /* namespace Kratos.*/
238 
Base class for all Conditions.
Definition: condition.h:59
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
virtual void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:574
virtual void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const
Definition: condition.h:260
virtual void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:440
virtual void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:408
An adapter scheme for obtaining mass and stiffness matrices for dynamic eigenvalue problems.
Definition: eigensolver_dynamic_scheme.hpp:60
Scheme< TSparseSpace, TDenseSpace > BaseType
Definition: eigensolver_dynamic_scheme.hpp:67
KRATOS_CLASS_POINTER_DEFINITION(EigensolverDynamicScheme)
BaseType::LocalSystemVectorType LocalSystemVectorType
Definition: eigensolver_dynamic_scheme.hpp:69
void CalculateSystemContributions(Element &rCurrentElement, LocalSystemMatrixType &LHS_Contribution, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &EquationId, const ProcessInfo &CurrentProcessInfo) override
This function is designed to be called in the builder and solver to introduce the selected time integ...
Definition: eigensolver_dynamic_scheme.hpp:91
void CalculateRHSContribution(Element &rCurrentElement, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &EquationId, const ProcessInfo &CurrentProcessInfo) override
This function is designed to calculate just the RHS contribution.
Definition: eigensolver_dynamic_scheme.hpp:143
~EigensolverDynamicScheme() override
Destructor.
Definition: eigensolver_dynamic_scheme.hpp:81
void CalculateSystemContributions(Condition &rCurrentCondition, LocalSystemMatrixType &LHS_Contribution, LocalSystemVectorType &RHS_Contribution, Condition::EquationIdVectorType &EquationId, const ProcessInfo &CurrentProcessInfo) override
Functions totally analogous to the precedent but applied to the "condition" objects.
Definition: eigensolver_dynamic_scheme.hpp:158
void CalculateRHSContribution(Condition &rCurrentCondition, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &EquationId, const ProcessInfo &CurrentProcessInfo) override
Functions totally analogous to the precedent but applied to the "condition" objects.
Definition: eigensolver_dynamic_scheme.hpp:211
BaseType::LocalSystemMatrixType LocalSystemMatrixType
Definition: eigensolver_dynamic_scheme.hpp:71
void CalculateLHSContribution(Element &rCurrentElement, LocalSystemMatrixType &LHS_Contribution, Element::EquationIdVectorType &EquationId, const ProcessInfo &CurrentProcessInfo) override
This function is designed to calculate just the LHS contribution.
Definition: eigensolver_dynamic_scheme.hpp:123
void CalculateLHSContribution(Condition &rCurrentCondition, LocalSystemMatrixType &LHS_Contribution, Condition::EquationIdVectorType &EquationId, const ProcessInfo &CurrentProcessInfo) override
Functions totally analogous to the precedent but applied to the "condition" objects.
Definition: eigensolver_dynamic_scheme.hpp:191
EigensolverDynamicScheme()
Constructor.
Definition: eigensolver_dynamic_scheme.hpp:78
Base class for all Elements.
Definition: element.h:60
virtual void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:570
virtual void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:437
virtual void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const
Definition: element.h:258
virtual void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:405
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
This class provides the implementation of the basic tasks that are needed by the solution strategy.
Definition: scheme.h:56
virtual void EquationId(const Element &rElement, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo)
This method gets the eqaution id corresponding to the current element.
Definition: scheme.h:636
typename TDenseSpace::VectorType LocalSystemVectorType
Local system vector type definition.
Definition: scheme.h:80
typename TDenseSpace::MatrixType LocalSystemMatrixType
Local system matrix type definition.
Definition: scheme.h:77
#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