MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
structural_element_base.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_element_base__
21 #define __mast__structural_element_base__
22 
23 // C++ includes
24 #include <memory>
25 #include <map>
26 
27 // MAST includes
28 #include "base/elem_base.h"
29 
30 
31 namespace MAST {
32 
33  // Forward declerations
34  class GeomElem;
35  class ElementPropertyCardBase;
36  class BoundaryConditionBase;
37  class FEMOperatorMatrix;
38  class StressStrainOutputBase;
39  template <typename ValType> class FieldFunction;
40 
41 
43  public MAST::ElementBase
44  {
45  public:
50  const MAST::GeomElem& elem,
52 
53  virtual ~StructuralElementBase();
54 
55 
60  virtual void set_solution(const RealVectorX& vec,
61  bool if_sens = false);
62 
63 
68  virtual void set_perturbed_solution(const RealVectorX& vec,
69  bool if_sens = false);
70 
71 
76  virtual void set_velocity(const RealVectorX& vec,
77  bool if_sens = false);
78 
83  virtual void set_perturbed_velocity(const RealVectorX& vec,
84  bool if_sens = false);
85 
86 
91  virtual void set_acceleration(const RealVectorX& vec,
92  bool if_sens = false);
93 
94 
99  virtual void set_perturbed_acceleration(const RealVectorX& vec,
100  bool if_sens = false);
101 
102 
108  const RealVectorX& local_solution(bool if_sens = false) const {
109 
110  if (!if_sens)
111  return _local_sol;
112  else
113  return _local_sol_sens;
114  }
115 
120  return _property;
121  }
122 
123 
127  virtual bool internal_residual (bool request_jacobian,
128  RealVectorX& f,
129  RealMatrixX& jac) = 0;
130 
135  virtual bool
136  linearized_internal_residual (bool request_jacobian,
137  RealVectorX& f,
138  RealMatrixX& jac);
139 
145  virtual void
147  const unsigned int s,
149  bool request_jacobian,
150  RealVectorX& f,
151  RealMatrixX& jac) = 0;
152 
156  virtual bool
158 
159 
160 
164  virtual bool damping_residual (bool request_jacobian,
165  RealVectorX& f,
166  RealMatrixX& jac_xdot,
167  RealMatrixX& jac) {
168 
169  // to be implemented nothing is done yet
170  return false;
171  }
172 
176  virtual bool inertial_residual (bool request_jacobian,
177  RealVectorX& f,
178  RealMatrixX& jac_xddot,
179  RealMatrixX& jac_xdot,
180  RealMatrixX& jac);
181 
186  virtual bool linearized_inertial_residual (bool request_jacobian,
187  RealVectorX& f,
188  RealMatrixX& jac_xddot,
189  RealMatrixX& jac_xdot,
190  RealMatrixX& jac);
191 
192 
199  bool side_external_residual (bool request_jacobian,
200  RealVectorX& f,
201  RealMatrixX& jac_xdot,
202  RealMatrixX& jac,
203  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
204 
211  bool
213  (bool request_jacobian,
214  RealVectorX& f,
215  RealMatrixX& jac_xdot,
216  RealMatrixX& jac,
217  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
218 
219 
220 
227  bool
229  (bool request_jacobian,
230  ComplexVectorX& f,
231  ComplexMatrixX& jac,
232  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
233 
234 
238  virtual bool prestress_residual (bool request_jacobian,
239  RealVectorX& f,
240  RealMatrixX& jac) = 0;
241 
248  bool volume_external_residual (bool request_jacobian,
249  RealVectorX& f,
250  RealMatrixX& jac_xdot,
251  RealMatrixX& jac,
252  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
253 
254 
259  (const MAST::FunctionBase& p,
260  const unsigned int s,
262  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc,
263  bool request_jacobian,
264  RealVectorX& f,
265  RealMatrixX& jac);
266 
267 
274  bool
276  (bool request_jacobian,
277  RealVectorX& f,
278  RealMatrixX& jac_xdot,
279  RealMatrixX& jac,
280  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
281 
282 
288  bool
290  (bool request_jacobian,
291  ComplexVectorX& f,
292  ComplexMatrixX& jac,
293  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
294 
295 
296 
300  virtual bool internal_residual_sensitivity (const MAST::FunctionBase& p,
301  bool request_jacobian,
302  RealVectorX& f,
303  RealMatrixX& jac) = 0;
304 
309  bool request_jacobian,
310  RealVectorX& f,
311  RealMatrixX& jac_xdot,
312  RealMatrixX& jac) {
313  return false;
314  }
315 
319  virtual bool inertial_residual_sensitivity (const MAST::FunctionBase& p,
320  bool request_jacobian,
321  RealVectorX& f,
322  RealMatrixX& jac_xddot,
323  RealMatrixX& jac_xdot,
324  RealMatrixX& jac);
325 
329  virtual void
331  const unsigned int s,
333  bool request_jacobian,
334  RealVectorX& f,
335  RealMatrixX& jac_xddot,
336  RealMatrixX& jac_xdot,
337  RealMatrixX& jac);
338 
343  bool request_jacobian,
344  RealVectorX& f,
345  RealMatrixX& jac_xdot,
346  RealMatrixX& jac,
347  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
348 
352  virtual bool prestress_residual_sensitivity (const MAST::FunctionBase& p,
353  bool request_jacobian,
354  RealVectorX& f,
355  RealMatrixX& jac) = 0;
356 
361  bool request_jacobian,
362  RealVectorX& f,
363  RealMatrixX& jac_xdot,
364  RealMatrixX& jac,
365  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
366 
367 
368 
373  virtual bool if_incompatible_modes() const = 0;
374 
375 
379  virtual unsigned int incompatible_mode_size() const {
380  // needs to be implemented only in inherited elements
381  // that use incompatible modes
382  libmesh_error();
383  return 0;
384  }
385 
386 
393  _incompatible_sol = &vec;
394  }
395 
396 
402  virtual void update_incompatible_mode_solution(const RealVectorX& dsol) {
403  // should be implemented in derived classes, if relevant
404  libmesh_error();
405  }
406 
407 
413  virtual bool calculate_stress(bool request_derivative,
414  const MAST::FunctionBase* f,
415  MAST::StressStrainOutputBase& output) = 0;
416 
417 
422  virtual void
425  const unsigned int s,
426  const MAST::FieldFunction<RealVectorX>& vel_f) = 0;
427 
428  virtual void
430  MAST::StressStrainOutputBase& output) = 0;
431 
432  virtual void
434  RealMatrixX& m) = 0;
435 
440 
441 
442  template <typename ValType>
443  void transform_matrix_to_global_system(const ValType& local_mat,
444  ValType& global_mat) const;
445 
446  template <typename ValType>
447  void transform_vector_to_local_system(const ValType& global_vec,
448  ValType& local_vec) const;
449 
450  template <typename ValType>
451  void transform_vector_to_global_system(const ValType& local_vec,
452  ValType& global_vec) const;
453 
454  protected:
455 
459  virtual bool
460  surface_traction_residual(bool request_jacobian,
461  RealVectorX& f,
462  RealMatrixX& jac,
463  const unsigned int side,
465 
470  virtual bool
472  bool request_jacobian,
473  RealVectorX& f,
474  RealMatrixX& jac,
475  const unsigned int side,
477 
478 
483  virtual bool
484  surface_traction_residual_shifted_boundary(bool request_jacobian,
485  RealVectorX& f,
486  RealMatrixX& jac,
487  const unsigned int side,
489 
494  virtual bool
496  bool request_jacobian,
497  RealVectorX& f,
498  RealMatrixX& jac,
499  const unsigned int side,
501 
505  virtual bool
506  surface_pressure_residual(bool request_jacobian,
507  RealVectorX& f,
508  RealMatrixX& jac,
509  const unsigned int side,
511 
512 
518  virtual bool
519  surface_pressure_residual(bool request_jacobian,
520  RealVectorX& f,
521  RealMatrixX& jac,
523 
529  virtual void
531  const unsigned int s,
534  bool request_jacobian,
535  RealVectorX& f,
536  RealMatrixX& jac);
537 
543  virtual bool
544  linearized_surface_pressure_residual(bool request_jacobian,
545  RealVectorX& f,
546  RealMatrixX& jac,
548 
549 
554  Real piston_theory_cp(const unsigned int order,
555  const Real vel_U,
556  const Real gamma,
557  const Real mach);
558 
559 
560 
565  Real piston_theory_dcp_dv(const unsigned int order,
566  const Real vel_U,
567  const Real gamma,
568  const Real mach);
569 
570 
578  virtual bool piston_theory_residual(bool request_jacobian,
579  RealVectorX &f,
580  RealMatrixX& jac_xdot,
581  RealMatrixX& jac,
583 
590  virtual bool piston_theory_residual(bool request_jacobian,
591  RealVectorX &f,
592  RealMatrixX& jac_xdot,
593  RealMatrixX& jac,
594  const unsigned int side,
596 
597 
598 
607  bool request_jacobian,
608  RealVectorX &f,
609  RealMatrixX& jac_xdot,
610  RealMatrixX& jac,
612 
620  bool request_jacobian,
621  RealVectorX &f,
622  RealMatrixX& jac_xdot,
623  RealMatrixX& jac,
624  const unsigned int side,
626 
627 
632  virtual bool
634  (bool request_jacobian,
635  ComplexVectorX& f,
636  ComplexMatrixX& jac,
637  const unsigned int side,
639 
640 
649  virtual bool
651  (bool request_jacobian,
652  ComplexVectorX& f,
653  ComplexMatrixX& jac,
655 
656 
661  virtual bool
663  (const MAST::FunctionBase& p,
664  bool request_jacobian,
665  ComplexVectorX& f,
666  ComplexMatrixX& jac,
667  const unsigned int side,
669 
670 
676  virtual bool
679  bool request_jacobian,
680  ComplexVectorX& f,
681  ComplexMatrixX& jac,
683 
684  libmesh_error(); // to be implemented
685  }
686 
687 
691  virtual bool
693  bool request_jacobian,
694  RealVectorX& f,
695  RealMatrixX& jac,
696  const unsigned int side,
698 
703  virtual bool
705  bool request_jacobian,
706  RealVectorX& f,
707  RealMatrixX& jac,
709 
710 
711 
716  virtual bool thermal_residual(bool request_jacobian,
717  RealVectorX& f,
718  RealMatrixX& jac,
720 
725  virtual bool thermal_residual_sensitivity(const MAST::FunctionBase& p,
726  bool request_jacobian,
727  RealVectorX& f,
728  RealMatrixX& jac,
730 
736  const unsigned int s,
739  bool request_jacobian,
740  RealVectorX& f,
741  RealMatrixX& jac) = 0;
742 
747 
748 
753 
754 
759 
760 
765 
766 
771 
772 
777 
778 
783 
784 
789 
790 
795 
796 
801 
802 
807 
808 
813 
814 
819 
820 
825 
826  };
827 
828 
832  std::unique_ptr<MAST::StructuralElementBase>
834  const MAST::GeomElem& elem,
836 
837 }
838 
839 
840 #endif // __mast__structural_element_base__
RealVectorX _local_accel_sens
local acceleration sensitivity
virtual void internal_residual_boundary_velocity(const MAST::FunctionBase &p, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
calculates the term on side s: .
bool linearized_frequency_domain_volume_external_residual(bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc)
Calculates the frequency domain volume external force contribution to system residual.
virtual void surface_pressure_boundary_velocity(const MAST::FunctionBase &p, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
Calculates the force vector and Jacobian due to surface pressure applied on the entire element domain...
virtual bool thermal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc)=0
Calculates the force vector and Jacobian due to thermal stresses.
virtual bool linearized_frequency_domain_surface_pressure_residual(bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0
Calculates the force vector and Jacobian due to small perturbation surface pressure.
virtual bool thermal_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc)=0
Calculates the sensitivity of force vector and Jacobian due to thermal stresses.
virtual void update_incompatible_mode_solution(const RealVectorX &dsol)
updates the incompatible solution for this element.
Data structure provides the mechanism to store stress and strain output from a structural analysis...
RealVectorX _local_delta_vel_sens
local perturbed velocity sensitivity
virtual bool inertial_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xddot, RealMatrixX &jac_xdot, RealMatrixX &jac)
sensitivity of the inertial force contribution to system residual
virtual bool piston_theory_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, MAST::BoundaryConditionBase &bc)=0
Calculates the force vector and Jacobian due to piston-theory based surface pressure on the entire el...
virtual bool linearized_inertial_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xddot, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual of linerized problem
virtual void thermal_residual_boundary_velocity(const MAST::FunctionBase &p, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
Calculates the sensitivity of force vector and Jacobian due to thermal stresses.
virtual void set_perturbed_acceleration(const RealVectorX &vec, bool if_sens=false)
stores vec as perturbed acceleration for element level calculations, or its sensitivity if if_sens is...
virtual void set_acceleration(const RealVectorX &vec, bool if_sens=false)
stores vec as acceleration for element level calculations, or its sensitivity if if_sens is true...
Matrix< Complex, Dynamic, 1 > ComplexVectorX
RealVectorX _local_delta_sol_sens
local perturbed solution sensitivity
virtual unsigned int incompatible_mode_size() const
virtual bool inertial_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xddot, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual
RealVectorX _local_delta_accel_sens
local perturbed acceleration sensitivity
virtual bool damping_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
damping force contribution to system residual
RealVectorX _local_accel
local acceleration
bool linearized_volume_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc)
volume external force contribution to system residual.
virtual bool internal_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
sensitivity of the internal force contribution to system residual
libMesh::Real Real
RealVectorX _local_sol_sens
local solution sensitivity
virtual bool if_incompatible_modes() const =0
void transform_vector_to_local_system(const ValType &global_vec, ValType &local_vec) const
virtual void thermal_residual_temperature_derivative(const MAST::FEBase &fe_thermal, RealMatrixX &m)=0
Real piston_theory_cp(const unsigned int order, const Real vel_U, const Real gamma, const Real mach)
RealVectorX _local_delta_vel
local perturbed velocity
virtual void set_velocity(const RealVectorX &vec, bool if_sens=false)
stores vec as velocity for element level calculations, or its sensitivity if if_sens is true...
RealVectorX _local_delta_sol
local perturbed solution
bool follower_forces
flag for follower forces
const MAST::ElementPropertyCardBase & elem_property()
returns a constant reference to the finite element object
Matrix< Real, Dynamic, Dynamic > RealMatrixX
Matrix< Complex, Dynamic, Dynamic > ComplexMatrixX
virtual bool prestress_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
sensitivity of the prestress force contribution to system residual
Real piston_theory_dcp_dv(const unsigned int order, const Real vel_U, const Real gamma, const Real mach)
virtual bool surface_traction_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0
Calculates the sensitivity of element vector and matrix quantities for surface traction boundary cond...
bool side_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &bc)
side external force contribution to system residual.
virtual void set_solution(const RealVectorX &vec, bool if_sens=false)
stores vec as solution for element level calculations, or its sensitivity if if_sens is true...
virtual bool calculate_stress(bool request_derivative, const MAST::FunctionBase *f, MAST::StressStrainOutputBase &output)=0
Calculates the stress tensor.
Matrix< Real, Dynamic, 1 > RealVectorX
virtual bool linearized_frequency_domain_surface_pressure_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0
Calculates the sensitivity of force vector and Jacobian due to small is applicable for perturbation s...
bool volume_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc)
volume external force contribution to system residual.
std::unique_ptr< MAST::StructuralElementBase > build_structural_element(MAST::SystemInitialization &sys, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p)
builds the structural element for the specified element type
virtual bool piston_theory_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, MAST::BoundaryConditionBase &bc)=0
Calculates the force vector and Jacobian due to piston-theory based surface pressure on the entire el...
const MAST::ElementPropertyCardBase & _property
element property
virtual bool surface_traction_residual_shifted_boundary_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0
Calculates the sensitivity of force vector and Jacobian due to surface traction and sensitiity due to...
StructuralElementBase(MAST::SystemInitialization &sys, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p)
Constructor.
virtual void calculate_stress_boundary_velocity(const MAST::FunctionBase &p, MAST::StressStrainOutputBase &output, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f)=0
Calculates the boundary velocity term contributions to the sensitivity of stress at the specified bou...
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 set_perturbed_velocity(const RealVectorX &vec, bool if_sens=false)
stores vec as perturbed velocity for element level calculations, or its sensitivity if if_sens is tru...
virtual bool surface_traction_residual_shifted_boundary(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0
Calculates the sensitivity of force vector and Jacobian due to surface traction and sensitiity due to...
const RealVectorX & local_solution(bool if_sens=false) const
RealVectorX _local_vel_sens
local velocity sensitivity
virtual bool linearized_surface_pressure_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc)
Calculates the force vector and Jacobian due to surface pressure applied on the entire element domain...
const MAST::GeomElem & elem() const
Definition: elem_base.h:109
virtual void calculate_stress_temperature_derivative(MAST::FEBase &fe_thermal, MAST::StressStrainOutputBase &output)=0
virtual void set_perturbed_solution(const RealVectorX &vec, bool if_sens=false)
stores vec as perturbed solution for element level calculations, or its sensitivity if if_sens is tru...
virtual bool damping_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
sensitivity of the damping force contribution to system residual
bool side_external_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &bc)
sensitivity of the side external force contribution to system residual
virtual bool prestress_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
prestress force contribution to system residual
bool linearized_side_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &bc)
side external force contribution to system residual.
virtual bool internal_residual_jac_dot_state_sensitivity(RealMatrixX &jac)=0
calculates d[J]/d{x} .
void transform_matrix_to_global_system(const ValType &local_mat, ValType &global_mat) const
bool linearized_frequency_domain_side_external_residual(bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &bc)
Calculates the external force due to frequency domain side external force contribution to system resi...
virtual bool linearized_internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual of the linearized problem
RealVectorX * _incompatible_sol
incompatible mode solution vector
RealVectorX _local_delta_accel
local perturbed acceleration
RealVectorX _local_sol
local solution
virtual bool surface_pressure_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0
Calculates the force vector and Jacobian due to surface pressure.
bool volume_external_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc)
sensitivity of the volume external force contribution to system residual
RealVectorX _local_vel
local velocity
void set_incompatible_mode_solution(RealVectorX &vec)
sets the pointer to the incompatible mode solution vector.
virtual bool surface_traction_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0
Calculates the force vector and Jacobian due to surface traction.
virtual bool internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
internal force contribution to system residual
virtual void inertial_residual_boundary_velocity(const MAST::FunctionBase &p, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xddot, RealMatrixX &jac_xdot, RealMatrixX &jac)
sensitivity of the inertial force contribution to system residual
void volume_external_residual_boundary_velocity(const MAST::FunctionBase &p, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
boundary velocity contribution of volume external force.
void transform_vector_to_global_system(const ValType &local_vec, ValType &global_vec) const
virtual bool surface_pressure_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0
Calculates the force vector and Jacobian due to surface pressure.
This is the base class for elements that implement calculation of finite element quantities over the ...
Definition: elem_base.h:72