libMesh
Public Member Functions | Public Attributes | List of all members
LaplaceQoI Class Reference

#include <L-qoi.h>

Inheritance diagram for LaplaceQoI:
[legend]

Public Member Functions

 LaplaceQoI ()
 
virtual ~LaplaceQoI ()
 
virtual void init_qoi (std::vector< Number > &sys_qoi)
 
virtual void init_context (DiffContext &context)
 Prepares the result of a build_context() call for use. More...
 
virtual void postprocess ()
 
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, outputting to elem_qoi_derivative. More...
 
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. More...
 
virtual UniquePtr< DifferentiableQoIclone ()
 Copy of this object. More...
 
virtual void init_qoi (std::vector< Number > &)
 Initialize system qoi. More...
 
virtual void clear_qoi ()
 Clear all the data structures associated with the QoI. More...
 
virtual void side_qoi (DiffContext &, const QoISet &)
 Does any work that needs to be done on side of elem in a quantity of interest assembly loop, outputting to elem_qoi. More...
 
virtual void side_qoi_derivative (DiffContext &, const QoISet &)
 Does any work that needs to be done on side of elem in a quantity of interest derivative assembly loop, outputting to elem_qoi_derivative. More...
 
virtual void thread_join (std::vector< Number > &qoi, const std::vector< Number > &other_qoi, const QoISet &qoi_indices)
 Method to combine thread-local qois. More...
 
virtual void parallel_op (const Parallel::Communicator &communicator, std::vector< Number > &sys_qoi, std::vector< Number > &local_qoi, const QoISet &qoi_indices)
 Method to populate system qoi data structure with process-local qoi. More...
 
virtual void finalize_derivative (NumericVector< Number > &derivatives, std::size_t qoi_index)
 Method to finalize qoi derivatives which require more than just a simple sum of element contributions. More...
 

Public Attributes

bool assemble_qoi_sides
 If assemble_qoi_sides is true (it is false by default), the assembly loop for a quantity of interest or its derivatives will loop over domain boundary sides. More...
 
bool assemble_qoi_internal_sides
 If assemble_qoi_internal_sides is true (it is false by default), the assembly loop for a quantity of interest or its derivatives will loop over element sides which do not fall on domain boundaries. More...
 
bool assemble_qoi_elements
 If assemble_qoi_elements is false (it is true by default), the assembly loop for a quantity of interest or its derivatives will skip computing on mesh elements, and will only compute on mesh sides. More...
 

Detailed Description

Definition at line 16 of file L-qoi.h.

Constructor & Destructor Documentation

LaplaceQoI::LaplaceQoI ( )

Definition at line 19 of file L-qoi.h.

19 {}
virtual LaplaceQoI::~LaplaceQoI ( )
virtual

Definition at line 20 of file L-qoi.h.

20 {}

Member Function Documentation

virtual void libMesh::DifferentiableQoI::clear_qoi ( )
virtualinherited

Clear all the data structures associated with the QoI.

Definition at line 74 of file diff_qoi.h.

Referenced by libMesh::DifferentiableSystem::clear().

74 {}
virtual UniquePtr<DifferentiableQoI> LaplaceQoI::clone ( )
virtual

Copy of this object.

User should override to copy any needed state.

Implements libMesh::DifferentiableQoI.

Definition at line 33 of file L-qoi.h.

34  {
35  return UniquePtr<DifferentiableQoI> (new LaplaceQoI(*this));
36  }
std::unique_ptr< T > UniquePtr
Definition: auto_ptr.h:46
LaplaceQoI()
Definition: L-qoi.h:19
void LaplaceQoI::element_qoi ( DiffContext ,
const QoISet  
)
virtual

Does any work that needs to be done on elem in a quantity of interest assembly loop, outputting to elem_qoi.

Only qois included in the supplied QoISet need to be assembled.

Reimplemented from libMesh::DifferentiableQoI.

Definition at line 27 of file L-qoi.C.

References std::abs(), libMesh::FEMContext::get_element_fe(), libMesh::FEMContext::get_element_qrule(), libMesh::DiffContext::get_qois(), libMesh::FEMContext::interior_value(), libmesh_nullptr, libMesh::QBase::n_points(), and libMesh::Real.

29 {
30  FEMContext & c = cast_ref<FEMContext &>(context);
31 
32  FEBase * elem_fe = libmesh_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 }
double abs(double a)
const class libmesh_nullptr_t libmesh_nullptr
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:262
Number interior_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:326
const std::vector< Number > & get_qois() const
Const accessor for QoI vector.
Definition: diff_context.h:318
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:61
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
Definition: fem_context.h:765
unsigned int n_points() const
Definition: quadrature.h:113
This class forms the foundation from which generic finite elements may be derived.
void LaplaceQoI::element_qoi_derivative ( DiffContext ,
const QoISet  
)
virtual

Does any work that needs to be done on elem in a quantity of interest derivative assembly loop, outputting to elem_qoi_derivative.

Only qois included in the supplied QoISet need their derivatives assembled.

Reimplemented from libMesh::DifferentiableQoI.

Definition at line 68 of file L-qoi.C.

