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

#include <H-qoi.h>

Inheritance diagram for CoupledSystemQoI:
[legend]

Public Member Functions

 CoupledSystemQoI ()=default
 
 CoupledSystemQoI (const CoupledSystemQoI &)=default
 
virtual ~CoupledSystemQoI ()=default
 
virtual void init_qoi_count (System &sys)
 Initialize system qoi. More...
 
virtual void init_context (DiffContext &)
 Prepares the result of a build_context() call for use. More...
 
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...
 
void init_qoi (std::vector< Number > &)
 Non-virtual, to try to help deprecated user code catch this change at compile time (if they specified override) 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 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...
 

Protected Attributes

unsigned int u_var
 
unsigned int p_var
 

Detailed Description

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

Constructor & Destructor Documentation

◆ CoupledSystemQoI() [1/2]

CoupledSystemQoI::CoupledSystemQoI ( )
default

◆ CoupledSystemQoI() [2/2]

CoupledSystemQoI::CoupledSystemQoI ( const CoupledSystemQoI )
default

◆ ~CoupledSystemQoI()

virtual CoupledSystemQoI::~CoupledSystemQoI ( )
virtualdefault

Member Function Documentation

◆ clear_qoi()

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

Clear all the data structures associated with the QoI.

Definition at line 91 of file diff_qoi.h.

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

91 {}

◆ 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 33 of file H-qoi.h.

34  {
35  // Calls defaulted copy constructor for this class
36  return std::make_unique<CoupledSystemQoI>(*this);
37  }

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

Definition at line 124 of file diff_qoi.h.

126  {}

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

Definition at line 136 of file diff_qoi.h.

138  {}

◆ 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.

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

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

◆ init_context()

void CoupledSystemQoI::init_context ( DiffContext )
virtual

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

FEMSystem-based problems will need to reimplement this in order to call FE::get_*() as their particular QoI requires. Trying to evaluate a QoI without overriding init_context is both inefficient and deprecated.

Reimplemented from libMesh::DifferentiableQoI.

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

References libMesh::FEMContext::get_element_fe(), libMesh::FEAbstract::get_JxW(), libMesh::FEAbstract::get_nothing(), libMesh::FEGenericBase< OutputType >::get_phi(), libMesh::FEMContext::get_side_fe(), libMesh::DiffContext::get_system(), libMesh::FEAbstract::get_xyz(), and libMesh::System::variable_number().

17 {
18  FEMContext & c = cast_ref<FEMContext &>(context);
19 
20  u_var = c.get_system().variable_number("u");
21  p_var = c.get_system().variable_number("p");
22 
23  // Note that the concentration and velocity components
24  // use the same basis.
25  FEBase * u_elem_fe = nullptr;
26  c.get_element_fe(u_var, u_elem_fe);
27  u_elem_fe->get_nothing();
28 
29  FEBase * p_elem_fe = nullptr;
30  c.get_element_fe(p_var, p_elem_fe);
31  p_elem_fe->get_nothing();
32 
33  FEBase * side_fe = nullptr;
34  c.get_side_fe(u_var, side_fe);
35 
36  side_fe->get_JxW();
37  side_fe->get_phi();
38  side_fe->get_xyz();
39 
40  // Don't waste time on side computations for p
41  FEBase * p_side_fe = nullptr;
42  c.get_side_fe(p_var, p_side_fe);
43  p_side_fe->get_nothing();
44 }
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
unsigned int u_var
Definition: H-qoi.h:40
unsigned int p_var
Definition: H-qoi.h:40
unsigned int variable_number(std::string_view var) const
Definition: system.C:1557
const System & get_system() const
Accessor for associated system.
Definition: diff_context.h:106
virtual_for_inffe const std::vector< Real > & get_JxW() const
Definition: fe_abstract.h:295
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62
virtual_for_inffe const std::vector< Point > & get_xyz() const
Definition: fe_abstract.h:272
void get_nothing() const
Definition: fe_abstract.h:261
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:277
This class forms the foundation from which generic finite elements may be derived.
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:207

◆ init_qoi() [1/2]

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

Initialize system qoi.

This version of the function required direct vector access, and is now deprecated.

Definition at line 72 of file diff_qoi.h.

72 {}

◆ init_qoi() [2/2]

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

Non-virtual, to try to help deprecated user code catch this change at compile time (if they specified override)

Definition at line 78 of file diff_qoi.h.

78 {}

◆ init_qoi_count()

void CoupledSystemQoI::init_qoi_count ( System )
virtual

Initialize system qoi.

Often this will just call sys.init_qois(some_desired_number_of_qois)

Reimplemented from libMesh::DifferentiableQoI.

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

References libMesh::System::init_qois().

11 {
12  //Only 1 qoi to worry about
13  sys.init_qois(1);
14 }

◆ 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.

References communicator.

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

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 }
DIE A HORRIBLE DEATH HERE typedef MPI_Comm communicator

◆ postprocess()

