libMesh
Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
SigmaPhysics Class Reference

#include <sigma_physics.h>

Inheritance diagram for SigmaPhysics:
[legend]

Public Member Functions

 SigmaPhysics ()
 
 SigmaPhysics (SigmaPhysics &&)=default
 
 SigmaPhysics (const SigmaPhysics &)=default
 
SigmaPhysicsoperator= (const SigmaPhysics &)=default
 
SigmaPhysicsoperator= (SigmaPhysics &&)=default
 
virtual ~SigmaPhysics ()=default
 
Real & k ()
 
bool & analytic_jacobians ()
 
virtual void init_data (System &sys)
 
virtual std::unique_ptr< DifferentiablePhysicsclone_physics () override
 Copy of this object. More...
 
virtual bool eulerian_residual (bool request_jacobian, DiffContext &context) override
 Adds a pseudo-convection contribution on elem to elem_residual, if the nodes of elem are being translated by a moving mesh. More...
 
virtual bool mass_residual (bool request_jacobian, DiffContext &) override
 Subtracts a mass vector contribution on elem from elem_residual. More...
 
virtual void clear_physics ()
 Clear any data structures associated with the physics. More...
 
virtual void init_physics (const System &sys)
 Initialize any data structures associated with the physics. More...
 
virtual bool element_constraint (bool request_jacobian, DiffContext &)
 Adds the constraint contribution on elem to elem_residual. More...
 
virtual bool side_time_derivative (bool request_jacobian, DiffContext &)
 Adds the time derivative contribution on side of elem to elem_residual. More...
 
virtual bool side_constraint (bool request_jacobian, DiffContext &)
 Adds the constraint contribution on side of elem to elem_residual. More...
 
