forked from smdogroup/tacs
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathTACSBuckling.h
192 lines (156 loc) · 6.16 KB
/
TACSBuckling.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
/*
This file is part of TACS: The Toolkit for the Analysis of Composite
Structures, a parallel finite-element code for structural and
multidisciplinary design optimization.
Copyright (C) 2010 University of Toronto
Copyright (C) 2012 University of Michigan
Copyright (C) 2014 Georgia Tech Research Corporation
Additional copyright (C) 2010 Graeme J. Kennedy and Joaquim
R.R.A. Martins All rights reserved.
TACS is licensed under the Apache License, Version 2.0 (the
"License"); you may not use this software except in compliance with
the License. You may obtain a copy of the License at
http://www.apache.org/licenses/LICENSE-2.0
*/
#ifndef TACS_BUCKLING_H
#define TACS_BUCKLING_H
/*
Implementation of buckling and frequency analysis and sensitivity
analysis of eigenvalues.
*/
#include "TACSAssembler.h"
#include "GSEP.h"
#include "TACSMg.h"
#include "JacobiDavidson.h"
/*
Linearized buckling analysis code.
The efficient solution of generalized eigenvalue problems requires a
shift and invert operation that involves among other things, copying
and axpy operations on matrices. These operations are not supported
with the TACSMat interface because it would be difficult to do this
for matrices of different alternate types. This intermediate object
maintains consistency between matrix types involved in the operation
without exposing the underlying matrix type.
*/
class TACSLinearBuckling : public TACSObject {
public:
TACSLinearBuckling( TACSAssembler *_tacs,
TacsScalar _sigma,
TACSMat *_gmat, TACSMat *_kmat,
TACSMat *_aux_mat, TACSKsm *_solver,
int _max_lanczos_vecs,
int _num_eigvals, double _eig_tol );
~TACSLinearBuckling();
// Retrieve the instance of TACSAssembler
// --------------------------------------
TACSAssembler* getTACS(){ return tacs; }
// Functions to set the shift value
// --------------------------------
TacsScalar getSigma();
void setSigma( TacsScalar sigma );
// Solve the eigenvalue problem
// ----------------------------
void solve( TACSVec *rhs=NULL, KSMPrint *ksm_print=NULL );
void evalEigenDVSens( int n, TacsScalar fdvSens[], int numDVs );
// Extract the eigenvalue or check the solution
// --------------------------------------------
TacsScalar extractEigenvalue( int n, TacsScalar *error );
TacsScalar extractEigenvector( int n, TACSBVec *ans, TacsScalar *error );
void checkEigenvector( int n );
TacsScalar checkOrthogonality();
void printOrthogonality();
private:
// Data for the eigenvalue analysis
TacsScalar sigma;
EPBucklingShiftInvert *ep_op;
SEP *sep;
// The tacs object
TACSAssembler *tacs;
// Tolerances/required number of eigenvalues
int max_lanczos_vecs;
int num_eigvals;
double eig_tol;
// These are used by the eigenvalue solver and to solve the linear systems
// for the path determination problem
TACSPc *pc;
TACSKsm *solver;
TACSMat *aux_mat, *kmat, *gmat;
// Vectors used in the analysis
TACSBVec *path; // The solution path
TACSBVec *res, *update, *eigvec;
// The multigrid object -- only defined if a multigrid
// preconditioner is used
TACSMg *mg;
};
/*!
The following class performs frequency analysis and gradient
evaluation of a TACS finite-element model.
The code computes eigenvalues and eigenvectors of the generalized
eigenproblem:
K u = lambda M u
The natural frequency of vibration are determined where lambda =
omega^2.
The code uses a Lanczos eigenproblem solver with full
orthogonalization. The full orthogonalization ensures that the
Lanczos basis is linearly independent to the required precision. The
derivatives of the eigenvalues are obtained using an efficient
method for computing the derivative of the inner product of two
vectors and the corresponding matrix.
*/
class TACSFrequencyAnalysis : public TACSObject {
public:
TACSFrequencyAnalysis( TACSAssembler *_tacs,
TacsScalar _sigma,
TACSMat *_mmat, TACSMat *_kmat,
TACSKsm *_solver, int max_lanczos,
int num_eigvals, double _eig_tol );
TACSFrequencyAnalysis( TACSAssembler *_tacs,
TacsScalar _sigma,
TACSMat *_mmat, TACSMat *_kmat,
TACSMat *_pcmat, TACSPc *_pc,
int max_jd_size,
int fgmres_size, int num_eigvals,
double eigtol=1e-9,
double eig_rtol=1e-9,
double eig_atol=1e-30,
int num_recycle=0,
JDRecycleType recycle_type=JD_NUM_RECYCLE );
~TACSFrequencyAnalysis();
// Retrieve the instance of TACSAssembler
// --------------------------------------
TACSAssembler* getTACS(){ return tacs; }
// Solve the generalized eigenvalue problem
// ----------------------------------------
TacsScalar getSigma();
void setSigma( TacsScalar _sigma );
void solve( KSMPrint *ksm_print=NULL, KSMPrint *ksm_file=NULL );
void evalEigenDVSens( int n, TacsScalar fdvSens[], int numDVs );
// Extract and check the solution
// ------------------------------
TacsScalar extractEigenvalue( int n, TacsScalar *error );
TacsScalar extractEigenvector( int n, TACSBVec *ans, TacsScalar *error );
void checkEigenvector( int n );
TacsScalar checkOrthogonality();
private:
// The TACS assembler object
TACSAssembler *tacs;
// The matrices used in the analysis
TACSMat *mmat; // The mass matrix
TACSMat *kmat; // The stiffness matrix
TACSMat *pcmat; // Matrix associated with the preconditioner (JD)
TACSKsm *solver; // Associated with kmat
TACSPc *pc; // The preconditioner
// The multigrid object -- only defined if a multigrid
// preconditioner is used
TACSMg *mg;
// The eigen solver
TacsScalar sigma;
EPGeneralizedShiftInvert *ep_op;
SEP *sep;
// Objects associated with the Jacobi-Davidson method
TACSJDFrequencyOperator *jd_op;
TACSJacobiDavidson *jd;
// Vectors required for eigen-sensitivity analysis
TACSBVec *eigvec, *res;
};
#endif // TACS_BUCKLING_H