libMesh
L-qoi.C
Go to the documentation of this file.
1 #include "L-qoi.h"
2 
3 using namespace libMesh;
4 
5 void LaplaceQoI::init_qoi(std::vector<Number> & sys_qoi)
6 {
7  // Only 1 qoi to worry about
8  sys_qoi.resize(1);
9 }
10 
12 {
13  FEMContext & c = cast_ref<FEMContext &>(context);
14 
15  // Now make sure we have requested all the data
16  // we need to build the linear system.
17  FEBase * elem_fe = nullptr;
18  c.get_element_fe(0, elem_fe);
19  elem_fe->get_JxW();
20  elem_fe->get_phi();
21  elem_fe->get_xyz();
22 }
23 
24 
25 // We only have one QoI, so we don't bother checking the qois argument
26 // to see if it was requested from us
28  const QoISet & /* qois */)
29 {
30  FEMContext & c = cast_ref<FEMContext &>(context);
31 
32  FEBase * elem_fe = nullptr;
33  c.get_element_fe(0, elem_fe);
34 
35  // Element Jacobian * quadrature weights for interior integration
36  const std::vector<Real> & JxW = elem_fe->get_JxW();
37 
38  const std::vector<Point> & xyz = elem_fe->get_xyz();
39 
40  unsigned int n_qpoints = c.get_element_qrule().n_points();
41 
42  Number dQoI_0 = 0.;
43 
44  // Loop over quadrature points
45  for (unsigned int qp = 0; qp != n_qpoints; qp++)
46  {
47  // Get co-ordinate locations of the current quadrature point
48  const Real xf = xyz[qp](0);
49  const Real yf = xyz[qp](1);
50 
51  // If in the sub-domain omega, add the contribution to the integral R
52  if (std::abs(xf - 0.875) <= 0.125 && std::abs(yf - 0.125) <= 0.125)
53  {
54  // Get the solution value at the quadrature point
55  Number T = c.interior_value(0, qp);
56 
57  // Update the elemental increment dR for each qp
58  dQoI_0 += JxW[qp] * T;
59  }
60  }
61 
62  // Update the computed value of the global functional R, by adding the contribution from this element
63  c.get_qois()[0] = c.get_qois()[0] + dQoI_0;
64 }
65 
66 // We only have one QoI, so we don't bother checking the qois argument
67 // to see if it was requested from us
69  const QoISet & /* qois */)
70 {
71  FEMContext & c = cast_ref<FEMContext &>(context);
72 
73  // First we get some references to cell-specific data that
74  // will be used to assemble the linear system.
75  FEBase * elem_fe = nullptr;
76  c.get_element_fe(0, elem_fe);
77 
78  // Element Jacobian * quadrature weights for interior integration
79  const std::vector<Real> & JxW = elem_fe->get_JxW();
80 
81  // The basis functions for the element
82  const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
83 
84  // The element quadrature points
85  const std::vector<Point > & q_point = elem_fe->get_xyz();
86 
87  // The number of local degrees of freedom in each variable
88  const unsigned int n_T_dofs = c.n_dof_indices(0);
89  unsigned int n_qpoints = c.get_element_qrule().n_points();
90 
91  // Fill the QoI RHS corresponding to this QoI. Since this is the 0th QoI
92  // we fill in the [0][i] subderivatives, i corresponding to the variable index.
93  // Our system has only one variable, so we only have to fill the [0][0] subderivative
95 
96  // Loop over the qps
97  for (unsigned int qp=0; qp != n_qpoints; qp++)
98  {
99  const Real x = q_point[qp](0);
100  const Real y = q_point[qp](1);
101 
102  // If in the sub-domain over which QoI 0 is supported, add contributions
103  // to the adjoint rhs
104  if (std::abs(x - 0.875) <= 0.125 && std::abs(y - 0.125) <= 0.125)
105  for (unsigned int i=0; i != n_T_dofs; i++)
106  Q(i) += JxW[qp]*phi[i][qp];
107  } // end of the quadrature point qp-loop
108 }
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::FEMContext::get_element_qrule
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
Definition: fem_context.h:790
libMesh::FEAbstract::get_xyz
const std::vector< Point > & get_xyz() const
Definition: fe_abstract.h:237
libMesh::QBase::n_points
unsigned int n_points() const
Definition: quadrature.h:126
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::FEGenericBase
This class forms the foundation from which generic finite elements may be derived.
Definition: exact_error_estimator.h:39
libMesh::FEAbstract::get_JxW
const std::vector< Real > & get_JxW() const
Definition: fe_abstract.h:244
LaplaceQoI::element_qoi_derivative
virtual void element_qoi_derivative(DiffContext &context, const QoISet &qois)
Does any work that needs to be done on elem in a quantity of interest derivative assembly loop,...
Definition: L-qoi.C:68
libMesh::DiffContext::get_qois
const std::vector< Number > & get_qois() const
Const accessor for QoI vector.
Definition: diff_context.h:319
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::DiffContext::get_qoi_derivatives
const std::vector< DenseVector< Number > > & get_qoi_derivatives() const
Const accessor for QoI derivatives.
Definition: diff_context.h:331
libMesh::QoISet
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
libMesh::DiffContext::n_dof_indices
unsigned int n_dof_indices() const
Total number of dof indices on the element.
Definition: diff_context.h:399
LaplaceQoI::init_qoi
virtual void init_qoi(std::vector< Number > &sys_qoi)
Definition: L-qoi.C:5
libMesh::DiffContext
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
libMesh::DenseSubVector< Number >
libMesh::FEMContext::get_element_fe
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:275
L-qoi.h
LaplaceQoI::element_qoi
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,...
Definition: L-qoi.C:27
LaplaceQoI::init_context
virtual void init_context(DiffContext &context)
Prepares the result of a build_context() call for use.
Definition: L-qoi.C:11
libMesh::FEGenericBase::get_phi
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:206
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::FEMContext::interior_value
Number interior_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:371
libMesh::FEMContext
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62