virtual bool nonlocal_time_derivative (bool request_jacobian, DiffContext &)
 Adds any nonlocal time derivative contributions (e.g. More...
 
virtual bool nonlocal_constraint (bool request_jacobian, DiffContext &)
 Adds any nonlocal constraint contributions (e.g. More...
 
virtual void time_evolving (unsigned int var, unsigned int order)
 Tells the DiffSystem that variable var is evolving with respect to time. More...
 
bool is_time_evolving (unsigned int var) const
 
virtual bool side_mass_residual (bool request_jacobian, DiffContext &)
 Subtracts a mass vector contribution on side of elem from elem_residual. More...
 
virtual bool nonlocal_mass_residual (bool request_jacobian, DiffContext &c)
 Subtracts any nonlocal mass vector contributions (e.g. More...
 
virtual bool damping_residual (bool request_jacobian, DiffContext &)
 Subtracts a damping vector contribution on elem from elem_residual. More...
 
virtual bool side_damping_residual (bool request_jacobian, DiffContext &)
 Subtracts a damping vector contribution on side of elem from elem_residual. More...
 
virtual bool nonlocal_damping_residual (bool request_jacobian, DiffContext &)
 Subtracts any nonlocal damping vector contributions (e.g. More...
 
virtual void set_mesh_system (System *sys)
 Tells the DifferentiablePhysics that system sys contains the isoparametric Lagrangian variables which correspond to the coordinates of mesh nodes, in problems where the mesh itself is expected to move in time. More...
 
const System * get_mesh_system () const
 
System * get_mesh_system ()
 
virtual void set_mesh_x_var (unsigned int var)
 Tells the DifferentiablePhysics that variable var from the mesh system should be used to update the x coordinate of mesh nodes, in problems where the mesh itself is expected to move in time. More...
 
unsigned int get_mesh_x_var () const
 
virtual void set_mesh_y_var (unsigned int var)
 Tells the DifferentiablePhysics that variable var from the mesh system should be used to update the y coordinate of mesh nodes. More...
 
unsigned int get_mesh_y_var () const
 
virtual void set_mesh_z_var (unsigned int var)
 Tells the DifferentiablePhysics that variable var from the mesh system should be used to update the z coordinate of mesh nodes. More...
 
unsigned int get_mesh_z_var () const
 
bool _eulerian_time_deriv (bool request_jacobian, DiffContext &)
 This method simply combines element_time_derivative() and eulerian_residual(), which makes its address useful as a pointer-to-member-function when refactoring. More...
 
bool have_first_order_vars () const
 
const std::set< unsigned int > & get_first_order_vars () const
 
bool is_first_order_var (unsigned int var) const
 
bool have_second_order_vars () const
 
const std::set< unsigned int > & get_second_order_vars () const
 
bool is_second_order_var (unsigned int var) const
 

Public Attributes

bool compute_internal_sides
 compute_internal_sides is false by default, indicating that side_* computations will only be done on boundary sides. More...
 

Protected Member Functions

virtual void init_context (DiffContext &context) override
 
virtual bool element_time_derivative (bool request_jacobian, DiffContext &context) override
 Adds the time derivative contribution on elem to elem_residual. More...
 

Protected Attributes

Real _k
 
Number computed_QoI [1]
 
std::string _fe_family
 
unsigned int _fe_order
 
unsigned int T_var
 
bool _analytic_jacobians
 
System * _mesh_sys
 System from which to acquire moving mesh information. More...
 
unsigned int _mesh_x_var
 Variables from which to acquire moving mesh information. More...
 
unsigned int _mesh_y_var
 
unsigned int _mesh_z_var
 
std::vector< unsigned int_time_evolving
 Stores unsigned int to tell us which variables are evolving as first order in time (1), second order in time (2), or are not time evolving (0). More...
 
std::set< unsigned int_first_order_vars
 Variable indices for those variables that are first order in time. More...
 
std::set< unsigned int_second_order_vars
 Variable indices for those variables that are second order in time. More...
 
std::map< unsigned int, unsigned int_second_order_dot_vars
 If the user adds any second order variables, then we need to also cache the map to their corresponding dot variable that will be added by this TimeSolver class. More...
 

Detailed Description

Definition at line 12 of file sigma_physics.h.

Constructor & Destructor Documentation

◆ SigmaPhysics() [1/3]

SigmaPhysics::SigmaPhysics ( )
inline

Definition at line 16 of file sigma_physics.h.

16 : FEMPhysics() {}
FEMPhysics()
Constructor.
Definition: fem_physics.h:51

◆ SigmaPhysics() [2/3]

SigmaPhysics::SigmaPhysics ( SigmaPhysics &&  )
default

◆ SigmaPhysics() [3/3]

SigmaPhysics::SigmaPhysics ( const SigmaPhysics )
default

◆ ~SigmaPhysics()

virtual SigmaPhysics::~SigmaPhysics ( )
virtualdefault

Member Function Documentation

◆ _eulerian_time_deriv()

bool libMesh::DifferentiablePhysics::_eulerian_time_deriv ( bool  request_jacobian,
DiffContext context 
)
inherited

This method simply combines element_time_derivative() and eulerian_residual(), which makes its address useful as a pointer-to-member-function when refactoring.

Definition at line 96 of file diff_physics.C.

References libMesh::DifferentiablePhysics::element_time_derivative(), and libMesh::DifferentiablePhysics::eulerian_residual().

Referenced by libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), and libMesh::NewmarkSolver::element_residual().

98 {
99  // For any problem we need time derivative terms
100  request_jacobian =
101  this->element_time_derivative(request_jacobian, context);
102 
103  // For a moving mesh problem we may need the pseudoconvection term too
104  return this->eulerian_residual(request_jacobian, context) &&
105  request_jacobian;
106 }
virtual bool element_time_derivative(bool request_jacobian, DiffContext &)
Adds the time derivative contribution on elem to elem_residual.
Definition: diff_physics.h:125
virtual bool eulerian_residual(bool request_jacobian, DiffContext &)
Adds a pseudo-convection contribution on elem to elem_residual, if the nodes of elem are being transl...
Definition: diff_physics.h:277

◆ analytic_jacobians()

bool& SigmaPhysics::analytic_jacobians ( )
inline

Definition at line 26 of file sigma_physics.h.

26 { return _analytic_jacobians; }
bool _analytic_jacobians
Definition: sigma_physics.h:57

◆ clear_physics()

void libMesh::DifferentiablePhysics::clear_physics ( )
virtualinherited

Clear any data structures associated with the physics.

Definition at line 28 of file diff_physics.C.

References libMesh::DifferentiablePhysics::_time_evolving.

Referenced by libMesh::DifferentiableSystem::clear().

29 {
30  _time_evolving.resize(0);
31 }
std::vector< unsigned int > _time_evolving
Stores unsigned int to tell us which variables are evolving as first order in time (1)...
Definition: diff_physics.h:543

◆ clone_physics()

std::unique_ptr< DifferentiablePhysics > SigmaPhysics::clone_physics ( )
overridevirtual

Copy of this object.

User should override to copy any needed state.

Implements libMesh::DifferentiablePhysics.

Definition at line 174 of file sigma_physics.C.

175 {
176  return std::make_unique<SigmaPhysics>(*this);
177 }

◆ damping_residual()

virtual bool libMesh::DifferentiablePhysics::damping_residual ( bool  request_jacobian,
DiffContext  
)
inlinevirtualinherited

Subtracts a damping vector contribution on elem from elem_residual.

This method is not used in first-order-in-time problems. For second-order-in-time problems, this is the \( C(u,\ddot{u})\ddot{u} \) term. This method is only called for UnsteadySolver-based TimeSolvers.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

If the problem has no damping, the default "do-nothing" is correct. Otherwise, this must be reimplemented.

Reimplemented in SecondOrderScalarSystemFirstOrderTimeSolverBase, and SecondOrderScalarSystemSecondOrderTimeSolverBase.

Definition at line 360 of file diff_physics.h.

Referenced by libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), and libMesh::NewmarkSolver::element_residual().

362  {
363  return request_jacobian;
364  }

◆ element_constraint()

virtual bool libMesh::DifferentiablePhysics::element_constraint ( bool  request_jacobian,
DiffContext  
)
inlinevirtualinherited

Adds the constraint contribution on elem to elem_residual.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

Users may need to reimplement this for their particular PDE.

To implement the constraint 0 = G(u), the user should examine u = elem_solution and add (G(u), phi_i) to elem_residual in elem_constraint().

Reimplemented in CoupledSystem, and NavierSystem.

Definition at line 144 of file diff_physics.h.

Referenced by libMesh::EulerSolver::element_residual(), libMesh::SteadySolver::element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::NewmarkSolver::element_residual(), and libMesh::EigenTimeSolver::element_residual().

146  {
147  return request_jacobian;
148  }

◆ element_time_derivative()

bool SigmaPhysics::element_time_derivative ( bool  request_jacobian,
DiffContext &   
)
overrideprotectedvirtual

Adds the time derivative contribution on elem to elem_residual.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

Users need to reimplement this for their particular PDE.

To implement the physics model du/dt = F(u), the user should examine u = elem_solution and add (F(u), phi_i) to elem_residual in elem_time_derivative().

Reimplemented from libMesh::DifferentiablePhysics.

Definition at line 106 of file sigma_physics.C.

References _analytic_jacobians, compute_jacobian(), libMesh::Real, and T_var.

108 {
109  bool compute_jacobian = request_jacobian && _analytic_jacobians;
110 
111  FEMContext & c = cast_ref<FEMContext &>(context);
112 
113  // First we get some references to cell-specific data that
114  // will be used to assemble the linear system.
115  FEBase * elem_fe = nullptr;
116  c.get_element_fe(T_var, elem_fe);
117 
118  // Element Jacobian * quadrature weights for interior integration
119  const std::vector<Real> & JxW = elem_fe->get_JxW();
120 
121  // Element basis functions
122  const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
123  const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
124 
125  // The number of local degrees of freedom in each variable
126  const unsigned int n_u_dofs = c.n_dof_indices(T_var);
127 
128  // The subvectors and submatrices we need to fill:
129  DenseSubMatrix<Number> & K = c.get_elem_jacobian(T_var, T_var);
130  DenseSubVector<Number> & F = c.get_elem_residual(T_var);
131 
132  // Quadrature point locations
133  const std::vector<Point > & q_point = elem_fe->get_xyz();
134 
135  // Now we will build the element Jacobian and residual.
136  // Constructing the residual requires the solution and its
137  // gradient from the previous timestep. This must be
138  // calculated at each quadrature point by summing the
139  // solution degree-of-freedom values by the appropriate
140  // weight functions.
141  unsigned int n_qpoints = c.get_element_qrule().n_points();
142 
143  // Conductivity
144  Real sigma = 1.0;
145 
146  // Forcing function
147  Real f = 1.0;
148 
149  for (unsigned int qp=0; qp != n_qpoints; qp++)
150  {
151  // Compute the solution gradient at the Newton iterate
152  Gradient grad_T = c.interior_gradient(T_var, qp);
153 
154  // Location of the current qp
155  const Real x = q_point[qp](0);
156 
157  // Spatially varying conductivity
158  sigma = 0.001 + x;
159 
160  for (unsigned int i=0; i != n_u_dofs; i++)
161  {
162 F(i) += JxW[qp] * ( ( -sigma * (grad_T * dphi[i][qp]) ) + (f * phi[i][qp]) );
163  }
164 
165  if (compute_jacobian)
166  for (unsigned int i=0; i != n_u_dofs; i++)
167  for (unsigned int j=0; j != n_u_dofs; ++j)
168  K(i,j) += JxW[qp] * -sigma * (dphi[i][qp] * dphi[j][qp]);
169  } // end of the quadrature point qp-loop
170 
171  return compute_jacobian;
172 }
bool _analytic_jacobians
Definition: sigma_physics.h:57
void compute_jacobian(const NumericVector< Number > &, SparseMatrix< Number > &J, NonlinearImplicitSystem &system)
Definition: assembly.C:315
NumberVectorValue Gradient
FEGenericBase< Real > FEBase
unsigned int T_var
Definition: sigma_physics.h:54
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ eulerian_residual()

bool libMesh::FEMPhysics::eulerian_residual ( bool  request_jacobian,
DiffContext context 
)
overridevirtualinherited

Adds a pseudo-convection contribution on elem to elem_residual, if the nodes of elem are being translated by a moving mesh.

This function assumes that the user's time derivative equations (except for any equations involving unknown mesh xyz coordinates themselves) are expressed in an Eulerian frame of reference, and that the user is satisfied with an unstabilized convection term. Lagrangian equations will probably require overriding eulerian_residual() with a blank function; ALE or stabilized formulations will require reimplementing eulerian_residual() entirely.

Reimplemented from libMesh::DifferentiablePhysics.

Reimplemented in SolidSystem.

Definition at line 38 of file fem_physics.C.

References libMesh::DifferentiablePhysics::_mesh_sys, libMesh::DifferentiablePhysics::_mesh_x_var, libMesh::DifferentiablePhysics::_mesh_y_var, libMesh::DifferentiablePhysics::_mesh_z_var, libMesh::FEMContext::get_element_qrule(), libMesh::FEMContext::interior_gradient(), libMesh::invalid_uint, libMesh::DifferentiablePhysics::is_time_evolving(), libMesh::libmesh_assert(), libMesh::libmesh_real(), libMesh::make_range(), libMesh::DiffContext::n_vars(), and libMesh::UnsteadySolver::old_nonlinear_solution().

40 {
41  // Only calculate a mesh movement residual if it's necessary
42  if (!_mesh_sys)
43  return request_jacobian;
44 
45  libmesh_not_implemented();
46 
47 #if 0
48  FEMContext & context = cast_ref<FEMContext &>(c);
49 
50  // This function only supports fully coupled mesh motion for now
51  libmesh_assert_equal_to (_mesh_sys, this);
52 
53  unsigned int n_qpoints = (context.get_element_qrule())->n_points();
54 
55  const unsigned int n_x_dofs = (_mesh_x_var == libMesh::invalid_uint) ?
56  0 : context.dof_indices_var[_mesh_x_var].size();
57  const unsigned int n_y_dofs = (_mesh_y_var == libMesh::invalid_uint) ?
58  0 : context.dof_indices_var[_mesh_y_var].size();
59  const unsigned int n_z_dofs = (_mesh_z_var == libMesh::invalid_uint) ?
60  0 : context.dof_indices_var[_mesh_z_var].size();
61 
62  const unsigned int mesh_xyz_var = n_x_dofs ? _mesh_x_var :
63  (n_y_dofs ? _mesh_y_var :
64  (n_z_dofs ? _mesh_z_var :
66 
67  // If we're our own _mesh_sys, we'd better be in charge of
68  // at least one coordinate, and we'd better have the same
69  // FE type for all coordinates we are in charge of
70  libmesh_assert_not_equal_to (mesh_xyz_var, libMesh::invalid_uint);
71  libmesh_assert(!n_x_dofs || context.element_fe_var[_mesh_x_var] ==
72  context.element_fe_var[mesh_xyz_var]);
73  libmesh_assert(!n_y_dofs || context.element_fe_var[_mesh_y_var] ==
74  context.element_fe_var[mesh_xyz_var]);
75  libmesh_assert(!n_z_dofs || context.element_fe_var[_mesh_z_var] ==
76  context.element_fe_var[mesh_xyz_var]);
77 
78  const std::vector<std::vector<Real>> & psi =
79  context.element_fe_var[mesh_xyz_var]->get_phi();
80 
81  for (auto var : make_range(context.n_vars()))
82  {
83  // Mesh motion only affects time-evolving variables
84  if (this->is_time_evolving(var))
85  continue;
86 
87  // The mesh coordinate variables themselves are Lagrangian,
88  // not Eulerian, and no convective term is desired.
89  if (/*_mesh_sys == this && */
90  (var == _mesh_x_var ||
91  var == _mesh_y_var ||
92  var == _mesh_z_var))
93  continue;
94 
95  // Some of this code currently relies on the assumption that
96  // we can pull mesh coordinate data from our own system
97  if (_mesh_sys != this)
98  libmesh_not_implemented();
99 
100  // This residual should only be called by unsteady solvers:
101  // if the mesh is steady, there's no mesh convection term!
102  UnsteadySolver * unsteady;
103  if (this->time_solver->is_steady())
104  return request_jacobian;
105  else
106  unsteady = cast_ptr<UnsteadySolver*>(this->time_solver.get());
107 
108  const std::vector<Real> & JxW =
109  context.element_fe_var[var]->get_JxW();
110 
111  const std::vector<std::vector<Real>> & phi =
112  context.element_fe_var[var]->get_phi();
113 
114  const std::vector<std::vector<RealGradient>> & dphi =
115  context.element_fe_var[var]->get_dphi();
116 
117  const unsigned int n_u_dofs = context.dof_indices_var[var].size();
118 
119  DenseSubVector<Number> & Fu = *context.elem_subresiduals[var];
120  DenseSubMatrix<Number> & Kuu = *context.elem_subjacobians[var][var];
121 
122  DenseSubMatrix<Number> * Kux = n_x_dofs ?
123  context.elem_subjacobians[var][_mesh_x_var] : nullptr;
124  DenseSubMatrix<Number> * Kuy = n_y_dofs ?
125  context.elem_subjacobians[var][_mesh_y_var] : nullptr;
126  DenseSubMatrix<Number> * Kuz = n_z_dofs ?
127  context.elem_subjacobians[var][_mesh_z_var] : nullptr;
128 
129  std::vector<Real> delta_x(n_x_dofs, 0.);
130  std::vector<Real> delta_y(n_y_dofs, 0.);
131  std::vector<Real> delta_z(n_z_dofs, 0.);
132 
133  for (unsigned int i = 0; i != n_x_dofs; ++i)
134  {
135  unsigned int j = context.dof_indices_var[_mesh_x_var][i];
136  delta_x[i] = libmesh_real(this->current_solution(j)) -
137  libmesh_real(unsteady->old_nonlinear_solution(j));
138  }
139 
140  for (unsigned int i = 0; i != n_y_dofs; ++i)
141  {
142  unsigned int j = context.dof_indices_var[_mesh_y_var][i];
143  delta_y[i] = libmesh_real(this->current_solution(j)) -
144  libmesh_real(unsteady->old_nonlinear_solution(j));
145  }
146 
147  for (unsigned int i = 0; i != n_z_dofs; ++i)
148  {
149  unsigned int j = context.dof_indices_var[_mesh_z_var][i];
150  delta_z[i] = libmesh_real(this->current_solution(j)) -
151  libmesh_real(unsteady->old_nonlinear_solution(j));
152  }
153 
154  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
155  {
156  Gradient grad_u = context.interior_gradient(var, qp);
157  RealGradient convection(0.);
158 
159  for (unsigned int i = 0; i != n_x_dofs; ++i)
160  convection(0) += delta_x[i] * psi[i][qp];
161  for (unsigned int i = 0; i != n_y_dofs; ++i)
162  convection(1) += delta_y[i] * psi[i][qp];
163  for (unsigned int i = 0; i != n_z_dofs; ++i)
164  convection(2) += delta_z[i] * psi[i][qp];
165 
166  for (unsigned int i = 0; i != n_u_dofs; ++i)
167  {
168  Number JxWxPhiI = JxW[qp] * phi[i][qp];
169  Fu(i) += (convection * grad_u) * JxWxPhiI;
170  if (request_jacobian)
171  {
172  Number JxWxPhiI = JxW[qp] * phi[i][qp];
173  for (unsigned int j = 0; j != n_u_dofs; ++j)
174  Kuu(i,j) += JxWxPhiI * (convection * dphi[j][qp]);
175 
176  Number JxWxPhiIoverDT = JxWxPhiI/this->deltat;
177 
178  Number JxWxPhiIxDUDXoverDT = JxWxPhiIoverDT * grad_u(0);
179  for (unsigned int j = 0; j != n_x_dofs; ++j)
180  (*Kux)(i,j) += JxWxPhiIxDUDXoverDT * psi[j][qp];
181 
182  Number JxWxPhiIxDUDYoverDT = JxWxPhiIoverDT * grad_u(1);
183  for (unsigned int j = 0; j != n_y_dofs; ++j)
184  (*Kuy)(i,j) += JxWxPhiIxDUDYoverDT * psi[j][qp];
185 
186  Number JxWxPhiIxDUDZoverDT = JxWxPhiIoverDT * grad_u(2);
187  for (unsigned int j = 0; j != n_z_dofs; ++j)
188  (*Kuz)(i,j) += JxWxPhiIxDUDZoverDT * psi[j][qp];
189  }
190  }
191  }
192  }
193 #endif // 0
194 
195  return request_jacobian;
196 }
T libmesh_real(T a)
RealVectorValue RealGradient
unsigned int _mesh_x_var
Variables from which to acquire moving mesh information.
Definition: diff_physics.h:536
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
Definition: libmesh.h:310
System * _mesh_sys
System from which to acquire moving mesh information.
Definition: diff_physics.h:531
NumberVectorValue Gradient
libmesh_assert(ctx)
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
bool is_time_evolving(unsigned int var) const
Definition: diff_physics.h:260

◆ get_first_order_vars()

const std::set<unsigned int>& libMesh::DifferentiablePhysics::get_first_order_vars ( ) const
inlineinherited
Returns
The set of first order in time variable indices. May be empty.

Definition at line 506 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_first_order_vars.

Referenced by libMesh::DifferentiableSystem::have_first_order_scalar_vars().

507  { return _first_order_vars; }
std::set< unsigned int > _first_order_vars
Variable indices for those variables that are first order in time.
Definition: diff_physics.h:548

◆ get_mesh_system() [1/2]

const System * libMesh::DifferentiablePhysics::get_mesh_system ( ) const
inlineinherited
Returns
A const reference to the system with variables corresponding to mesh nodal coordinates, or nullptr if the mesh is fixed. Useful for ALE calculations.

Definition at line 613 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_sys.

Referenced by libMesh::FEMSystem::build_context().

614 {
615  return _mesh_sys;
616 }
System * _mesh_sys
System from which to acquire moving mesh information.
Definition: diff_physics.h:531

◆ get_mesh_system() [2/2]

System * libMesh::DifferentiablePhysics::get_mesh_system ( )
inlineinherited
Returns
A reference to the system with variables corresponding to mesh nodal coordinates, or nullptr if the mesh is fixed.

Definition at line 619 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_sys.

620 {
621  return _mesh_sys;
622 }
System * _mesh_sys
System from which to acquire moving mesh information.
Definition: diff_physics.h:531

◆ get_mesh_x_var()

unsigned int libMesh::DifferentiablePhysics::get_mesh_x_var ( ) const
inlineinherited
Returns
The variable number corresponding to the mesh x coordinate. Useful for ALE calculations.

Definition at line 625 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_x_var.

Referenced by libMesh::FEMSystem::build_context().

626 {
627  return _mesh_x_var;
628 }
unsigned int _mesh_x_var
Variables from which to acquire moving mesh information.
Definition: diff_physics.h:536

◆ get_mesh_y_var()

unsigned int libMesh::DifferentiablePhysics::get_mesh_y_var ( ) const
inlineinherited
Returns
The variable number corresponding to the mesh y coordinate. Useful for ALE calculations.

Definition at line 631 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_y_var.

Referenced by libMesh::FEMSystem::build_context().

632 {
633  return _mesh_y_var;
634 }

◆ get_mesh_z_var()

unsigned int libMesh::DifferentiablePhysics::get_mesh_z_var ( ) const
inlineinherited
Returns
The variable number corresponding to the mesh z coordinate. Useful for ALE calculations.

Definition at line 637 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_z_var.

Referenced by libMesh::FEMSystem::build_context().

638 {
639  return _mesh_z_var;
640 }

◆ get_second_order_vars()

const std::set<unsigned int>& libMesh::DifferentiablePhysics::get_second_order_vars ( ) const
inlineinherited
Returns
The set of second order in time variable indices. May be empty.

Definition at line 519 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_second_order_vars.

Referenced by libMesh::DifferentiableSystem::add_second_order_dot_vars(), libMesh::DiffContext::DiffContext(), libMesh::Euler2Solver::element_residual(), libMesh::DifferentiableSystem::have_second_order_scalar_vars(), and libMesh::FEMContext::pre_fe_reinit().

520  { return _second_order_vars; }
std::set< unsigned int > _second_order_vars
Variable indices for those variables that are second order in time.
Definition: diff_physics.h:553

◆ have_first_order_vars()

bool libMesh::DifferentiablePhysics::have_first_order_vars ( ) const
inlineinherited

Definition at line 500 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_first_order_vars.

Referenced by libMesh::DifferentiableSystem::have_first_order_scalar_vars().

501  { return !_first_order_vars.empty(); }
std::set< unsigned int > _first_order_vars
Variable indices for those variables that are first order in time.
Definition: diff_physics.h:548

◆ have_second_order_vars()

bool libMesh::DifferentiablePhysics::have_second_order_vars ( ) const
inlineinherited

Definition at line 513 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_second_order_vars.

Referenced by libMesh::EulerSolver::element_residual(), and libMesh::DifferentiableSystem::have_second_order_scalar_vars().

514  { return !_second_order_vars.empty(); }
std::set< unsigned int > _second_order_vars
Variable indices for those variables that are second order in time.
Definition: diff_physics.h:553

◆ init_context()

void SigmaPhysics::init_context ( DiffContext &  context)
overrideprotectedvirtual

Reimplemented from libMesh::DifferentiablePhysics.

Definition at line 63 of file sigma_physics.C.

References libMesh::NumericVector< Number >, and T_var.

64 {
65  FEMContext & c = cast_ref<FEMContext &>(context);
66 
67  FEBase * elem_fe = nullptr;
68  c.get_element_fe(T_var, elem_fe);
69 
70  // Now make sure we have requested all the data
71  // we need to build the linear system.
72  elem_fe->get_JxW();
73  elem_fe->get_dphi();
74  elem_fe->get_phi();
75  elem_fe->get_xyz();
76 
77  // Don't waste time on side computations for T
78  FEBase * side_fe = nullptr;
79  c.get_side_fe(T_var, side_fe);
80  side_fe->get_nothing();
81 
82  // We'll have a more automatic solution to preparing adjoint
83  // solutions for time integration, eventually...
84  if (c.is_adjoint())
85  {
86  // A reference to the system context is built with
87  const System & sys = c.get_system();
88 
89  // Get a pointer to the adjoint solution vector
90  NumericVector<Number> & adjoint_solution0 =
91  const_cast<System &>(sys).get_adjoint_solution(0);
92 
93  // Add this adjoint solution to the vectors that diff context should localize
94  c.add_localized_vector(adjoint_solution0, sys);
95 
96  // Get a pointer to the adjoint solution vector
97  NumericVector<Number> & adjoint_solution1 =
98  const_cast<System &>(sys).get_adjoint_solution(1);
99 
100  // Add this adjoint solution to the vectors that diff context should localize
101  c.add_localized_vector(adjoint_solution1, sys);
102  }
103 
104 }
FEGenericBase< Real > FEBase
unsigned int T_var
Definition: sigma_physics.h:54
template class LIBMESH_EXPORT NumericVector< Number >

◆ init_data()

void SigmaPhysics::init_data ( System &  sys)
virtual

Definition at line 39 of file sigma_physics.C.

References _analytic_jacobians, _k, libMesh::LAGRANGE, T_var, libMesh::DifferentiablePhysics::time_evolving(), and libMesh::zero.

40 {
41  T_var = dynamic_cast<FEMSystem&>(sys).add_variable ("T", static_cast<Order>(1), LAGRANGE);
42 
43  GetPot infile("heat.in");
44  _k = infile("k", 1.0);
45  _analytic_jacobians = infile("analytic_jacobians", true);
46 
47  // The temperature is evolving, with a first order time derivative
48  this->time_evolving(T_var, 1);
49 
50  ZeroFunction<Number> zero;
51 
52  ConstFunction<Number> one(1.0);
53 
54  // Dirichlet boundary conditions for the primal.
55  sys.get_dof_map().add_dirichlet_boundary(DirichletBoundary ({0,1,2,3}, {T_var}, &one));
56 
57  // Set adjoint boundary conditions.
58  sys.get_dof_map().add_adjoint_dirichlet_boundary(DirichletBoundary ({0,1,2,3}, {T_var}, &zero), 0);
59  sys.get_dof_map().add_adjoint_dirichlet_boundary(DirichletBoundary ({0,1,2,3}, {T_var}, &zero), 1);
60 
61 }
bool _analytic_jacobians
Definition: sigma_physics.h:57
virtual void time_evolving(unsigned int var, unsigned int order)
Tells the DiffSystem that variable var is evolving with respect to time.
Definition: diff_physics.C:41
const Number zero
.
Definition: libmesh.h:304
unsigned int T_var
Definition: sigma_physics.h:54

◆ init_physics()

void libMesh::DifferentiablePhysics::init_physics ( const System sys)
virtualinherited

Initialize any data structures associated with the physics.

Definition at line 35 of file diff_physics.C.

References libMesh::DifferentiablePhysics::_time_evolving, and libMesh::System::n_vars().

Referenced by libMesh::DifferentiableSystem::init_data().

36 {
37  // give us flags for every variable that might be time evolving
38  _time_evolving.resize(sys.n_vars(), false);
39 }
std::vector< unsigned int > _time_evolving
Stores unsigned int to tell us which variables are evolving as first order in time (1)...
Definition: diff_physics.h:543

◆ is_first_order_var()

bool libMesh::DifferentiablePhysics::is_first_order_var ( unsigned int  var) const
inlineinherited

Definition at line 509 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_first_order_vars.

510  { return _first_order_vars.find(var) != _first_order_vars.end(); }
std::set< unsigned int > _first_order_vars
Variable indices for those variables that are first order in time.
Definition: diff_physics.h:548

◆ is_second_order_var()

bool libMesh::DifferentiablePhysics::is_second_order_var ( unsigned int  var) const
inlineinherited

Definition at line 522 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_second_order_vars.

Referenced by libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns().

523  { return _second_order_vars.find(var) != _second_order_vars.end(); }
std::set< unsigned int > _second_order_vars
Variable indices for those variables that are second order in time.
Definition: diff_physics.h:553

◆ is_time_evolving()

bool libMesh::DifferentiablePhysics::is_time_evolving ( unsigned int  var) const
inlineinherited
Returns
true iff variable var is evolving with respect to time. In general, the user's init() function should have set time_evolving() for any variables which behave like du/dt = F(u), and should not call time_evolving() for any variables which behave like 0 = G(u).

Definition at line 260 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_time_evolving, and libMesh::libmesh_assert().

Referenced by libMesh::FEMPhysics::eulerian_residual(), libMesh::FEMSystem::init_context(), libMesh::FEMPhysics::mass_residual(), and libMesh::DifferentiablePhysics::nonlocal_mass_residual().

261  {
262  libmesh_assert_less(var,_time_evolving.size());
263  libmesh_assert( _time_evolving[var] == 0 ||
264  _time_evolving[var] == 1 ||
265  _time_evolving[var] == 2 );
266  return _time_evolving[var];
267  }
libmesh_assert(ctx)
std::vector< unsigned int > _time_evolving
Stores unsigned int to tell us which variables are evolving as first order in time (1)...
Definition: diff_physics.h:543

◆ k()

Real& SigmaPhysics::k ( )
inline

Definition at line 25 of file sigma_physics.h.

25 { return _k; }

◆ mass_residual()

bool libMesh::FEMPhysics::mass_residual ( bool  request_jacobian,
DiffContext c 
)
overridevirtualinherited

Subtracts a mass vector contribution on elem from elem_residual.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

Many problems can use the reimplementation in FEMPhysics::mass_residual which subtracts (du/dt,v) for each transient variable u; users with more complicated transient problems will need to reimplement this themselves.

Reimplemented from libMesh::DifferentiablePhysics.

Reimplemented in SecondOrderScalarSystemFirstOrderTimeSolverBase, SecondOrderScalarSystemSecondOrderTimeSolverBase, FirstOrderScalarSystemBase, ElasticitySystem, ElasticitySystem, and NavierSystem.

Definition at line 200 of file fem_physics.C.

References libMesh::DiffContext::elem_solution_rate_derivative, libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_jacobian(), libMesh::DiffContext::get_elem_residual(), libMesh::FEMContext::get_element_fe(), libMesh::FEMContext::get_element_qrule(), libMesh::FEAbstract::get_JxW(), libMesh::FEGenericBase< OutputType >::get_phi(), libMesh::FEMContext::interior_rate(), libMesh::DifferentiablePhysics::is_time_evolving(), libMesh::make_range(), libMesh::QBase::n_points(), and libMesh::DiffContext::n_vars().

202 {
203  FEMContext & context = cast_ref<FEMContext &>(c);
204 
205  unsigned int n_qpoints = context.get_element_qrule().n_points();
206 
207  for (auto var : make_range(context.n_vars()))
208  {
209  if (!this->is_time_evolving(var))
210  continue;
211 
212  FEBase * elem_fe = nullptr;
213  context.get_element_fe( var, elem_fe );
214 
215  const std::vector<Real> & JxW = elem_fe->get_JxW();
216 
217  const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
218 
219  const unsigned int n_dofs = cast_int<unsigned int>
220  (context.get_dof_indices(var).size());
221 
222  DenseSubVector<Number> & Fu = context.get_elem_residual(var);
223  DenseSubMatrix<Number> & Kuu = context.get_elem_jacobian( var, var );
224 
225  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
226  {
227  Number uprime;
228  context.interior_rate(var, qp, uprime);
229  const Number JxWxU = JxW[qp] * uprime;
230  for (unsigned int i = 0; i != n_dofs; ++i)
231  {
232  Fu(i) -= JxWxU * phi[i][qp];
233  if (request_jacobian && context.elem_solution_rate_derivative)
234  {
235  const Number JxWxPhiIxDeriv = JxW[qp] * phi[i][qp] *
236  context.elem_solution_rate_derivative;
237  Kuu(i,i) -= JxWxPhiIxDeriv * phi[i][qp];
238  for (unsigned int j = i+1; j < n_dofs; ++j)
239  {
240  const Number Kij = JxWxPhiIxDeriv * phi[j][qp];
241  Kuu(i,j) -= Kij;
242  Kuu(j,i) -= Kij;
243  }
244  }
245  }
246  }
247  }
248 
249  return request_jacobian;
250 }
FEGenericBase< Real > FEBase
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
bool is_time_evolving(unsigned int var) const
Definition: diff_physics.h:260

◆ nonlocal_constraint()

virtual bool libMesh::DifferentiablePhysics::nonlocal_constraint ( bool  request_jacobian,
DiffContext  
)
inlinevirtualinherited

Adds any nonlocal constraint contributions (e.g.

some components of constraints in scalar variable equations) to elem_residual

If this method receives request_jacobian = true, then it should also modify elem_jacobian and return true if possible. If the Jacobian changes have not been computed then the method should return false.

Users may need to reimplement this for PDEs on systems to which SCALAR variables with non-transient equations have been added.

Definition at line 233 of file diff_physics.h.

Referenced by libMesh::EulerSolver::nonlocal_residual(), libMesh::SteadySolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), libMesh::EigenTimeSolver::nonlocal_residual(), and libMesh::NewmarkSolver::nonlocal_residual().

235  {
236  return request_jacobian;
237  }

◆ nonlocal_damping_residual()

virtual bool libMesh::DifferentiablePhysics::nonlocal_damping_residual ( bool  request_jacobian,
DiffContext  
)
inlinevirtualinherited

Subtracts any nonlocal damping vector contributions (e.g.

any first time derivative coefficients in scalar variable equations) from elem_residual

If this method receives request_jacobian = true, then it should also modify elem_jacobian and return true if possible. If the Jacobian changes have not been computed then the method should return false.

Definition at line 394 of file diff_physics.h.

Referenced by libMesh::EulerSolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), and libMesh::NewmarkSolver::nonlocal_residual().

396  {
397  return request_jacobian;
398  }

◆ nonlocal_mass_residual()

bool libMesh::DifferentiablePhysics::nonlocal_mass_residual ( bool  request_jacobian,
DiffContext c 
)
virtualinherited

Subtracts any nonlocal mass vector contributions (e.g.

any time derivative coefficients in scalar variable equations) from elem_residual

If this method receives request_jacobian = true, then it should also modify elem_jacobian and return true if possible. If the Jacobian changes have not been computed then the method should return false.

Many problems can use the reimplementation in FEMPhysics::mass_residual which subtracts (du/dt,v) for each transient scalar variable u; users with more complicated transient scalar variable equations will need to reimplement this themselves.

Definition at line 57 of file diff_physics.C.

References libMesh::DiffContext::elem_solution_rate_derivative, libMesh::FEType::family, libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_jacobian(), libMesh::DiffContext::get_elem_residual(), libMesh::DiffContext::get_elem_solution(), libMesh::DiffContext::get_system(), libMesh::DifferentiablePhysics::is_time_evolving(), libMesh::make_range(), libMesh::DiffContext::n_vars(), libMesh::SCALAR, libMesh::Variable::type(), and libMesh::System::variable().

Referenced by libMesh::EulerSolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), libMesh::EigenTimeSolver::nonlocal_residual(), and libMesh::NewmarkSolver::nonlocal_residual().

59 {
60  FEMContext & context = cast_ref<FEMContext &>(c);
61 
62  for (auto var : make_range(context.n_vars()))
63  {
64  if (!this->is_time_evolving(var))
65  continue;
66 
67  if (c.get_system().variable(var).type().family != SCALAR)
68  continue;
69 
70  const std::vector<dof_id_type> & dof_indices =
71  context.get_dof_indices(var);
72 
73  const unsigned int n_dofs = cast_int<unsigned int>
74  (dof_indices.size());
75 
76  DenseSubVector<Number> & Fs = context.get_elem_residual(var);
77  DenseSubMatrix<Number> & Kss = context.get_elem_jacobian( var, var );
78 
80  context.get_elem_solution(var);
81 
82  for (unsigned int i=0; i != n_dofs; ++i)
83  {
84  Fs(i) -= Us(i);
85 
86  if (request_jacobian)
87  Kss(i,i) -= context.elem_solution_rate_derivative;
88  }
89  }
90 
91  return request_jacobian;
92 }
Defines a dense subvector for use in finite element computations.
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
bool is_time_evolving(unsigned int var) const
Definition: diff_physics.h:260

◆ nonlocal_time_derivative()

virtual bool libMesh::DifferentiablePhysics::nonlocal_time_derivative ( bool  request_jacobian,
DiffContext  
)
inlinevirtualinherited

Adds any nonlocal time derivative contributions (e.g.

some components of time derivatives in scalar variable equations) to elem_residual

If this method receives request_jacobian = true, then it should also modify elem_jacobian and return true if possible. If the Jacobian changes have not been computed then the method should return false.

Users may need to reimplement this for PDEs on systems to which SCALAR variables have been added.

Definition at line 214 of file diff_physics.h.

Referenced by libMesh::EulerSolver::nonlocal_residual(), libMesh::SteadySolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), libMesh::EigenTimeSolver::nonlocal_residual(), and libMesh::NewmarkSolver::nonlocal_residual().

