libMesh
Public Member Functions | Public Attributes | List of all members
libMesh::DifferentiableQoI Class Referenceabstract

This class provides a specific system class. More...

#include <diff_qoi.h>

Inheritance diagram for libMesh::DifferentiableQoI:
[legend]

Public Member Functions

 DifferentiableQoI ()
 Constructor. More...
 
virtual ~DifferentiableQoI ()
 Destructor. 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 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 init_context (DiffContext &)
 Prepares the result of a build_context() call for use. More...
 
virtual std::unique_ptr< DifferentiableQoIclone ()=0
 Copy of this object. 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

This class provides a specific system class.

It aims to generalize any system, linear or nonlinear, which provides both a residual and a Jacobian.

This class is part of the new DifferentiableSystem framework, which is still experimental. Users of this framework should beware of bugs and future API changes.

Author
Roy H. Stogner
Date
2006

Definition at line 52 of file diff_qoi.h.

Constructor & Destructor Documentation

◆ DifferentiableQoI()

libMesh::DifferentiableQoI::DifferentiableQoI ( )

Constructor.

Optionally initializes required data structures.

Definition at line 26 of file diff_qoi.C.

26  :
27  assemble_qoi_sides(false),
30 {
31 }

◆ ~DifferentiableQoI()

virtual libMesh::DifferentiableQoI::~DifferentiableQoI ( )
inlinevirtual

Destructor.

Definition at line 65 of file diff_qoi.h.

65 {}

Member Function Documentation

◆ clear_qoi()

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

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> libMesh::DifferentiableQoI::clone ( )
pure virtual

Copy of this object.

User should override to copy any needed state.

Implemented in libMesh::DifferentiableSystem, LaplaceQoI, and CoupledSystemQoI.

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

◆ element_qoi()

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

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  
)
inlinevirtual

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 
)
virtual

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 )
inlinevirtual

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()

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

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().

◆ 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 
)
virtual

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().

◆ side_qoi()

virtual void libMesh::DifferentiableQoI::side_qoi ( DiffContext ,
const QoISet  
)
inlinevirtual

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

135  {}

◆ side_qoi_derivative()

virtual void libMesh::DifferentiableQoI::side_qoi_derivative ( DiffContext ,
const QoISet  
)
inlinevirtual

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

147  {}

◆ thread_join()

void libMesh::DifferentiableQoI::thread_join ( std::vector< Number > &  qoi,
const std::vector< Number > &  other_qoi,
const QoISet qoi_indices 
)
virtual

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

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

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

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::DifferentiableQoI::assemble_qoi_elements
bool assemble_qoi_elements
If assemble_qoi_elements is false (it is true by default), the assembly loop for a quantity of intere...
Definition: diff_qoi.h:101
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::DifferentiableQoI::assemble_qoi_internal_sides
bool assemble_qoi_internal_sides
If assemble_qoi_internal_sides is true (it is false by default), the assembly loop for a quantity of ...
Definition: diff_qoi.h:93
libMesh::DifferentiableQoI::assemble_qoi_sides
bool assemble_qoi_sides
If assemble_qoi_sides is true (it is false by default), the assembly loop for a quantity of interest ...
Definition: diff_qoi.h:85