libMesh
element_qoi_derivative.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 "L-shaped.h"
11 
12 // Bring in everything from the libMesh namespace
13 using namespace libMesh;
14 
15 // We only have one QoI, so we don't bother checking the qois argument
16 // to see if it was requested from us
18  const QoISet & /* qois */)
19 {
20  FEMContext & c = cast_ref<FEMContext &>(context);
21 
22  // First we get some references to cell-specific data that
23  // will be used to assemble the linear system.
24  FEBase * elem_fe = nullptr;
25  c.get_element_fe(0, elem_fe);
26 
27  // Element Jacobian * quadrature weights for interior integration
28  const std::vector<Real> & JxW = elem_fe->get_JxW();
29 
30  // The basis functions for the element
31  const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
32 
33  // The element quadrature points
34  const std::vector<Point > & q_point = elem_fe->get_xyz();
35 
36  // The number of local degrees of freedom in each variable
37  const unsigned int n_T_dofs = c.n_dof_indices(0);
38  unsigned int n_qpoints = c.get_element_qrule().n_points();
39 
40  // Fill the QoI RHS corresponding to this QoI. Since this is the 0th QoI
41  // we fill in the [0][i] subderivatives, i corresponding to the variable index.
42  // Our system has only one variable, so we only have to fill the [0][0] subderivative
44 
45  // Loop over the qps
46  for (unsigned int qp=0; qp != n_qpoints; qp++)
47  {
48  const Real x = q_point[qp](0);
49  const Real y = q_point[qp](1);
50 
51  // If in the sub-domain over which QoI 0 is supported, add contributions
52  // to the adjoint rhs
53  if (std::abs(x - 0.875) <= 0.125 && std::abs(y - 0.125) <= 0.125)
54  for (unsigned int i=0; i != n_T_dofs; i++)
55  Q(i) += JxW[qp] *phi[i][qp];
56  } // end of the quadrature point qp-loop
57 }
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
LaplaceSystem::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: element_qoi_derivative.C:17
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
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
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
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62