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

#include <H-qoi.h>

Inheritance diagram for CoupledSystemQoI:
[legend]

Public Member Functions

 CoupledSystemQoI ()
 
virtual ~CoupledSystemQoI ()
 
virtual void init_qoi (std::vector< Number > &sys_qoi)
 
virtual void postprocess ()
 
virtual void side_qoi_derivative (DiffContext &context, const QoISet &qois)
 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 side_qoi (DiffContext &context, const QoISet &qois)
 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 std::unique_ptr< 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 element_qoi (DiffContext &, const QoISet &)
 Does any work that needs to be done on elem in a quantity of interest assembly loop, outputting to elem_qoi. More...
 
virtual void element_qoi_derivative (DiffContext &, const QoISet &)
 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 init_context (DiffContext &)
 Prepares the result of a build_context() call for use. 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 15 of file H-qoi.h.

Constructor & Destructor Documentation

◆ CoupledSystemQoI()

CoupledSystemQoI::CoupledSystemQoI ( )
inline

Definition at line 18 of file H-qoi.h.

18 {}

◆ ~CoupledSystemQoI()

virtual CoupledSystemQoI::~CoupledSystemQoI ( )
inlinevirtual

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

19 {}

Member Function Documentation

◆ clear_qoi()

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

Clear all the data structures associated with the QoI.

Definition at line 77 of file diff_qoi.h.

77 {}

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

◆ clone()

virtual std::unique_ptr<DifferentiableQoI> CoupledSystemQoI::clone ( )
inlinevirtual

Copy of this object.

User should override to copy any needed state.

Implements libMesh::DifferentiableQoI.

Definition at line 30 of file H-qoi.h.

31  {
32  DifferentiableQoI * my_clone = new CoupledSystemQoI;
33  *my_clone = *this;
34  return std::unique_ptr<DifferentiableQoI>(my_clone);
35  }

◆ element_qoi()

virtual void libMesh::DifferentiableQoI::element_qoi ( DiffContext ,
const QoISet  
)
inlinevirtualinherited

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 in LaplaceQoI.

Definition at line 110 of file diff_qoi.h.

112  {}

◆ element_qoi_derivative()

virtual void libMesh::DifferentiableQoI::element_qoi_derivative ( DiffContext ,
const QoISet  
)
inlinevirtualinherited

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 in LaplaceSystem, LaplaceSystem, LaplaceQoI, and HeatSystem.

Definition at line 122 of file diff_qoi.h.

124  {}

◆ finalize_derivative()

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 53 of file diff_qoi.C.

54 {
55  // by default, do nothing
56 }

Referenced by libMesh::FEMSystem::assemble_qoi_derivative().

◆ init_context()

virtual void libMesh::DifferentiableQoI::init_context ( DiffContext )
inlinevirtualinherited

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 in L2System, HeatSystem, HeatSystem, CoupledSystem, ElasticitySystem, LaplaceSystem, CurlCurlSystem, LaplaceSystem, PoissonSystem, LaplaceSystem, LaplaceSystem, CurlCurlSystem, SolidSystem, NavierSystem, LaplaceQoI, and libMesh::FEMSystem.

Definition at line 155 of file diff_qoi.h.

155 {}

◆ init_qoi() [1/2]

virtual void libMesh::DifferentiableQoI::init_qoi ( std::vector< Number > &  )
inlinevirtualinherited

Initialize system qoi.

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

Definition at line 71 of file diff_qoi.h.

71 {}

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

◆ init_qoi() [2/2]

void CoupledSystemQoI::init_qoi ( std::vector< Number > &  sys_qoi)
virtual

Definition at line 8 of file H-qoi.C.

9 {
10  //Only 1 qoi to worry about
11  sys_qoi.resize(1);
12 }

◆ parallel_op()

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 41 of file diff_qoi.C.

45 {
46  // Sum everything into local_qoi
47  communicator.sum(local_qoi);
48 
49  // Now put into system qoi
50  sys_qoi = local_qoi;
51 }

Referenced by libMesh::FEMSystem::assemble_qoi().

◆ postprocess()

virtual void CoupledSystemQoI::postprocess ( void  )
inlinevirtual

Definition at line 22 of file H-qoi.h.

22 {}

◆ side_qoi()

void CoupledSystemQoI::side_qoi ( DiffContext ,
const QoISet  
)
virtual

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 from libMesh::DifferentiableQoI.

Definition at line 78 of file H-qoi.C.

80 {
81 
82  FEMContext & c = cast_ref<FEMContext &>(context);
83 
84  // First we get some references to cell-specific data that
85  // will be used to assemble the linear system.
86  FEBase * side_fe = nullptr;
87  c.get_side_fe(0, side_fe);
88 
89  // Element Jacobian * quadrature weights for interior integration
90  const std::vector<Real> & JxW = side_fe->get_JxW();
91  const std::vector<Point> & q_point = side_fe->get_xyz();
92 
93  // Loop over qp's, compute the function at each qp and add
94  // to get the QoI
95  unsigned int n_qpoints = c.get_side_qrule().n_points();
96 
97  Number dQoI_0 = 0.;
98  Number u = 0.;
99  Number C = 0.;
100 
101  // If side is on the left outlet
102  if (c.has_side_boundary_id(2))
103  {
104  //Loop over all the qps on this side
105  for (unsigned int qp=0; qp != n_qpoints; qp++)
106  {
107  Real x = q_point[qp](0);
108 
109  if (x < 0.)
110  {
111  // Get u and C at the qp
112  u = c.side_value(0, qp);
113  C = c.side_value(3, qp);
114 
115  dQoI_0 += JxW[qp] * -u * C;
116  } // end if
117  } // end quadrature loop
118  } // end if on bdry
119 
120  c.get_qois()[0] += dQoI_0;
121 }

