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-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__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:
51  const MAST::GeomElem& elem,
53 
54  virtual ~StructuralElementBase();
55 
56 
61  virtual void set_solution(const RealVectorX& vec,
62  bool if_sens = false);
63 
64 
69  virtual void set_perturbed_solution(const RealVectorX& vec,
70  bool if_sens = false);
71 
72 
77  virtual void set_velocity(const RealVectorX& vec,
78  bool if_sens = false);
79 
84  virtual void set_perturbed_velocity(const RealVectorX& vec,
85  bool if_sens = false);
86 
87 
92  virtual void set_acceleration(const RealVectorX& vec,
93  bool if_sens = false);
94 
95 
100  virtual void set_perturbed_acceleration(const RealVectorX& vec,
101  bool if_sens = false);
102 
103 
109  const RealVectorX& local_solution(bool if_sens = false) const {
110 
111  if (!if_sens)
112  return _local_sol;
113  else
114  return _local_sol_sens;
115  }
116 
121  return _property;
122  }
123 
124 
128  virtual bool internal_residual (bool request_jacobian,
129  RealVectorX& f,
130  RealMatrixX& jac) = 0;
131 
136  virtual bool
137  linearized_internal_residual (bool request_jacobian,
138  RealVectorX& f,
139  RealMatrixX& jac);
140 
146  virtual void
148  const unsigned int s,
150  bool request_jacobian,
151  RealVectorX& f,
152  RealMatrixX& jac) = 0;
153 
157  virtual bool
159 
160 
161 
165  virtual bool damping_residual (bool request_jacobian,
166  RealVectorX& f,
167  RealMatrixX& jac_xdot,
168  RealMatrixX& jac) {
169 
170  // to be implemented nothing is done yet
171  return false;
172  }
173 
177  virtual bool inertial_residual (bool request_jacobian,
178  RealVectorX& f,
179  RealMatrixX& jac_xddot,
180  RealMatrixX& jac_xdot,
181  RealMatrixX& jac);
182 
187  virtual bool linearized_inertial_residual (bool request_jacobian,
188  RealVectorX& f,
189  RealMatrixX& jac_xddot,
190  RealMatrixX& jac_xdot,
191  RealMatrixX& jac);
192 
193 
200  bool side_external_residual (bool request_jacobian,
201  RealVectorX& f,
202  RealMatrixX& jac_xdot,
203  RealMatrixX& jac,
204  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
205 
212  bool
214  (bool request_jacobian,
215  RealVectorX& f,
216  RealMatrixX& jac_xdot,
217  RealMatrixX& jac,
218  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
219 
220 
221 
228  bool
230  (bool request_jacobian,
231  ComplexVectorX& f,
232  ComplexMatrixX& jac,
233  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
234 
235 
239  virtual bool prestress_residual (bool request_jacobian,
240  RealVectorX& f,
241  RealMatrixX& jac) = 0;
242 
249  bool volume_external_residual (bool request_jacobian,
250  RealVectorX& f,
251  RealMatrixX& jac_xdot,
252  RealMatrixX& jac,
253  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
254 
255 
260  (const MAST::FunctionBase& p,
261  const unsigned int s,
263  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc,
264  bool request_jacobian,
265  RealVectorX& f,
266  RealMatrixX& jac);
267 
268 
275  bool
277  (bool request_jacobian,
278  RealVectorX& f,
279  RealMatrixX& jac_xdot,
280  RealMatrixX& jac,
281  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
282 
283 
289  bool
291  (bool request_jacobian,
292  ComplexVectorX& f,
293  ComplexMatrixX& jac,
294  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
295 
296 
297 
301  virtual bool internal_residual_sensitivity (const MAST::FunctionBase& p,
302  bool request_jacobian,
303  RealVectorX& f,
304  RealMatrixX& jac) = 0;
305 
310  bool request_jacobian,
311  RealVectorX& f,
312  RealMatrixX& jac_xdot,
313  RealMatrixX& jac) {
314  return false;
315  }
316 
320  virtual bool inertial_residual_sensitivity (const MAST::FunctionBase& p,
321  bool request_jacobian,
322  RealVectorX& f,
323  RealMatrixX& jac_xddot,
324  RealMatrixX& jac_xdot,
325  RealMatrixX& jac);
326 
330  virtual void
332  const unsigned int s,
334  bool request_jacobian,
335  RealVectorX& f,
336  RealMatrixX& jac_xddot,
337  RealMatrixX& jac_xdot,
338  RealMatrixX& jac);
339 
344  bool request_jacobian,
345  RealVectorX& f,
346  RealMatrixX& jac_xdot,
347  RealMatrixX& jac,
348  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
349 
353  virtual bool prestress_residual_sensitivity (const MAST::FunctionBase& p,
354  bool request_jacobian,
355  RealVectorX& f,
356  RealMatrixX& jac) = 0;
357 
362  bool request_jacobian,
363  RealVectorX& f,
364  RealMatrixX& jac_xdot,
365  RealMatrixX& jac,
366  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
367 
368 
369 
374  virtual bool if_incompatible_modes() const = 0;
375 
376 
380  virtual unsigned int incompatible_mode_size() const {
381  // needs to be implemented only in inherited elements
382  // that use incompatible modes
383  libmesh_error();
384  return 0;
385  }
386 
387 
394  _incompatible_sol = &vec;
395  }
396 
397 
403  virtual void update_incompatible_mode_solution(const RealVectorX& dsol) {
404  // should be implemented in derived classes, if relevant
405  libmesh_error();
406  }
407 
408 
414  virtual bool calculate_stress(bool request_derivative,
415  const MAST::FunctionBase* f,
416  MAST::StressStrainOutputBase& output) = 0;
417 
418 
423  virtual void
426  const unsigned int s,
427  const MAST::FieldFunction<RealVectorX>& vel_f) = 0;
428 
429  virtual void
431  MAST::StressStrainOutputBase& output) = 0;
432 
433  virtual void
435  RealMatrixX& m) = 0;
436 
441 
442 
443  template <typename ValType>
444  void transform_matrix_to_global_system(const ValType& local_mat,
445  ValType& global_mat) const;
446 
447  template <typename ValType>
448  void transform_vector_to_local_system(const ValType& global_vec,
449  ValType& local_vec) const;
450 
451  template <typename ValType>
452  void transform_vector_to_global_system(const ValType& local_vec,
453  ValType& global_vec) const;
454 
455  protected:
456 
460  virtual bool
461  surface_pressure_residual(bool request_jacobian,
462  RealVectorX& f,
463  RealMatrixX& jac,
464  const unsigned int side,
466 
467 
473  virtual bool
474  surface_pressure_residual(bool request_jacobian,
475  RealVectorX& f,
476  RealMatrixX& jac,
478 
484  virtual void
486  const unsigned int s,
489  bool request_jacobian,
490  RealVectorX& f,
491  RealMatrixX& jac);
492 
498  virtual bool
499  linearized_surface_pressure_residual(bool request_jacobian,
500  RealVectorX& f,
501  RealMatrixX& jac,
503 
504 
509  Real piston_theory_cp(const unsigned int order,
510  const Real vel_U,
511  const Real gamma,
512  const Real mach);
513 
514 
515 
520  Real piston_theory_dcp_dv(const unsigned int order,
521  const Real vel_U,
522  const Real gamma,
523  const Real mach);
524 
525 
533  virtual bool piston_theory_residual(bool request_jacobian,
534  RealVectorX &f,
535  RealMatrixX& jac_xdot,
536  RealMatrixX& jac,
538 
545  virtual bool piston_theory_residual(bool request_jacobian,
546  RealVectorX &f,
547  RealMatrixX& jac_xdot,
548  RealMatrixX& jac,
549  const unsigned int side,
551 
552 
553 
562  bool request_jacobian,
563  RealVectorX &f,
564  RealMatrixX& jac_xdot,
565  RealMatrixX& jac,
567 
575  bool request_jacobian,
576  RealVectorX &f,
577  RealMatrixX& jac_xdot,
578  RealMatrixX& jac,
579  const unsigned int side,
581 
582 
587  virtual bool
589  (bool request_jacobian,
590  ComplexVectorX& f,
591  ComplexMatrixX& jac,
592  const unsigned int side,
594 
595 
604  virtual bool
606  (bool request_jacobian,
607  ComplexVectorX& f,
608  ComplexMatrixX& jac,
610 
611 
616  virtual bool
618  (const MAST::FunctionBase& p,
619  bool request_jacobian,
620  ComplexVectorX& f,
621  ComplexMatrixX& jac,
622  const unsigned int side,
624 
625 
631  virtual bool
634  bool request_jacobian,
635  ComplexVectorX& f,
636  ComplexMatrixX& jac,
638 
639  libmesh_error(); // to be implemented
640  }
641 
642 
646  virtual bool
648  bool request_jacobian,
649  RealVectorX& f,
650  RealMatrixX& jac,
651  const unsigned int side,
653 
658  virtual bool
660  bool request_jacobian,
661  RealVectorX& f,
662  RealMatrixX& jac,
664 
665 
666 
671  virtual bool thermal_residual(bool request_jacobian,
672  RealVectorX& f,
673  RealMatrixX& jac,
675 
680  virtual bool thermal_residual_sensitivity(const MAST::FunctionBase& p,
681  bool request_jacobian,
682  RealVectorX& f,
683  RealMatrixX& jac,
685 
691  const unsigned int s,
694  bool request_jacobian,
695  RealVectorX& f,
696  RealMatrixX& jac) = 0;
697 
702 
703 
708 
709 
714 
715 
720 
721 
726 
727 
732 
733 
738 
739 
744 
745 
750 
751 
756 
757 
762 
763 
768 
769 
774 
775 
780 
781  };
782 
783 
787  std::unique_ptr<MAST::StructuralElementBase>
789  MAST::AssemblyBase& assembly,
790  const MAST::GeomElem& elem,
792 
793 }
794 
795 
796 #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: .
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.
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.
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
std::unique_ptr< MAST::StructuralElementBase > build_structural_element(MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p)
builds the structural element for the specified element type
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...
const RealVectorX & local_solution(bool if_sens=false) const
const MAST::GeomElem & elem() const
Definition: elem_base.h:117
Matrix< Complex, Dynamic, 1 > ComplexVectorX
void transform_vector_to_local_system(const ValType &global_vec, ValType &local_vec) const
RealVectorX _local_delta_sol_sens
local perturbed solution sensitivity
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 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_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
MAST::AssemblyBase & assembly()
Definition: elem_base.h:103
virtual bool if_incompatible_modes() const =0
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...
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
void transform_matrix_to_global_system(const ValType &local_mat, ValType &global_mat) const
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)
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.
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...
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...
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.
const MAST::ElementPropertyCardBase & _property
element property
virtual unsigned int incompatible_mode_size() const
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...
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...
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 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...
StructuralElementBase(MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p)
Constructor.
virtual bool prestress_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
prestress force contribution to system residual
virtual bool internal_residual_jac_dot_state_sensitivity(RealMatrixX &jac)=0
calculates d[J]/d{x} .
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
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.
RealVectorX _local_vel
local velocity
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.
void set_incompatible_mode_solution(RealVectorX &vec)
sets the pointer to the incompatible mode solution vector.
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 transform_vector_to_global_system(const ValType &local_vec, ValType &global_vec) const
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 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