216  {
217  return request_jacobian;
218  }

◆ operator=() [1/2]

SigmaPhysics& SigmaPhysics::operator= ( const SigmaPhysics )
default

◆ operator=() [2/2]

SigmaPhysics& SigmaPhysics::operator= ( SigmaPhysics &&  )
default

◆ set_mesh_system()

void libMesh::DifferentiablePhysics::set_mesh_system ( System sys)
inlinevirtualinherited

Tells the DifferentiablePhysics that system sys contains the isoparametric Lagrangian variables which correspond to the coordinates of mesh nodes, in problems where the mesh itself is expected to move in time.

The system with mesh coordinate data (which may be this system itself, for fully coupled moving mesh problems) is currently assumed to have new (end of time step) mesh coordinates stored in solution, old (beginning of time step) mesh coordinates stored in _old_nonlinear_solution, and constant velocity motion during each time step.

Activating this function ensures that local (but not neighbor!) element geometry is correctly repositioned when evaluating element residuals.

Currently sys must be *this for a tightly coupled moving mesh problem or nullptr to stop mesh movement; loosely coupled moving mesh problems are not implemented.

This code is experimental. "Trust but verify, and not in that order"

Definition at line 569 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_sys.

Referenced by SolidSystem::init_data().

570 {
571  // For now we assume that we're doing fully coupled mesh motion
572  // if (sys && sys != this)
573  // libmesh_not_implemented();
574 
575  // For the foreseeable future we'll assume that we keep these
576  // Systems in the same EquationSystems
577  // libmesh_assert_equal_to (&this->get_equation_systems(),
578  // &sys->get_equation_systems());
579 
580  // And for the immediate future this code may not even work
581  libmesh_experimental();
582 
583  _mesh_sys = sys;
584 }
System * _mesh_sys
System from which to acquire moving mesh information.
Definition: diff_physics.h:531

◆ set_mesh_x_var()

void libMesh::DifferentiablePhysics::set_mesh_x_var ( unsigned int  var)
inlinevirtualinherited

Tells the DifferentiablePhysics that variable var from the mesh system should be used to update the x coordinate of mesh nodes, in problems where the mesh itself is expected to move in time.

The system with mesh coordinate data (which may be this system itself, for fully coupled moving mesh problems) is currently assumed to have new (end of time step) mesh coordinates stored in solution, old (beginning of time step) mesh coordinates stored in _old_nonlinear_solution, and constant velocity motion during each time step.

Activating this function ensures that local (but not neighbor!) element geometry is correctly repositioned when evaluating element residuals.

Definition at line 589 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_x_var.

Referenced by SolidSystem::init_data().

590 {
591  _mesh_x_var = var;
592 }
unsigned int _mesh_x_var
Variables from which to acquire moving mesh information.
Definition: diff_physics.h:536

◆ set_mesh_y_var()

void libMesh::DifferentiablePhysics::set_mesh_y_var ( unsigned int  var)
inlinevirtualinherited

Tells the DifferentiablePhysics that variable var from the mesh system should be used to update the y coordinate of mesh nodes.

Definition at line 597 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_y_var.

Referenced by SolidSystem::init_data().

598 {
599  _mesh_y_var = var;
600 }