References libMesh::FEAbstract::get_JxW(), libMesh::DiffContext::get_qois(), libMesh::FEMContext::get_side_fe(), libMesh::FEMContext::get_side_qrule(), libMesh::FEAbstract::get_xyz(), libMesh::FEMContext::has_side_boundary_id(), libMesh::QBase::n_points(), libMesh::Real, and libMesh::FEMContext::side_value().

◆ side_qoi_derivative()

void CoupledSystemQoI::side_qoi_derivative ( DiffContext ,
const QoISet  
)
virtual

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 from libMesh::DifferentiableQoI.

Definition at line 17 of file H-qoi.C.

19 {
20  FEMContext & c = cast_ref<FEMContext &>(context);
21 
22  FEBase * side_fe = nullptr;
23  c.get_side_fe(0, side_fe);
24 
25  // Element Jacobian * quadrature weights for interior integration
26  const std::vector<Real> & JxW = side_fe->get_JxW();
27 
28  // Side basis functions
29  const std::vector<std::vector<Real>> & phi = side_fe->get_phi();
30  const std::vector<Point > & q_point = side_fe->get_xyz();
31 
32  // The number of local degrees of freedom in each variable
33  const unsigned int n_u_dofs = c.n_dof_indices(1);
34 
37 
38  // Now we will build the element Jacobian and residual.
39  // Constructing the residual requires the solution and its
40  // gradient from the previous timestep. This must be
41  // calculated at each quadrature point by summing the
42  // solution degree-of-freedom values by the appropriate
43  // weight functions.
44  unsigned int n_qpoints = c.get_side_qrule().n_points();
45 
46  Number u = 0.;
47  Number C = 0.;
48 
49  // // If side is on outlet
50  if (c.has_side_boundary_id(2))
51  {
52  // Loop over all the qps on this side
53  for (unsigned int qp=0; qp != n_qpoints; qp++)
54  {
55  Real x = q_point[qp](0);
56 
57  // If side is on left outlet
58  if (x < 0.)
59  {
60  // Get u at the qp
61  u = c.side_value(0, qp);
62  C = c.side_value(3, qp);
63 
64  // Add the contribution from each basis function
65  for (unsigned int i=0; i != n_u_dofs; i++)
66  {
67  Qu(i) += JxW[qp] * -phi[i][qp] * C;
68  QC(i) += JxW[qp] * phi[i][qp] * -u;
69  }
70  } // end if
71 
72  } // end quadrature loop
73 
74  } // end if on outlet
75 }

References libMesh::FEAbstract::get_JxW(), libMesh::FEGenericBase< OutputType >::get_phi(), libMesh::DiffContext::get_qoi_derivatives(), libMesh::FEMContext::get_side_fe(), libMesh::FEMContext::get_side_qrule(), libMesh::FEAbstract::get_xyz(), libMesh::FEMContext::has_side_boundary_id(), libMesh::DiffContext::n_dof_indices(), libMesh::QBase::n_points(), libMesh::Real, and libMesh::FEMContext::side_value().

◆ thread_join()

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 33 of file diff_qoi.C.

36 {
37  for (auto i : index_range(qoi))
38  qoi[i] += other_qoi[i];
39 }

References libMesh::index_range().

Member Data Documentation

◆ assemble_qoi_elements

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 101 of file diff_qoi.h.

◆ assemble_qoi_internal_sides

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 93 of file diff_qoi.h.

◆ assemble_qoi_sides

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 85 of file diff_qoi.h.

Referenced by main().


The documentation for this class was generated from the following files:
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::FEAbstract::get_xyz
const std::vector< Point > & get_xyz() const
Definition: fe_abstract.h:237
CoupledSystemQoI::CoupledSystemQoI
CoupledSystemQoI()
Definition: H-qoi.h:18
libMesh::QBase::n_points
unsigned int n_points() const
Definition: quadrature.h:126
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
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
libMesh::DiffContext::get_qois
const std::vector< Number > & get_qois() const
Const accessor for QoI vector.
Definition: diff_context.h:319
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::DifferentiableQoI
This class provides a specific system class.
Definition: diff_qoi.h:52
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::DenseSubVector< Number >
libMesh::FEMContext::has_side_boundary_id
bool has_side_boundary_id(boundary_id_type id) const
Reports if the boundary id is found on the current side.
Definition: fem_context.C:254
libMesh::FEMContext::side_value
Number side_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:620
libMesh::FEMContext::get_side_qrule
const QBase & get_side_qrule() const
Accessor for element side quadrature rule for the dimension of the current _elem.
Definition: fem_context.h:797
libMesh::FEMContext::get_side_fe
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:312
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