MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
structural_buckling_eigenproblem_assembly.h
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2020 Manav Bhatia and MAST authors
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 #ifndef __mast__structural_buckling_eigenproblem_assembly__
21 #define __mast__structural_buckling_eigenproblem_assembly__
22 
23 
24 // MAST includes
26 
27 namespace MAST {
28 
29  // Forward declerations
30  class Parameter;
31 
32 
35  public:
36 
41 
47 
48 
53  virtual void
54  eigenproblem_assemble(libMesh::SparseMatrix<Real>* A,
55  libMesh::SparseMatrix<Real>* B);
56 
57 
76  void set_buckling_data (bool use_linearized_approach,
77  MAST::Parameter& p,
78  const Real lambda1,
79  const Real lambda2,
80  libMesh::NumericVector<Real>& x1,
81  libMesh::NumericVector<Real>& x2);
82 
83 
88 
89 
93  virtual void clear_discipline_and_system();
94 
95 
96  protected:
97 
101  std::map<const libMesh::Elem*, RealVectorX> _incompatible_sol;
102 
103 
108 
113 
114 
120 
125  libMesh::NumericVector<Real> *_sol1, *_sol2;
126 
127 
131  libMesh::NumericVector<Real> *_sol1_sens, *_sol2_sens;
132  };
133 
134 }
135 
136 
137 #endif // __mast__structural_buckling_eigenproblem_assembly__
StructuralBucklingEigenproblemAssembly()
constructor associates the eigen system with this assembly object
libMesh::NumericVector< Real > * _sol1_sens
the equilibrium solution sensitivity
This is a scalar function whose value can be changed and one that can be used as a design variable in...
Definition: parameter.h:35
libMesh::Real Real
virtual void clear_discipline_and_system()
clears the states and load factors for buckling eigenproblem
bool _use_linearized_formulation
whether or not to use the linearized formulation
Assembles the system of equations for an eigenproblem of type .
Real _lambda1
values of load factors for which the two stiffness matrices are calculated for the buckling eigenvalu...
virtual ~StructuralBucklingEigenproblemAssembly()
destructor resets the association with the eigen system from this assembly object ...
Real critical_point_estimate_from_eigenproblem(Real v) const
calculates the critical load factor based on the eigensolution
std::map< const libMesh::Elem *, RealVectorX > _incompatible_sol
map of local incompatible mode solution per 3D elements
void set_buckling_data(bool use_linearized_approach, MAST::Parameter &p, const Real lambda1, const Real lambda2, libMesh::NumericVector< Real > &x1, libMesh::NumericVector< Real > &x2)
set the states and load factors for buckling eigenproblem.
virtual void eigenproblem_assemble(libMesh::SparseMatrix< Real > *A, libMesh::SparseMatrix< Real > *B)
assembles the global A and B matrices for the modal eigenvalue problem
MAST::Parameter * _load_param
load parameter used to define the two stiffness matrices
libMesh::NumericVector< Real > * _sol1
the equilibrium solutions associated with _lambda1 and _lambda2 load factors.