◆ set_mesh_z_var()

void libMesh::DifferentiablePhysics::set_mesh_z_var ( unsigned int  var)
inlinevirtualinherited

Tells the DifferentiablePhysics that variable var from the mesh system should be used to update the z coordinate of mesh nodes.

Definition at line 605 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_z_var.

Referenced by SolidSystem::init_data().

606 {
607  _mesh_z_var = var;
608 }

◆ side_constraint()

virtual bool libMesh::DifferentiablePhysics::side_constraint ( bool  request_jacobian,
DiffContext  
)
inlinevirtualinherited

Adds the constraint contribution on side of elem to elem_residual.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

Users may need to reimplement this for their particular PDE depending on the boundary conditions.

To implement a weak form of the constraint 0 = G(u), the user should examine u = elem_solution and add (G(u), phi_i) boundary integral contributions to elem_residual in side_constraint().

Reimplemented in LaplaceSystem, LaplaceSystem, LaplaceSystem, and NavierSystem.

Definition at line 195 of file diff_physics.h.

Referenced by libMesh::EulerSolver::side_residual(), libMesh::SteadySolver::side_residual(), libMesh::Euler2Solver::side_residual(), libMesh::EigenTimeSolver::side_residual(), and libMesh::NewmarkSolver::side_residual().

197  {
198  return request_jacobian;
199  }

◆ side_damping_residual()

virtual bool libMesh::DifferentiablePhysics::side_damping_residual ( bool  request_jacobian,
DiffContext  
)
inlinevirtualinherited

Subtracts a damping vector contribution on side of elem from elem_residual.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

For most problems, the default implementation of "do nothing" is correct; users with boundary conditions including first time derivatives may need to reimplement this themselves.

Definition at line 378 of file diff_physics.h.

Referenced by libMesh::EulerSolver::side_residual(), libMesh::Euler2Solver::side_residual(), and libMesh::NewmarkSolver::side_residual().

380  {
381  return request_jacobian;
382  }

◆ side_mass_residual()

virtual bool libMesh::DifferentiablePhysics::side_mass_residual ( bool  request_jacobian,
DiffContext  
)
inlinevirtualinherited

Subtracts a mass vector contribution on side of elem from elem_residual.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

For most problems, the default implementation of "do nothing" is correct; users with boundary conditions including time derivatives may need to reimplement this themselves.

Definition at line 320 of file diff_physics.h.

Referenced by libMesh::EulerSolver::side_residual(), libMesh::Euler2Solver::side_residual(), libMesh::EigenTimeSolver::side_residual(), and libMesh::NewmarkSolver::side_residual().

