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

Constructor & Destructor Documentation

libMesh::DifferentiableQoI::DifferentiableQoI ( )

Constructor.

Optionally initializes required data structures.

Definition at line 24 of file diff_qoi.C.

24  :
25  assemble_qoi_sides(false),
28 {
29 }
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:98
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:90
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:82
virtual libMesh::DifferentiableQoI::~DifferentiableQoI ( )
virtual

Destructor.

Definition at line 62 of file diff_qoi.h.

62 {}

Member Function Documentation

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

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> 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(), and init_context().

virtual void libMesh::DifferentiableQoI::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 in LaplaceQoI.

Definition at line 107 of file diff_qoi.h.

109  {}
virtual void libMesh::DifferentiableQoI::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 in HeatSystem, LaplaceSystem, LaplaceSystem, and LaplaceQoI.

Definition at line 119 of file diff_qoi.h.

121  {}
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 51 of file diff_qoi.C.

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

52 {
53  // by default, do nothing
54 }
virtual void libMesh::DifferentiableQoI::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 in libMesh::FEMSystem, HeatSystem, CoupledSystem, ElasticitySystem, LaplaceSystem, CurlCurlSystem, LaplaceSystem, PoissonSystem, LaplaceSystem, LaplaceSystem, CurlCurlSystem, L2System, SolidSystem, NavierSystem, HeatSystem, and LaplaceQoI.

Definition at line 152 of file diff_qoi.h.

References clone(), finalize_derivative(), parallel_op(), and thread_join().

152 {}
virtual void libMesh::DifferentiableQoI::init_qoi ( std::vector< Number > &  )
virtual

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

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

Definition at line 130 of file diff_qoi.h.

132  {}
virtual void libMesh::DifferentiableQoI::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 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 
)
virtual

Method to combine thread-local qois.

By default, simply sums thread qois.

Definition at line 31 of file diff_qoi.C.

Referenced by 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

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

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

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: