libMesh
side_postprocess.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 #include "libmesh/threads.h"
9 
10 // Local includes
11 #include "L-shaped.h"
12 
13 // Bring in everything from the libMesh namespace
14 using namespace libMesh;
15 
16 // Define the postprocess function to compute QoI 1, the integral of the the normal
17 // derivative of the solution over part of the boundary
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 * side_fe = nullptr;
25  c.get_side_fe(0, side_fe);
26 
27  // Element Jacobian * quadrature weights for interior integration
28  const std::vector<Real> & JxW = side_fe->get_JxW();
29 
30  const std::vector<Point > & q_point = side_fe->get_xyz();
31 
32  const std::vector<Point> & face_normals = side_fe->get_normals();
33 
34  unsigned int n_qpoints = c.get_side_qrule().n_points();
35 
36  Number dQoI_1 = 0.;
37 
38  // Loop over qp's, compute the function at each qp and add
39  // to get the QoI
40  for (unsigned int qp=0; qp != n_qpoints; qp++)
41  {
42  const Real x = q_point[qp](0);
43  const Real y = q_point[qp](1);
44 
45  const Real TOL = 1.e-5;
46 
47  // If on the bottom horizontal bdry (y = -1)
48  if (std::abs(y - 1.0) <= TOL && x > 0.0)
49  {
50  // Get the value of the gradient at this point
51  const Gradient grad_T = c.side_gradient(0, qp);
52 
53  // Add the contribution of this qp to the integral QoI
54  dQoI_1 += JxW[qp] * (grad_T * face_normals[qp] * x * (x - 1.));
55  }
56 
57  } // end of the quadrature point qp-loop
58 
59  static Threads::spin_mutex local_mtx;
60  Threads::spin_mutex::scoped_lock lock(local_mtx);
61  computed_QoI[1] = computed_QoI[1] + dQoI_1;
62 }
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
void get_side_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for edge/face (2D/3D) finite element object for variable var for the largest dimension in th...
Definition: fem_context.h:317
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
The libMesh namespace provides an interface to certain functionality in the library.
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62
virtual void side_postprocess(DiffContext &context)
Does any work that needs to be done on side of elem in a postprocessing loop.
unsigned int n_points() const
Definition: quadrature.h:131
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Gradient side_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:719
This class forms the foundation from which generic finite elements may be derived.
const QBase & get_side_qrule() const
Accessor for element side quadrature rule for the dimension of the current _elem. ...
Definition: fem_context.h:809