322  {
323  return request_jacobian;
324  }

◆ side_time_derivative()

virtual bool libMesh::DifferentiablePhysics::side_time_derivative ( bool  request_jacobian,
DiffContext  
)
inlinevirtualinherited

Adds the time derivative contribution on side of elem to elem_residual.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

Users may need to reimplement this for their particular PDE depending on the boundary conditions.

To implement a weak form of the source term du/dt = F(u) on sides, such as might arise in a flux boundary condition, the user should examine u = elem_solution and add (F(u), phi_i) boundary integral contributions to elem_residual in side_constraint().

Reimplemented in ElasticitySystem, ElasticitySystem, CurlCurlSystem, CurlCurlSystem, and SolidSystem.

Definition at line 174 of file diff_physics.h.

Referenced by libMesh::EulerSolver::side_residual(), libMesh::SteadySolver::side_residual(), libMesh::Euler2Solver::side_residual(), libMesh::EigenTimeSolver::side_residual(), and libMesh::NewmarkSolver::side_residual().

176  {
177  return request_jacobian;
178  }

◆ time_evolving()

void libMesh::DifferentiablePhysics::time_evolving ( unsigned int  var,
unsigned int  order 
)
virtualinherited

Tells the DiffSystem that variable var is evolving with respect to time.

In general, the user's init() function should call time_evolving() with order 1 for any variables which behave like du/dt = F(u), with order 2 for any variables that behave like d^2u/dt^2 = F(u), and should not call time_evolving() for any variables which behave like 0 = G(u).

