libMesh
steady_solver.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 #include "libmesh/diff_solver.h"
20 #include "libmesh/diff_system.h"
21 #include "libmesh/steady_solver.h"
22 
23 #include "libmesh/adjoint_refinement_estimator.h"
24 #include "libmesh/error_vector.h"
25 
26 namespace libMesh
27 {
28 
29 
30 
31 SteadySolver::~SteadySolver () = default;
32 
33 
34 
35 bool SteadySolver::element_residual(bool request_jacobian,
36  DiffContext & context)
37 {
38  return this->_general_residual(request_jacobian,
39  context,
42 }
43 
44 
45 
46 bool SteadySolver::side_residual(bool request_jacobian,
47  DiffContext & context)
48 {
49  return this->_general_residual(request_jacobian,
50  context,
53 }
54 
55 
56 
57 bool SteadySolver::nonlocal_residual(bool request_jacobian,
58  DiffContext & context)
59 {
60  return this->_general_residual(request_jacobian,
61  context,
64 }
65 
66 
67 
68 bool SteadySolver::_general_residual(bool request_jacobian,
69  DiffContext & context,
70  ResFuncType time_deriv,
71  ResFuncType constraint)
72 {
73  // If a fixed solution is requested, it will just be the current
74  // solution
76  {
77  context.get_elem_fixed_solution() = context.get_elem_solution();
78  context.fixed_solution_derivative = 1.0;
79  }
80 
81  bool jacobian_computed =
82  (_system.get_physics()->*time_deriv)(request_jacobian, context);
83 
84  // The user shouldn't compute a jacobian unless requested
85  libmesh_assert (request_jacobian || !jacobian_computed);
86 
87  bool jacobian_computed2 =
88  (_system.get_physics()->*constraint)(jacobian_computed, context);
89 
90  // The user shouldn't compute a jacobian unless requested
91  libmesh_assert (jacobian_computed || !jacobian_computed2);
92 
93  return jacobian_computed2;
94 }
95 
97 {
98  this->_system.assemble_qoi();
99 
100  return;
101 }
102 
103 void SteadySolver::integrate_adjoint_sensitivity(const QoISet & qois, const ParameterVector & parameter_vector, SensitivityData & sensitivities)
104 {
105  this->_system.ImplicitSystem::adjoint_qoi_parameter_sensitivity(qois, parameter_vector, sensitivities);
106 
107  return;
108 }
109 
110 #ifdef LIBMESH_ENABLE_AMR
112  (AdjointRefinementEstimator & adjoint_refinement_error_estimator,
113  ErrorVector & QoI_elementwise_error)
114 {
115  // Base class assumes a direct steady state error estimate
116  adjoint_refinement_error_estimator.estimate_error(_system, QoI_elementwise_error);
117 
118  // Also get the spatially integrated errors for all the QoIs in the QoI set
119  for (auto j : make_range(_system.n_qois()))
120  {
121  // Skip this QoI if not in the QoI Set
122  if (adjoint_refinement_error_estimator.qoi_set().has_index(j))
123  {
124  _system.set_qoi_error_estimate(j, adjoint_refinement_error_estimator.get_global_QoI_error_estimate(j));
125  }
126  }
127 
128  return;
129 }
130 #endif // LIBMESH_ENABLE_AMR
131 
132 } // namespace libMesh
virtual bool side_constraint(bool request_jacobian, DiffContext &)
Adds the constraint contribution on side of elem to elem_residual.
Definition: diff_physics.h:195
QoISet & qoi_set()
Access to the QoISet (default: weight all QoIs equally) to use when computing errors.
Real fixed_solution_derivative
The derivative of elem_fixed_solution with respect to the nonlinear solution, for use by systems cons...
Definition: diff_context.h:517
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
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 ~SteadySolver()
Destructor.
Data structure for specifying which Parameters should be independent variables in a parameter sensiti...
const DenseVector< Number > & get_elem_fixed_solution() const
Accessor for element fixed solution.
Definition: diff_context.h:210
This class implements a "brute force" goal-oriented error estimator which computes an estimate of err...
virtual void integrate_qoi_timestep() override
A method to integrate the system::QoI functionals.
Definition: steady_solver.C:96
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
Definition: error_vector.h:50
virtual void integrate_adjoint_refinement_error_estimate(AdjointRefinementEstimator &adjoint_refinement_error_estimator, ErrorVector &QoI_elementwise_error) override
A method to compute the adjoint refinement error estimate at the current timestep.
virtual bool nonlocal_time_derivative(bool request_jacobian, DiffContext &)
Adds any nonlocal time derivative contributions (e.g.
Definition: diff_physics.h:214
The libMesh namespace provides an interface to certain functionality in the library.
virtual bool nonlocal_constraint(bool request_jacobian, DiffContext &)
Adds any nonlocal constraint contributions (e.g.
Definition: diff_physics.h:233
bool has_index(std::size_t) const
Return whether or not this index is in the set to be calculated.
Definition: qoi_set.h:224
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
const DenseVector< Number > & get_elem_solution() const
Accessor for element solution.
Definition: diff_context.h:112
Data structure for holding completed parameter sensitivity calculations.
bool use_fixed_solution
A boolean to be set to true by systems using elem_fixed_solution, for optional use by e...
Definition: system.h:1563
libmesh_assert(ctx)
Number & get_global_QoI_error_estimate(unsigned int qoi_index)
This is an accessor function to access the computed global QoI error estimates.
virtual bool _general_residual(bool request_jacobian, DiffContext &, ResFuncType time_deriv, ResFuncType constraint)
This method is the underlying implementation of the public residual methods.
Definition: steady_solver.C:68
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false) override
This function does uniform refinements and an adjoint solve to get an adjoint solution on each cell...
const DifferentiablePhysics * get_physics() const
Definition: diff_system.h:181
virtual void assemble_qoi(const QoISet &qoi_indices=QoISet()) override
Prepares qoi for quantity of interest assembly, then calls user qoi function.
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
virtual bool nonlocal_residual(bool request_jacobian, DiffContext &) override
This method uses the DifferentiablePhysics&#39; nonlocal_time_derivative() and nonlocal_constraint() to b...
Definition: steady_solver.C:57
virtual bool element_residual(bool request_jacobian, DiffContext &) override
This method uses the DifferentiablePhysics&#39; element_time_derivative() and element_constraint() to bui...
Definition: steady_solver.C:35
virtual bool side_residual(bool request_jacobian, DiffContext &) override
This method uses the DifferentiablePhysics&#39; side_time_derivative() and side_constraint() to build a f...
Definition: steady_solver.C:46
virtual void integrate_adjoint_sensitivity(const QoISet &qois, const ParameterVector &parameter_vector, SensitivityData &sensitivities) override
A method to integrate the adjoint sensitivity w.r.t a given parameter vector.
virtual bool side_time_derivative(bool request_jacobian, DiffContext &)
Adds the time derivative contribution on side of elem to elem_residual.
Definition: diff_physics.h:174
virtual bool element_constraint(bool request_jacobian, DiffContext &)
Adds the constraint contribution on elem to elem_residual.
Definition: diff_physics.h:144