MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
second_order_newmark_transient_solver.cpp
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 // MAST includes
23 #include "base/elem_base.h"
25 #include "base/nonlinear_system.h"
26 
27 
28 // libMesh includes
29 #include "libmesh/numeric_vector.h"
30 
31 
34 beta(0.25),
35 gamma(0.5)
36 { }
37 
38 
40 { }
41 
42 
43 
44 
45 void
47 set_element_data(const std::vector<libMesh::dof_id_type>& dof_indices,
48  const std::vector<libMesh::NumericVector<Real>*>& sols){
49 
50  libmesh_assert_equal_to(sols.size(), 3);
51 
52  const unsigned int n_dofs = (unsigned int)dof_indices.size();
53 
54  // get the current state and velocity estimates
55  // also get the current discrete velocity replacement
57  sol = RealVectorX::Zero(n_dofs),
58  vel = RealVectorX::Zero(n_dofs),
59  accel = RealVectorX::Zero(n_dofs);
60 
61 
62  // get the references to current and previous sol and velocity
63  const libMesh::NumericVector<Real>
64  &sol_global = *sols[0],
65  &vel_global = *sols[1],
66  &acc_global = *sols[2];
67 
68  for (unsigned int i=0; i<n_dofs; i++) {
69 
70  sol(i) = sol_global(dof_indices[i]);
71  vel(i) = vel_global(dof_indices[i]);
72  accel(i) = acc_global(dof_indices[i]);
73  }
74 
75  _assembly_ops->set_elem_solution(sol);
76  _assembly_ops->set_elem_velocity(vel);
77  _assembly_ops->set_elem_acceleration(accel);
78 }
79 
80 
81 
82 void
84 extract_element_sensitivity_data(const std::vector<libMesh::dof_id_type>& dof_indices,
85  const std::vector<libMesh::NumericVector<Real>*>& sols,
86  std::vector<RealVectorX>& local_sols) {
87 
88  libmesh_assert_equal_to(sols.size(), 3);
89 
90  const unsigned int n_dofs = (unsigned int)dof_indices.size();
91 
92  local_sols.resize(3);
93 
95  &sol = local_sols[0],
96  &vel = local_sols[1],
97  &accel = local_sols[2];
98 
99  sol = RealVectorX::Zero(n_dofs);
100  vel = RealVectorX::Zero(n_dofs);
101  accel = RealVectorX::Zero(n_dofs);
102 
103 
104  // get the references to current and previous sol and velocity
105  const libMesh::NumericVector<Real>
106  &sol_global = *sols[0],
107  &vel_global = *sols[1],
108  &acc_global = *sols[2];
109 
110  for (unsigned int i=0; i<n_dofs; i++) {
111 
112  sol(i) = sol_global(dof_indices[i]);
113  vel(i) = vel_global(dof_indices[i]);
114  accel(i) = acc_global(dof_indices[i]);
115  }
116 
117 
118 }
119 
120 
121 
122 void
124 set_element_perturbed_data(const std::vector<libMesh::dof_id_type>& dof_indices,
125  const std::vector<libMesh::NumericVector<Real>*>& sols){
126 
127  libmesh_assert_equal_to(sols.size(), 3);
128 
129  const unsigned int n_dofs = (unsigned int)dof_indices.size();
130 
131  // get the current state and velocity estimates
132  // also get the current discrete velocity replacement
134  sol = RealVectorX::Zero(n_dofs),
135  vel = RealVectorX::Zero(n_dofs),
136  accel = RealVectorX::Zero(n_dofs);
137 
138 
139  // get the references to current and previous sol and velocity
140  const libMesh::NumericVector<Real>
141  &sol_global = *sols[0],
142  &vel_global = *sols[1],
143  &acc_global = *sols[2];
144 
145  for (unsigned int i=0; i<n_dofs; i++) {
146 
147  sol(i) = sol_global(dof_indices[i]);
148  vel(i) = vel_global(dof_indices[i]);
149  accel(i) = acc_global(dof_indices[i]);
150  }
151 
152  _assembly_ops->set_elem_perturbed_solution(sol);
153  _assembly_ops->set_elem_perturbed_velocity(vel);
154  _assembly_ops->set_elem_perturbed_acceleration(accel);
155 }
156 
157 
158 
159 
160 void
162 update_velocity(libMesh::NumericVector<Real>& vec,
163  const libMesh::NumericVector<Real>& sol) {
164 
165  const libMesh::NumericVector<Real>
166  &prev_sol = this->solution(1),
167  &prev_vel = this->velocity(1),
168  &prev_acc = this->acceleration(1);
169 
170  // first calculate the acceleration
171  vec = sol;
172  vec.add( -1., prev_sol);
173  vec.scale(gamma/beta/dt);
174  vec.add( 1.-gamma/beta, prev_vel);
175  vec.add((1.-gamma/2./beta)*dt, prev_acc);
176  vec.close();
177 }
178 
179 
180 
181 
182 void
184 update_acceleration(libMesh::NumericVector<Real>& vec,
185  const libMesh::NumericVector<Real>& sol) {
186 
187  const libMesh::NumericVector<Real>
188  &prev_sol = this->solution(1),
189  &prev_vel = this->velocity(1),
190  &prev_acc = this->acceleration(1);
191 
192  // first calculate the acceleration
193  vec = sol;
194  vec.add( -1., prev_sol);
195  vec.scale(1./beta/dt/dt);
196  vec.add( -1./beta/dt, prev_vel);
197  vec.add(-(.5-beta)/beta, prev_acc);
198  vec.close();
199 }
200 
201 
202 
203 void
205 update_sensitivity_velocity(libMesh::NumericVector<Real>& vec,
206  const libMesh::NumericVector<Real>& sol) {
207 
208  const libMesh::NumericVector<Real>
209  &prev_sol = this->solution_sensitivity(1),
210  &prev_vel = this->velocity_sensitivity(1),
211  &prev_acc = this->acceleration_sensitivity(1);
212 
213  // first calculate the acceleration
214  vec = sol;
215  vec.add( -1., prev_sol);
216  vec.scale(gamma/beta/dt);
217  vec.add( 1.-gamma/beta, prev_vel);
218  vec.add((1.-gamma/2./beta)*dt, prev_acc);
219  vec.close();
220 }
221 
222 
223 
224 void
226 update_sensitivity_acceleration(libMesh::NumericVector<Real>& vec,
227  const libMesh::NumericVector<Real>& sol) {
228 
229  const libMesh::NumericVector<Real>
230  &prev_sol = this->solution_sensitivity(1),
231  &prev_vel = this->velocity_sensitivity(1),
232  &prev_acc = this->acceleration_sensitivity(1);
233 
234  // first calculate the acceleration
235  vec = sol;
236  vec.add( -1., prev_sol);
237  vec.scale(1./beta/dt/dt);
238  vec.add( -1./beta/dt, prev_vel);
239  vec.add(-(.5-beta)/beta, prev_acc);
240  vec.close();
241 }
242 
243 
244 
245 
246 void
248 update_delta_velocity(libMesh::NumericVector<Real>& vec,
249  const libMesh::NumericVector<Real>& sol) {
250 
251  // first calculate the acceleration
252  vec = sol;
253  vec.scale( gamma/beta/dt);
254  vec.close();
255 }
256 
257 
258 
259 
260 void
262 update_delta_acceleration(libMesh::NumericVector<Real>& vec,
263  const libMesh::NumericVector<Real>& sol) {
264 
265  // first calculate the acceleration
266  vec = sol;
267  vec.scale( 1./beta/dt/dt);
268  vec.close();
269 }
270 
271 
272 
273 
274 
275 void
277 elem_calculations(bool if_jac,
278  RealVectorX& vec,
279  RealMatrixX& mat) {
280  // make sure that the assembly object is provided
281  libmesh_assert(_assembly);
282  unsigned int n_dofs = (unsigned int)vec.size();
283 
285  f_x = RealVectorX::Zero(n_dofs),
286  f_m = RealVectorX::Zero(n_dofs);
287 
289  f_m_jac_xddot = RealMatrixX::Zero(n_dofs, n_dofs),
290  f_m_jac_xdot = RealMatrixX::Zero(n_dofs, n_dofs),
291  f_m_jac = RealMatrixX::Zero(n_dofs, n_dofs),
292  f_x_jac_xdot = RealMatrixX::Zero(n_dofs, n_dofs),
293  f_x_jac = RealMatrixX::Zero(n_dofs, n_dofs);
294 
295  // perform the element assembly
296  _assembly_ops->elem_calculations(if_jac,
297  f_m, // mass vector
298  f_x, // forcing vector
299  f_m_jac_xddot, // Jac of mass wrt x_dotdot
300  f_m_jac_xdot, // Jac of mass wrt x_dot
301  f_m_jac, // Jac of mass wrt x
302  f_x_jac_xdot, // Jac of forcing vector wrt x_dot
303  f_x_jac); // Jac of forcing vector wrt x
304 
305 
306  if (_if_highest_derivative_solution) {
307 
308  //
309  // The residual here is modeled as
310  // r(x, xdot, xddot) = f_m(x, xdot, xddot) + f_x(x, xdot)= 0
311  //
312  // then, given x = x0, and xdot = xdot0, the residual can be used to
313  // evaluate xddot at the same time step as
314  //
315  // r(x0, xdot0, xddot) + dr/dxddot dxddot = 0
316  //
317 
318  // system residual
319  vec = (f_m + f_x);
320 
321  // system Jacobian
322  if (if_jac)
323  mat = f_m_jac_xddot;
324  }
325  else {
326 
327  // N-R solver: r(x) = 0, and provides estimates for x (current solution)
328  // in an iterative fashion.
329  // Newmark solver uses x and evaluates x_ddot and x_dot (acc and vel at
330  // current x).
331  //
332  // N-R requires r(x) and dr/dx (residual and Jacobian)
333  // N-R asks Newmark to provide these quantities
334  //
335  // Newmark asks elements (which provide physics) for contributions to
336  // f_m and f_x (using the current estimates of x_ddot, x_dot, x)
337  //
338  //
339  // The residual here is modeled as
340  // r = (f_m + f_x )= 0
341  // where, (for example)
342  // f_m = int_Omega phi u_dot [typical mass vector in conduction, for example]
343  // f_x = int_Omega phi_i u_i - int_Gamma phi q_n [typical conductance and heat flux combination, for example]
344  //
345  // This method assumes
346  // x = x0 + dt x0_dot + (1/2-beta) dt^2 x0_ddot + beta dt^2 x_ddot
347  // or, x_ddot = (x-x0)/beta/dt^2 - 1/beta/dt x0_dot - (1/2-beta)/beta x0_ddot
348  //
349  // x_dot = x0_dot + (1-gamma) dt x0_ddot + gamma dt x_ddot
350  // or, x_dot = x0_dot + (1-gamma) dt x0_ddot +
351  // gamma dt [(x-x0)/beta/dt^2 - 1/beta/dt x0_dot - (1/2-beta)/beta x0_ddot]
352  // = x0_dot + (1-gamma) dt x0_ddot +
353  // gamma/beta/dt (x-x0) - gamma/beta x0_dot - gamma (1/2-beta)/beta dt x0_ddot
354  // = gamma/beta/dt (x-x0) + (1 - gamma/beta) x0_dot + (1 - gamma/2/beta) dt x0_ddot
355  //
356  // Both f_m and f_x can be functions of x_dot and x. Then, the
357  // Jacobian is
358  // dr/dx = df_m/dx + df_x/dx +
359  // df_m/dx_dot dx_dot/dx + df_x/dx_dot dx_dot/dx
360  // = (df_m/dx + df_x/dx) +
361  // (df_m/dx_dot + df_x/dx_dot) (gamma/beta/dt) +
362  // df_m/dx_ddot (1/beta/dt/dt) +
363  //
364 
365 
366  // system residual
367  vec = (f_m + f_x);
368 
369  // system Jacobian
370  if (if_jac)
371  mat = ((1./beta/dt/dt) * f_m_jac_xddot +
372  (gamma/beta/dt)* (f_m_jac_xdot+f_x_jac_xdot) +
373  (f_m_jac + f_x_jac));
374  }
375 }
376 
377 
378 
379 void
382 
383  // make sure that the assembly object is provided
384  libmesh_assert(_assembly_ops);
385 
386  // perform the element assembly
387  _assembly_ops->linearized_jacobian_solution_product(vec);
388 }
389 
390 
391 
392 
393 void
396  RealVectorX& vec) {
397 
398  // make sure that the assembly object is provided
399  libmesh_assert(_assembly_ops);
400  unsigned int n_dofs = (unsigned int)vec.size();
401 
403  f_x = RealVectorX::Zero(n_dofs),
404  f_m = RealVectorX::Zero(n_dofs);
405 
406  // perform the element assembly
407  _assembly_ops->elem_sensitivity_calculations(f,
408  f_m, // mass vector
409  f_x); // forcing vector
410 
411  // system residual
412  vec = (f_m + f_x);
413 }
414 
415 
416 void
418 elem_sensitivity_contribution_previous_timestep(const std::vector<RealVectorX>& prev_sols,
419  RealVectorX& vec) {
420 
421  // nothing to be done for highest derivative term
422  if (_if_highest_derivative_solution) return;
423 
424  // make sure that the assembly object is provided
425  libmesh_assert(_assembly_ops);
426  libmesh_assert_equal_to(prev_sols.size(), 3);
427 
428  unsigned int n_dofs = (unsigned int)vec.size();
429 
430  const RealVectorX
431  &sol = prev_sols[0],
432  &vel = prev_sols[1],
433  &accel = prev_sols[2];
435  dummy_vec = RealVectorX::Zero(n_dofs);
436 
438  f_m_jac_xdot = RealMatrixX::Zero(n_dofs, n_dofs),
439  dummy_mat = RealMatrixX::Zero(n_dofs, n_dofs);
440 
441  // perform the element assembly
442  _assembly_ops->elem_calculations(true,
443  dummy_vec, // mass vector
444  dummy_vec, // forcing vector
445  f_m_jac_xdot, // Jac of mass wrt x_dot
446  dummy_mat, // Jac of mass wrt x
447  dummy_mat); // Jac of forcing vector wrt x
448 
449  libmesh_assert(false); // this expression is in error and needs to be set
450  vec -= f_m_jac_xdot * ( (1./beta/dt)*sol + (1.-beta)/beta * vel);
451 }
452 
453 
454 void
457  RealVectorX& vec) {
458 
459  libmesh_assert(false); // to be implemented
460 }
461 
462 
463 
464 void
467  RealVectorX& vec) {
468  libmesh_assert(false); // to be implemented
469 }
470 
471 
472 
473 void
477  RealVectorX& vec) {
478  libmesh_assert(false); // to be implemented
479 }
480 
481 
482 
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
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 ...
virtual void update_acceleration(libMesh::NumericVector< Real > &vec, const libMesh::NumericVector< Real > &sol)
update the transient acceleration based on the current solution
virtual void update_sensitivity_velocity(libMesh::NumericVector< Real > &vel, const libMesh::NumericVector< Real > &sol)
update the transient sensitivity velocity based on the current sensitivity solution ...
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 ...
libMesh::NumericVector< Real > & solution(unsigned int prev_iter=0) const
virtual void elem_sensitivity_calculations(const MAST::FunctionBase &f, RealVectorX &vec)=0
performs the element sensitivity calculations over elem, and returns the element residual sensitivity...
libMesh::NumericVector< Real > & solution_sensitivity(unsigned int prev_iter=0) const
libMesh::NumericVector< Real > & acceleration_sensitivity(unsigned int prev_iter=0) const
virtual void set_element_perturbed_data(const std::vector< libMesh::dof_id_type > &dof_indices, const std::vector< libMesh::NumericVector< Real > *> &sols)
provides the element with the transient data for calculations
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 ...
Matrix< Real, Dynamic, Dynamic > RealMatrixX
libMesh::NumericVector< Real > & velocity(unsigned int prev_iter=0) const
virtual void elem_shape_sensitivity_calculations(const MAST::FunctionBase &f, RealVectorX &vec)=0
performs the element shape sensitivity calculations over elem, and returns the element residual sensi...
libMesh::NumericVector< Real > & acceleration(unsigned int prev_iter=0) const
virtual void elem_topology_sensitivity_calculations(const MAST::FunctionBase &f, RealVectorX &vec)=0
performs the element topology sensitivity calculations over elem, and returns the element residual se...
Matrix< Real, Dynamic, 1 > RealVectorX
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_velocity(libMesh::NumericVector< Real > &vec, const libMesh::NumericVector< Real > &sol)
update the transient velocity based on the current solution
virtual void set_element_data(const std::vector< libMesh::dof_id_type > &dof_indices, const std::vector< libMesh::NumericVector< Real > *> &sols)
provides the element with the transient data for calculations
libMesh::NumericVector< Real > & velocity_sensitivity(unsigned int prev_iter=0) const
virtual void elem_linearized_jacobian_solution_product(RealVectorX &vec)=0
performs the element calculations over elem, and returns the element vector quantity in vec...
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)
provides the element with the sensitivity of transient data for calculations