Most derived systems will not have to reimplement this function; however any system which reimplements mass_residual() may have to reimplement time_evolving() to prepare data structures.

Definition at line 41 of file diff_physics.C.

References libMesh::DifferentiablePhysics::_first_order_vars, libMesh::DifferentiablePhysics::_second_order_vars, and libMesh::DifferentiablePhysics::_time_evolving.

Referenced by libMesh::DifferentiableSystem::add_second_order_dot_vars(), init_data(), CurlCurlSystem::init_data(), and HeatSystem::init_data().

43 {
44  libmesh_error_msg_if(order != 1 && order != 2, "Input order must be 1 or 2!");
45 
46  if (_time_evolving.size() <= var)
47  _time_evolving.resize(var+1, 0);
48 
49  _time_evolving[var] = order;
50 
51  if (order == 1)
52  _first_order_vars.insert(var);
53  else
54  _second_order_vars.insert(var);
55 }
std::set< unsigned int > _second_order_vars
Variable indices for those variables that are second order in time.
Definition: diff_physics.h:553
std::set< unsigned int > _first_order_vars
Variable indices for those variables that are first order in time.
Definition: diff_physics.h:548
std::vector< unsigned int > _time_evolving
Stores unsigned int to tell us which variables are evolving as first order in time (1)...
Definition: diff_physics.h:543

Member Data Documentation

◆ _analytic_jacobians

bool SigmaPhysics::_analytic_jacobians
protected

Definition at line 57 of file sigma_physics.h.

Referenced by element_time_derivative(), and init_data().

◆ _fe_family

std::string SigmaPhysics::_fe_family
protected

Definition at line 50 of file sigma_physics.h.

◆ _fe_order

unsigned int SigmaPhysics::_fe_order
protected

Definition at line 51 of file sigma_physics.h.

◆ _first_order_vars

std::set<unsigned int> libMesh::DifferentiablePhysics::_first_order_vars
protectedinherited

◆ _k

Real SigmaPhysics::_k
protected

Definition at line 44 of file sigma_physics.h.

Referenced by init_data().

◆ _mesh_sys

System* libMesh::DifferentiablePhysics::_mesh_sys
protectedinherited

◆ _mesh_x_var

unsigned int libMesh::DifferentiablePhysics::_mesh_x_var
protectedinherited

◆ _mesh_y_var

unsigned int libMesh::DifferentiablePhysics::_mesh_y_var
protectedinherited

◆ _mesh_z_var

unsigned int libMesh::DifferentiablePhysics::_mesh_z_var
protectedinherited

◆ _second_order_dot_vars

std::map<unsigned int,unsigned int> libMesh::DifferentiablePhysics::_second_order_dot_vars
protectedinherited

If the user adds any second order variables, then we need to also cache the map to their corresponding dot variable that will be added by this TimeSolver class.

Definition at line 560 of file diff_physics.h.

Referenced by libMesh::DifferentiableSystem::add_second_order_dot_vars(), and libMesh::DifferentiableSystem::get_second_order_dot_var().

◆ _second_order_vars

std::set<unsigned int> libMesh::DifferentiablePhysics::_second_order_vars
protectedinherited

◆ _time_evolving

std::vector<unsigned int> libMesh::DifferentiablePhysics::_time_evolving
protectedinherited

Stores unsigned int to tell us which variables are evolving as first order in time (1), second order in time (2), or are not time evolving (0).

Definition at line 543 of file diff_physics.h.

Referenced by libMesh::DifferentiablePhysics::clear_physics(), libMesh::DifferentiablePhysics::init_physics(), libMesh::DifferentiablePhysics::is_time_evolving(), and libMesh::DifferentiablePhysics::time_evolving().

◆ compute_internal_sides

bool libMesh::DifferentiablePhysics::compute_internal_sides
inherited

compute_internal_sides is false by default, indicating that side_* computations will only be done on boundary sides.

If compute_internal_sides is true, computations will be done on sides between elements as well.

Definition at line 156 of file diff_physics.h.

◆ computed_QoI

Number SigmaPhysics::computed_QoI[1]
protected

Definition at line 47 of file sigma_physics.h.

◆ T_var

unsigned int SigmaPhysics::T_var
protected

Definition at line 54 of file sigma_physics.h.

Referenced by element_time_derivative(), init_context(), and init_data().


The documentation for this class was generated from the following files: