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 ()=default
 Destructor. 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 init_qoi_count (System &)
 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 }
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:115
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:107
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:99

◆ ~DifferentiableQoI()

virtual libMesh::DifferentiableQoI::~DifferentiableQoI ( )
virtualdefault

Destructor.

Member Function Documentation

◆ clear_qoi()

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

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

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

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

◆ init_context()

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

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

Definition at line 171 of file diff_qoi.h.

171 { libmesh_deprecated(); }

◆ init_qoi() [1/2]

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

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 > &  )
inline

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

virtual void libMesh::DifferentiableQoI::init_qoi_count ( System )
inlinevirtual

Initialize system qoi.

Often this will just call sys.init_qois(some_desired_number_of_qois)

Reimplemented in LaplaceQoI, and CoupledSystemQoI.

Definition at line 85 of file diff_qoi.h.

85 {}

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

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

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

149  {}

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

161  {}

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

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

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

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

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


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