libMesh
element_qoi.C
Go to the documentation of this file.
1 // General libMesh includes
2 #include "libmesh/libmesh_common.h"
3 #include "libmesh/elem.h"
4 #include "libmesh/fe_base.h"
5 #include "libmesh/fem_context.h"
6 #include "libmesh/point.h"
7 #include "libmesh/quadrature.h"
8 
9 // Local includes
10 #include "heatsystem.h"
11 
12 // Bring in everything from the libMesh namespace
13 using namespace libMesh;
14 
15 // Define the postprocess function to compute QoI 0, the weighted flux
16 
17 void HeatSystem::element_qoi (DiffContext & context, const QoISet & /* qois */)
18 {
19  FEMContext & c = cast_ref<FEMContext &>(context);
20 
21  FEBase * elem_fe = nullptr;
22  c.get_element_fe( 0, elem_fe );
23 
24  // Element Jacobian * quadrature weights for interior integration
25  const std::vector<Real> & JxW = elem_fe->get_JxW();
26 
27  // The number of local degrees of freedom in each variable
28 
29  unsigned int n_qpoints = c.get_element_qrule().n_points();
30 
31  // A reference to the system context is built with
32  const System & sys = c.get_system();
33 
34  Number dQoI_0 = 0.;
35  Number dQoI_1 = 0.;
36 
37  // Smooth QoI weight
38  Number q = 0.0;
39 
40  // Loop over quadrature points
41  for (unsigned int qp = 0; qp != n_qpoints; qp++)
42  {
43  Number u = c.interior_value(0, qp);
44 
45  // The final time QoI
46  if(std::abs(sys.time - (dynamic_cast<const HeatSystem &>(sys).QoI_time_instant)[0]) < TOLERANCE)
47  {
48  dQoI_0 += JxW[qp] * ( u * 1.0 );
49  }
50 
51  // Smooth QoI contribution
52  q = 1.0;
53 
54  dQoI_1 += (JxW[qp] * ( u * q ));
55  }
56  // End loop over quadrature points
57 
58  // Update the computed value of the global functional R, by adding the contribution from this element
59  c.get_qois()[0] = c.get_qois()[0] + dQoI_0;
60  c.get_qois()[1] = c.get_qois()[1] + dQoI_1;
61 
62 }
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Definition: system.h:1615
const std::vector< Number > & get_qois() const
Const accessor for QoI vector.
Definition: diff_context.h:317
Number interior_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:412
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
static constexpr Real TOLERANCE
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
The libMesh namespace provides an interface to certain functionality in the library.
const System & get_system() const
Accessor for associated system.
Definition: diff_context.h:106
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62
unsigned int n_points() const
Definition: quadrature.h:131
virtual void element_qoi(DiffContext &context, const QoISet &qois)
Does any work that needs to be done on elem in a quantity of interest assembly loop, outputting to elem_qoi.
Definition: element_qoi.C:17
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for interior finite element object for variable var for the largest dimension in the mesh...
Definition: fem_context.h:277
This class forms the foundation from which generic finite elements may be derived.
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
Definition: fem_context.h:802