20 #include "libmesh/getpot.h" 22 #include "heatsystem.h" 24 #include "libmesh/fe_base.h" 25 #include "libmesh/fe_interface.h" 26 #include "libmesh/fem_context.h" 27 #include "libmesh/mesh.h" 28 #include "libmesh/quadrature.h" 29 #include "libmesh/string_to_enum.h" 30 #include "libmesh/system.h" 31 #include "libmesh/equation_systems.h" 32 #include "libmesh/zero_function.h" 33 #include "libmesh/dirichlet_boundaries.h" 34 #include "libmesh/dof_map.h" 35 #include "libmesh/numeric_vector.h" 40 Utility::string_to_enum<FEFamily>(
_fe_family));
44 std::ifstream i(
"heat.in");
45 libmesh_error_msg_if(!i,
'[' << this->
processor_id() <<
"] Can't find heat.in; exiting early.");
47 GetPot infile(
"heat.in");
48 _k = infile(
"k", 1.0);
59 ZeroFunction<Number>
zero;
61 #ifdef LIBMESH_ENABLE_DIRICHLET 65 (DirichletBoundary ({0,1,2,3}, {
T_var},
zero,
69 FEMSystem::init_data();
76 FEMContext & c = cast_ref<FEMContext &>(context);
78 FEBase * elem_fe =
nullptr;
79 c.get_element_fe(0, elem_fe);
87 FEBase * side_fe =
nullptr;
88 c.get_side_fe(0, side_fe);
89 side_fe->get_nothing();
96 const System & sys = c.get_system();
103 c.add_localized_vector(adjoint_solution, sys);
106 FEMSystem::init_context(context);
109 #define optassert(X) {if (!(X)) libmesh_error_msg("Assertion " #X " failed.");} 113 DiffContext & context)
117 FEMContext & c = cast_ref<FEMContext &>(context);
121 FEBase * elem_fe =
nullptr;
122 c.get_element_fe(0, elem_fe);
125 const std::vector<Real> & JxW = elem_fe->get_JxW();
128 const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
131 const unsigned int n_u_dofs = c.n_dof_indices(0);
134 DenseSubMatrix<Number> & K = c.get_elem_jacobian(0, 0);
135 DenseSubVector<Number> & F = c.get_elem_residual(0);
143 unsigned int n_qpoints = c.get_element_qrule().n_points();
145 for (
unsigned int qp=0; qp != n_qpoints; qp++)
148 Gradient grad_T = c.interior_gradient(0, qp);
150 for (
unsigned int i=0; i != n_u_dofs; i++)
151 F(i) += JxW[qp] * -
parameters[0] * (grad_T * dphi[i][qp]);
153 for (
unsigned int i=0; i != n_u_dofs; i++)
154 for (
unsigned int j=0; j != n_u_dofs; ++j)
155 K(i,j) += JxW[qp] * -
parameters[0] * (dphi[i][qp] * dphi[j][qp]);
169 for (std::size_t j=0, Np = parameters_in.size(); j != Np; ++j)
171 Number old_parameter = *parameters_in[j];
173 *parameters_in[j] = old_parameter -
dp;
179 std::unique_ptr<NumericVector<Number>> R_minus = this->
rhs->
clone();
183 R_minus_dp.push_back(-R_minus->dot(this->get_adjoint_solution(0)));
185 *parameters_in[j] = old_parameter +
dp;
191 std::unique_ptr<NumericVector<Number>> R_plus = this->
rhs->
clone();
193 R_plus_dp.push_back(-R_plus->dot(this->get_adjoint_solution(0)));
195 *parameters_in[j] = old_parameter;
std::vector< Number > deltat_vector
virtual bool element_time_derivative(bool request_jacobian, DiffContext &context)
Adds the time derivative contribution on elem to elem_residual.
const EquationSystems & get_equation_systems() const
NumericVector< Number > * rhs
The system matrix.
virtual std::unique_ptr< NumericVector< T > > clone() const =0
virtual void time_evolving(unsigned int var, unsigned int order)
Tells the DiffSystem that variable var is evolving with respect to time.
void compute_jacobian(const NumericVector< Number > &, SparseMatrix< Number > &J, NonlinearImplicitSystem &system)
std::vector< Number > R_plus_dp
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time...
virtual void init_context(DiffContext &context)
NumberVectorValue Gradient
unsigned int add_variable(std::string_view var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
std::vector< Number > parameters
FEGenericBase< Real > FEBase
void perturb_accumulate_residuals(ParameterVector ¶meters)
virtual void close()=0
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
virtual void update()
Update the local values to reflect the solution on neighboring processors.
std::vector< Number > R_minus_dp
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
T & set(const std::string &)
virtual void init_data()
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used...
NumericVector< Number > & get_adjoint_solution(unsigned int i=0)
Parameters parameters
Data structure holding arbitrary parameters.
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false) override
Prepares matrix or rhs for matrix assembly.
processor_id_type processor_id() const
System(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
const DofMap & get_dof_map() const
template class LIBMESH_EXPORT NumericVector< Number >