MAST
stabilized_first_order_transient_sensitivity_solver.h
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2019 Manav Bhatia
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__stabilized_first_order_newmark_transient_sensitivity_solver__
21 #define __mast__stabilized_first_order_newmark_transient_sensitivity_solver__
22 
23 // MAST includes
25 
26 
27 namespace MAST {
28 
29  // Forward declerations
30  class OutputAssemblyElemOperations;
31 
38  public:
40 
42 
47 
49 
53  unsigned int max_index;
54 
58  void set_eigenvalue_stabilization(bool f);
59 
64  void set_nolinear_solution_location(std::string& file_root,
65  std::string& dir);
66 
67 
71  virtual void sensitivity_solve(MAST::AssemblyBase& assembly,
72  const MAST::FunctionBase& f);
73 
74  virtual Real
76  const MAST::FunctionBase& p,
78 
79  virtual void
80  set_element_data(const std::vector<libMesh::dof_id_type>& dof_indices,
81  const std::vector<libMesh::NumericVector<Real>*>& sols);
82 
83 
84  virtual void
85  extract_element_sensitivity_data(const std::vector<libMesh::dof_id_type>& dof_indices,
86  const std::vector<libMesh::NumericVector<Real>*>& sols,
87  std::vector<RealVectorX>& local_sols);
88 
89  virtual void
90  update_sensitivity_velocity(libMesh::NumericVector<Real>& vec,
91  const libMesh::NumericVector<Real>& sol);
92 
93  virtual void
94  elem_calculations(bool if_jac,
95  RealVectorX& vec,
96  RealMatrixX& mat);
97 
98  virtual void
100  RealVectorX& vec);
101 
102  virtual void
103  elem_sensitivity_contribution_previous_timestep(const std::vector<RealVectorX>& prev_sols,
104  RealVectorX& vec);
105 
106  virtual void
107  update_velocity(libMesh::NumericVector<Real>& vel,
108  const libMesh::NumericVector<Real>& sol);
109 
110  virtual void
111  update_acceleration(libMesh::NumericVector<Real>& acc,
112  const libMesh::NumericVector<Real>& sol){
113  libmesh_error();
114  }
115 
116  virtual void
117  update_sensitivity_acceleration(libMesh::NumericVector<Real>& acc,
118  const libMesh::NumericVector<Real>& sol){
119  libmesh_error();
120  }
121 
122  virtual void
123  update_delta_velocity(libMesh::NumericVector<Real>& vel,
124  const libMesh::NumericVector<Real>& sol){
125  libmesh_error();
126  }
127 
128  virtual void
129  update_delta_acceleration(libMesh::NumericVector<Real>& acc,
130  const libMesh::NumericVector<Real>& sol){
131  libmesh_error();
132  }
133 
134  virtual void
136  libmesh_error();
137  }
138 
139  virtual void
141  RealVectorX& vec) {
142  libmesh_error();
143  }
144 
145  virtual void
148  RealVectorX& vec) {
149  libmesh_error();
150  }
151 
152  virtual void
154  libmesh_error();
155  }
156 
157  virtual void
159  (const std::vector<libMesh::dof_id_type>& dof_indices,
160  const std::vector<libMesh::NumericVector<Real>*>& sols) {
161  libmesh_error();
162  }
163  protected:
164 
165  Real _compute_norm_amplification_factor(const libMesh::NumericVector<Real>& sol0,
166  const libMesh::NumericVector<Real>& sol1);
167 
168  Real _compute_eig_amplification_factor(libMesh::SparseMatrix<Real>& A,
169  libMesh::SparseMatrix<Real>& B);
170 
171 
175  unsigned int _index0, _index1;
176 
177  // nonlinear solutions are stored in this directory
178  std::string _sol_name_root, _sol_dir;
179  };
180 
181 }
182 
183 #endif // __mast__stabilized_first_order_newmark_transient_sensitivity_solver__
184 
virtual void elem_topology_sensitivity_calculations(const MAST::FunctionBase &f, const MAST::FieldFunction< RealVectorX > &vel, RealVectorX &vec)
performs the element topology sensitivity calculations over elem, and returns the element residual se...
void set_nolinear_solution_location(std::string &file_root, std::string &dir)
sets the directory where the nonlinear solutions are stored.
virtual void elem_sensitivity_calculations(const MAST::FunctionBase &f, RealVectorX &vec)
performs the element sensitivity calculations over elem, and returns the element residual sensitivity...
virtual void extract_element_sensitivity_data(const std::vector< libMesh::dof_id_type > &dof_indices, const std::vector< libMesh::NumericVector< Real > * > &sols, std::vector< RealVectorX > &local_sols)
void set_eigenvalue_stabilization(bool f)
sets if the eigenvalue-based stabilization will be used.
virtual void elem_linearized_jacobian_solution_product(RealVectorX &vec)
performs the element calculations over elem, and returns the element vector quantity in vec...
This provides the base class for definitin of element level contribution of output quantity in an ana...
libMesh::Real Real
virtual void elem_shape_sensitivity_calculations(const MAST::FunctionBase &f, RealVectorX &vec)
performs the element shape sensitivity calculations over elem, and returns the element residual sensi...
virtual void sensitivity_solve(MAST::AssemblyBase &assembly, const MAST::FunctionBase &f)
solvers the current time step for sensitivity wrt f
virtual void update_sensitivity_velocity(libMesh::NumericVector< Real > &vec, const libMesh::NumericVector< Real > &sol)
update the transient sensitivity velocity based on the current sensitivity solution ...
Real _compute_eig_amplification_factor(libMesh::SparseMatrix< Real > &A, libMesh::SparseMatrix< Real > &B)
Matrix< Real, Dynamic, Dynamic > RealMatrixX
virtual void set_element_data(const std::vector< libMesh::dof_id_type > &dof_indices, const std::vector< libMesh::NumericVector< Real > * > &sols)
Matrix< Real, Dynamic, 1 > RealVectorX
virtual void elem_sensitivity_contribution_previous_timestep(const std::vector< RealVectorX > &prev_sols, RealVectorX &vec)
Real _compute_norm_amplification_factor(const libMesh::NumericVector< Real > &sol0, const libMesh::NumericVector< Real > &sol1)
virtual void update_acceleration(libMesh::NumericVector< Real > &acc, const libMesh::NumericVector< Real > &sol)
update the transient acceleration based on the current solution
virtual void set_element_perturbed_data(const std::vector< libMesh::dof_id_type > &dof_indices, const std::vector< libMesh::NumericVector< Real > * > &sols)
This solver implements the stabilized sensitivity analysis solver for chaotic systems where the linea...
virtual void update_sensitivity_acceleration(libMesh::NumericVector< Real > &acc, const libMesh::NumericVector< Real > &sol)
update the transient sensitivity acceleration based on the current sensitivity solution ...
virtual void update_delta_acceleration(libMesh::NumericVector< Real > &acc, const libMesh::NumericVector< Real > &sol)
update the perturbation in transient acceleration based on the current perturbed solution ...
unsigned int max_index
index of solution that is used for current linearization
virtual void elem_calculations(bool if_jac, RealVectorX &vec, RealMatrixX &mat)
performs the element calculations over elem, and returns the element vector and matrix quantities in ...
virtual void update_delta_velocity(libMesh::NumericVector< Real > &vel, const libMesh::NumericVector< Real > &sol)
update the perturbation in transient velocity based on the current perturbed solution ...
virtual Real evaluate_q_sens_for_previous_interval(MAST::AssemblyBase &assembly, const MAST::FunctionBase &p, MAST::OutputAssemblyElemOperations &output)
virtual void update_velocity(libMesh::NumericVector< Real > &vel, const libMesh::NumericVector< Real > &sol)
update the transient velocity based on the current solution
virtual void elem_second_derivative_dot_solution_assembly(RealMatrixX &mat)
calculates over elem, and returns the matrix in vec .