virtual void CoupledSystemQoI::postprocess ( void  )
inlinevirtual

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

25 {}

◆ 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 110 of file H-qoi.C.

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

112 {
113 
114  FEMContext & c = cast_ref<FEMContext &>(context);
115 
116  // First we get some references to cell-specific data that
117  // will be used to assemble the linear system.
118  FEBase * side_fe = nullptr;
119  c.get_side_fe(u_var, side_fe);
120 
121  // Element Jacobian * quadrature weights for interior integration
122  const std::vector<Real> & JxW = side_fe->get_JxW();
123  const std::vector<Point> & q_point = side_fe->get_xyz();
124 
125  // Loop over qp's, compute the function at each qp and add
126  // to get the QoI
127  unsigned int n_qpoints = c.get_side_qrule().n_points();
128 
129  Number dQoI_0 = 0.;
130  Number u = 0.;
131  Number C = 0.;
132 
133  // If side is on the left outlet
134  if (c.has_side_boundary_id(2))
135  {
136  //Loop over all the qps on this side
137  for (unsigned int qp=0; qp != n_qpoints; qp++)
138  {
139  Real x = q_point[qp](0);
140 
141  if (x < 0.)
142  {
143  // Get u and C at the qp
144  u = c.side_value(0, qp);
145  C = c.side_value(3, qp);
146 
147  dQoI_0 += JxW[qp] * -u * C;
148  } // end if
149  } // end quadrature loop
150  } // end if on bdry
151 
152  c.get_qois()[0] += dQoI_0;
153 }
const std::vector< Number > & get_qois() const
Const accessor for QoI vector.
Definition: diff_context.h:317
Number side_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:656
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
unsigned int u_var
Definition: H-qoi.h:40
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:298
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62
unsigned int n_points() const
Definition: quadrature.h:123
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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

◆ 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 49 of file H-qoi.C.

References libMesh::DiffContext::get_qoi_derivatives(), libMesh::FEMContext::get_side_fe(), libMesh::FEMContext::get_side_qrule(), libMesh::FEMContext::has_side_boundary_id(), libMesh::DiffContext::n_dof_indices(), libMesh::QBase::n_points(), libMesh::Real, and libMesh::FEMContext::side_value().

51 {
52  FEMContext & c = cast_ref<FEMContext &>(context);
53 
54  FEBase * side_fe = nullptr;
55  c.get_side_fe(u_var, side_fe);
56 
57  // Element Jacobian * quadrature weights for interior integration
58  const std::vector<Real> & JxW = side_fe->get_JxW();
59 
60  // Side basis functions
61  const std::vector<std::vector<Real>> & phi = side_fe->get_phi();
62  const std::vector<Point > & q_point = side_fe->get_xyz();
63 
64  // The number of local degrees of freedom in each variable
65  const unsigned int n_u_dofs = c.n_dof_indices(1);
66 
69 
70  // Now we will build the element Jacobian and residual.
71  // Constructing the residual requires the solution and its
72  // gradient from the previous timestep. This must be
73  // calculated at each quadrature point by summing the
74  // solution degree-of-freedom values by the appropriate
75  // weight functions.
76  unsigned int n_qpoints = c.get_side_qrule().n_points();
77 
78  Number u = 0.;
79  Number C = 0.;
80 
81  // // If side is on outlet
82  if (c.has_side_boundary_id(2))
83  {
84  // Loop over all the qps on this side
85  for (unsigned int qp=0; qp != n_qpoints; qp++)
86  {
87  Real x = q_point[qp](0);
88 
89  // If side is on left outlet
90  if (x < 0.)
91  {
92  // Get u at the qp
93  u = c.side_value(0, qp);
94  C = c.side_value(3, qp);
95 
96  // Add the contribution from each basis function
97  for (unsigned int i=0; i != n_u_dofs; i++)
98  {
99  Qu(i) += JxW[qp] * -phi[i][qp] * C;
100  QC(i) += JxW[qp] * phi[i][qp] * -u;
101  }
102  } // end if
103 
104  } // end quadrature loop
105 
106  } // end if on outlet
107 }
Number side_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:656
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
unsigned int n_dof_indices() const
Total number of dof indices on the element.
Definition: diff_context.h:395
unsigned int u_var
Definition: H-qoi.h:40
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:298
Defines a dense subvector for use in finite element computations.
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62
unsigned int n_points() const
Definition: quadrature.h:123
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const std::vector< DenseVector< Number > > & get_qoi_derivatives() const
Const accessor for QoI derivatives.
Definition: diff_context.h:329
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

◆ 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.

References libMesh::index_range().

36 {
37  for (auto i : index_range(qoi))
38  qoi[i] += other_qoi[i];
39 }
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:111

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 115 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 107 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 99 of file diff_qoi.h.

Referenced by main().

◆ p_var

unsigned int CoupledSystemQoI::p_var
protected

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

◆ u_var

unsigned int CoupledSystemQoI::u_var
protected

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


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