MAST
transient_solver_base.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 
21 #ifndef __mast__transient_solver_base__
22 #define __mast__transient_solver_base__
23 
24 // MAST includes
25 #include "base/mast_data_types.h"
27 
28 // libMesh includes
29 #include "libmesh/numeric_vector.h"
30 
31 
32 namespace MAST {
33 
34  // Forward declerations
35  class TransientAssemblyElemOperations;
36  class ElementBase;
37  class NonlinearSystem;
38 
39 
42  public:
47  TransientSolverBase(unsigned int o,
48  unsigned int n);
49 
50  virtual ~TransientSolverBase();
51 
52 
53  virtual void
55 
56 
57  virtual void clear_assembly();
58 
64 
67 
68 
72  virtual void clear_elem_operation_object();
73 
78 
79 
87  libMesh::NumericVector<Real>&
88  solution(unsigned int prev_iter = 0) const;
89 
94  libMesh::NumericVector<Real>&
95  solution_sensitivity(unsigned int prev_iter = 0) const;
96 
104  libMesh::NumericVector<Real>&
105  velocity(unsigned int prev_iter = 0) const;
106 
107 
112  libMesh::NumericVector<Real>&
113  velocity_sensitivity(unsigned int prev_iter = 0) const;
114 
115 
123  libMesh::NumericVector<Real>&
124  acceleration(unsigned int prev_iter = 0) const;
125 
130  libMesh::NumericVector<Real>&
131  acceleration_sensitivity(unsigned int prev_iter = 0) const;
132 
133 
137  virtual void update_velocity(libMesh::NumericVector<Real>& vel,
138  const libMesh::NumericVector<Real>& sol) = 0;
139 
143  virtual void update_acceleration(libMesh::NumericVector<Real>& acc,
144  const libMesh::NumericVector<Real>& sol) = 0;
145 
150  virtual void update_sensitivity_velocity(libMesh::NumericVector<Real>& vel,
151  const libMesh::NumericVector<Real>& sol) = 0;
152 
157  virtual void update_sensitivity_acceleration(libMesh::NumericVector<Real>& acc,
158  const libMesh::NumericVector<Real>& sol) = 0;
159 
160 
165  virtual void
166  update_delta_velocity(libMesh::NumericVector<Real>& vel,
167  const libMesh::NumericVector<Real>& sol) = 0;
168 
173  virtual void
174  update_delta_acceleration(libMesh::NumericVector<Real>& acc,
175  const libMesh::NumericVector<Real>& sol) = 0;
176 
180  virtual void solve(MAST::AssemblyBase& assembly);
181 
185  virtual void sensitivity_solve(MAST::AssemblyBase& assembly,
186  const MAST::FunctionBase& f);
187 
188 
197 
202  void
204  const MAST::FunctionBase& f);
205 
206 
207 
212  virtual void advance_time_step();
213 
219  virtual void advance_time_step_with_sensitivity();
220 
221 
226  virtual void
227  build_local_quantities(const libMesh::NumericVector<Real>& current_sol,
228  std::vector<libMesh::NumericVector<Real>*>& qtys);
229 
237  virtual void
238  build_sensitivity_local_quantities(unsigned int prev_iter,
239  std::vector<libMesh::NumericVector<Real>*>& qtys);
240 
249  virtual void
250  build_perturbed_local_quantities
251  (const libMesh::NumericVector<Real>& current_sol,
252  std::vector<libMesh::NumericVector<Real>*>& qtys);
253 
254 
258  virtual void
259  set_element_data(const std::vector<libMesh::dof_id_type>& dof_indices,
260  const std::vector<libMesh::NumericVector<Real>*>& sols) = 0;
261 
262 
267  virtual void
268  extract_element_sensitivity_data(const std::vector<libMesh::dof_id_type>& dof_indices,
269  const std::vector<libMesh::NumericVector<Real>*>& sols,
270  std::vector<RealVectorX>& local_sols) = 0;
271 
276  virtual void
277  elem_sensitivity_contribution_previous_timestep(const std::vector<RealVectorX>& prev_sols,
278  RealVectorX& vec) = 0;
279 
283  virtual void
284  set_element_perturbed_data
285  (const std::vector<libMesh::dof_id_type>& dof_indices,
286  const std::vector<libMesh::NumericVector<Real>*>& sols) = 0;
287 
288 
292  virtual void
293  set_elem_data(unsigned int dim,
294  const libMesh::Elem& ref_elem,
295  MAST::GeomElem& elem) const;
296 
300  virtual void
301  init(const MAST::GeomElem& elem);
302 
306  virtual void clear_elem();
307 
308  protected:
309 
313  bool _first_step, _first_sensitivity_step;
314 
315 
320  const unsigned int _ode_order;
321 
326  const unsigned int _n_iters_to_store;
327 
328 
334 
339  bool _if_highest_derivative_solution;
340 
341  };
342 
343 }
344 
345 
346 #endif // __mast__transient_solver_base__
virtual void advance_time_step_with_sensitivity()
advances the time step and copies the current sensitivity solution to old sensitivity solution...
virtual void advance_time_step()
advances the time step and copies the current solution to old solution, and so on.
virtual void elem_sensitivity_contribution_previous_timestep(const std::vector< RealVectorX > &prev_sols, RealVectorX &vec)
computes the contribution for this element from previous time step
void solve_highest_derivative_and_advance_time_step(MAST::AssemblyBase &assembly)
To be used only for initial conditions.
virtual void sensitivity_solve(MAST::AssemblyBase &assembly, const MAST::FunctionBase &f)
solvers the current time step for sensitivity wrt f
virtual void update_sensitivity_acceleration(libMesh::NumericVector< Real > &acc, const libMesh::NumericVector< Real > &sol)=0
update the transient sensitivity acceleration based on the current sensitivity solution ...
virtual void set_assembly(MAST::AssemblyBase &assembly)
sets the assembly object
virtual MAST::TransientAssemblyElemOperations & get_elem_operation_object()
virtual void update_delta_acceleration(libMesh::NumericVector< Real > &acc, const libMesh::NumericVector< Real > &sol)=0
update the perturbation in transient acceleration based on the current perturbed solution ...
virtual void update_delta_velocity(libMesh::NumericVector< Real > &vel, const libMesh::NumericVector< Real > &sol)=0
update the perturbation in transient velocity based on the current perturbed solution ...
libMesh::Real Real
virtual void set_elem_operation_object(MAST::TransientAssemblyElemOperations &elem_ops)
Attaches the assembly elem operations object that provides the x_dot, M and J quantities for the elem...
libMesh::NumericVector< Real > & solution_sensitivity(unsigned int prev_iter=0) const
libMesh::NumericVector< Real > & velocity(unsigned int prev_iter=0) const
virtual void set_elem_data(unsigned int dim, const libMesh::Elem &ref_elem, MAST::GeomElem &elem) const =0
some analyses may want to set additional element data before initialization of the GeomElem...
virtual void update_velocity(libMesh::NumericVector< Real > &vel, const libMesh::NumericVector< Real > &sol)=0
update the transient velocity based on the current solution
virtual void build_local_quantities(const libMesh::NumericVector< Real > &current_sol, std::vector< libMesh::NumericVector< Real > * > &qtys)
localizes the relevant solutions for system assembly.
virtual void build_sensitivity_local_quantities(unsigned int prev_iter, std::vector< libMesh::NumericVector< Real > * > &qtys)
localizes the relevant solutions for system assembly.
libMesh::NumericVector< Real > & velocity_sensitivity(unsigned int prev_iter=0) const
Matrix< Real, Dynamic, 1 > RealVectorX
virtual void init(const MAST::GeomElem &elem)=0
initializes the object for calculation of element quantities for the specified elem.
virtual void update_acceleration(libMesh::NumericVector< Real > &acc, const libMesh::NumericVector< Real > &sol)=0
update the transient acceleration based on the current solution
libMesh::NumericVector< Real > & acceleration_sensitivity(unsigned int prev_iter=0) const
virtual void clear_assembly()
clears the assembly object
libMesh::NumericVector< Real > & solution(unsigned int prev_iter=0) const
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
virtual void solve(MAST::AssemblyBase &assembly)
solves the current time step for solution and velocity
virtual void clear_elem()
clears the element initialization
TransientSolverBase(unsigned int o, unsigned int n)
constructor requires the number of iterations to store for the derived solver.
virtual void update_sensitivity_velocity(libMesh::NumericVector< Real > &vel, const libMesh::NumericVector< Real > &sol)=0
update the transient sensitivity velocity based on the current sensitivity solution ...
void solve_highest_derivative_and_advance_time_step_with_sensitivity(MAST::AssemblyBase &assembly, const MAST::FunctionBase &f)
solves for the sensitivity of highest derivative and advances the time-step.
virtual void clear_elem_operation_object()
Clears the assembly elem operations object.
libMesh::NumericVector< Real > & acceleration(unsigned int prev_iter=0) const