References std::abs(), libMesh::DiffContext::get_dof_indices(), libMesh::FEMContext::get_element_fe(), libMesh::FEMContext::get_element_qrule(), libMesh::DiffContext::get_qoi_derivatives(), libmesh_nullptr, libMesh::QBase::n_points(), libMesh::Real, and libMesh::x.

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 = libmesh_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.get_dof_indices(0).size();
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 }
double abs(double a)
const std::vector< dof_id_type > & get_dof_indices() const
Accessor for element dof indices.
Definition: diff_context.h:366
const class libmesh_nullptr_t libmesh_nullptr
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:262
Defines a dense subvector for use in finite element computations.
PetscErrorCode Vec x
const std::vector< DenseVector< Number > > & get_qoi_derivatives() const
Const accessor for QoI derivatives.
Definition: diff_context.h:330
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:61
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
Definition: fem_context.h:765
unsigned int n_points() const
Definition: quadrature.h:113
This class forms the foundation from which generic finite elements may be derived.
void libMesh::DifferentiableQoI::finalize_derivative ( NumericVector< Number > &  derivatives,
std::size_t  qoi_index 
)
virtualinherited

Method to finalize qoi derivatives which require more than just a simple sum of element contributions.

Definition at line 51 of file diff_qoi.C.

Referenced by libMesh::FEMSystem::assemble_qoi_derivative(), and libMesh::DifferentiableQoI::init_context().

52 {
53  // by default, do nothing
54 }
void LaplaceQoI::init_context ( DiffContext )
virtual

Prepares the result of a build_context() call for use.

Most FEMSystem-based problems will need to reimplement this in order to call FE::get_*() as their particular QoI requires.

Reimplemented from libMesh::DifferentiableQoI.

Definition at line 11 of file L-qoi.C.

References libMesh::FEMContext::get_element_fe(), and libmesh_nullptr.

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 = libmesh_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 }
const class libmesh_nullptr_t libmesh_nullptr
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:262
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:61
This class forms the foundation from which generic finite elements may be derived.
void LaplaceQoI::init_qoi ( std::vector< Number > &  sys_qoi)
virtual

Definition at line 5 of file L-qoi.C.

6 {
7  // Only 1 qoi to worry about
8  sys_qoi.resize(1);
9 }
virtual void libMesh::DifferentiableQoI::init_qoi ( std::vector< Number > &  )
virtualinherited

Initialize system qoi.

By default, does nothing in order to maintain backward compatibility for FEMSystem applications that control qoi.

Definition at line 68 of file diff_qoi.h.

Referenced by libMesh::DifferentiableSystem::attach_qoi().

68 {}
void libMesh::DifferentiableQoI::parallel_op ( const Parallel::Communicator communicator,
std::vector< Number > &  sys_qoi,
std::vector< Number > &  local_qoi,
const QoISet qoi_indices 
)
virtualinherited

Method to populate system qoi data structure with process-local qoi.

By default, simply sums process qois into system qoi.

Definition at line 39 of file diff_qoi.C.

References libMesh::Parallel::Communicator::sum().

Referenced by libMesh::FEMSystem::assemble_qoi(), and libMesh::DifferentiableQoI::init_context().

43 {
44  // Sum everything into local_qoi
45  communicator.sum(local_qoi);
46 
47  // Now put into system qoi
48  sys_qoi = local_qoi;
49 }
MPI_Comm communicator
Communicator object for talking with subsets of processors.
Definition: parallel.h:181
virtual void LaplaceQoI::postprocess ( void  )
virtual

Definition at line 27 of file L-qoi.h.

27 {}
virtual void libMesh::DifferentiableQoI::side_qoi ( DiffContext ,
const QoISet  
)
virtualinherited

Does any work that needs to be done on side of elem in a quantity of interest assembly loop, outputting to elem_qoi.

Only qois included in the supplied QoISet need to be assembled.

Reimplemented in CoupledSystemQoI.

Definition at line 130 of file diff_qoi.h.

132  {}
virtual void libMesh::DifferentiableQoI::side_qoi_derivative ( DiffContext ,
const QoISet  
)
virtualinherited

Does any work that needs to be done on side of elem in a quantity of interest derivative assembly loop, outputting to elem_qoi_derivative.

Only qois included in the supplied QoISet need their derivatives assembled.

Reimplemented in LaplaceSystem, LaplaceSystem, and CoupledSystemQoI.

Definition at line 142 of file diff_qoi.h.

144  {}
void libMesh::DifferentiableQoI::thread_join ( std::vector< Number > &  qoi,
const std::vector< Number > &  other_qoi,
const QoISet qoi_indices 
)
virtualinherited

Method to combine thread-local qois.

By default, simply sums thread qois.

Definition at line 31 of file diff_qoi.C.

Referenced by libMesh::DifferentiableQoI::init_context().

34 {
35  for (std::size_t i=0; i != qoi.size(); ++i)
36  qoi[i] += other_qoi[i];
37 }

Member Data Documentation

bool libMesh::DifferentiableQoI::assemble_qoi_elements
inherited

If assemble_qoi_elements is false (it is true by default), the assembly loop for a quantity of interest or its derivatives will skip computing on mesh elements, and will only compute on mesh sides.

Definition at line 98 of file diff_qoi.h.

bool libMesh::DifferentiableQoI::assemble_qoi_internal_sides
inherited

If assemble_qoi_internal_sides is true (it is false by default), the assembly loop for a quantity of interest or its derivatives will loop over element sides which do not fall on domain boundaries.

Definition at line 90 of file diff_qoi.h.

bool libMesh::DifferentiableQoI::assemble_qoi_sides
inherited

If assemble_qoi_sides is true (it is false by default), the assembly loop for a quantity of interest or its derivatives will loop over domain boundary sides.

To add domain interior sides, also set assemble_qoi_internal_sides to true.

Definition at line 82 of file diff_qoi.h.

Referenced by main().


The documentation for this class was generated from the following files: