1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2020 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
17 
18 
19 
20 #ifndef LIBMESH_DIFF_QOI_H
21 #define LIBMESH_DIFF_QOI_H
22 
23 // Local Includes
24 #include "libmesh/diff_context.h"
25 
26 // C++ includes
27 #include <memory>
28 
29 namespace libMesh
30 {
31 
32 // Forward declarations
33 class DiffContext;
34 class QoISet;
35 
36 namespace Parallel {
37   class Communicator;
38 }
39 
40 /**
41  * This class provides a specific system class.  It aims
42  * to generalize any system, linear or nonlinear, which
43  * provides both a residual and a Jacobian.
44  *
45  * This class is part of the new DifferentiableSystem framework,
46  * which is still experimental.  Users of this framework should
47  * beware of bugs and future API changes.
48  *
49  * \author Roy H. Stogner
50  * \date 2006
51  */
52 class DifferentiableQoI
53 {
54 public:
55 
56   /**
57    * Constructor.  Optionally initializes required
58    * data structures.
59    */
60   DifferentiableQoI ();
61 
62   /**
63    * Destructor.
64    */
~DifferentiableQoI()65   virtual ~DifferentiableQoI () {}
66 
67   /**
68    * Initialize system qoi. By default, does nothing in order to maintain backward
69    * compatibility for FEMSystem applications that control qoi.
70    */
init_qoi(std::vector<Number> &)71   virtual void init_qoi( std::vector<Number> & /*sys_qoi*/){}
72 
73   /**
74    * Clear all the data structures associated with
75    * the QoI.
76    */
clear_qoi()77   virtual void clear_qoi () {}
78 
79   /**
80    * If \p assemble_qoi_sides is true (it is false by default), the
81    * assembly loop for a quantity of interest or its derivatives will
82    * loop over domain boundary sides.  To add domain interior sides,
83    * also set assemble_qoi_internal_sides to true.
84    */
85   bool assemble_qoi_sides;
86 
87   /**
88    * If \p assemble_qoi_internal_sides is true (it is false by
89    * default), the assembly loop for a quantity of interest or its
90    * derivatives will loop over element sides which do not fall on
91    * domain boundaries.
92    */
93   bool assemble_qoi_internal_sides;
94 
95   /**
96    * If \p assemble_qoi_elements is false (it is true by default), the
97    * assembly loop for a quantity of interest or its derivatives will
98    * skip computing on mesh elements, and will only compute on mesh
99    * sides.
100    */
101   bool assemble_qoi_elements;
102 
103   /**
104    * Does any work that needs to be done on \p elem in a quantity of
105    * interest assembly loop, outputting to elem_qoi.
106    *
107    * Only qois included in the supplied \p QoISet need to be
108    * assembled.
109    */
element_qoi(DiffContext &,const QoISet &)110   virtual void element_qoi (DiffContext &,
111                             const QoISet &)
112   {}
113 
114   /**
115    * Does any work that needs to be done on \p elem in a quantity of
116    * interest derivative assembly loop, outputting to
117    * elem_qoi_derivative
118    *
119    * Only qois included in the supplied \p QoISet need their
120    * derivatives assembled.
121    */
element_qoi_derivative(DiffContext &,const QoISet &)122   virtual void element_qoi_derivative (DiffContext &,
123                                        const QoISet &)
124   {}
125 
126   /**
127    * Does any work that needs to be done on \p side of \p elem in a
128    * quantity of interest assembly loop, outputting to elem_qoi.
129    *
130    * Only qois included in the supplied \p QoISet need to be
131    * assembled.
132    */
side_qoi(DiffContext &,const QoISet &)133   virtual void side_qoi (DiffContext &,
134                          const QoISet &)
135   {}
136 
137   /**
138    * Does any work that needs to be done on \p side of \p elem in a
139    * quantity of interest derivative assembly loop, outputting to
140    * elem_qoi_derivative.
141    *
142    * Only qois included in the supplied \p QoISet need their
143    * derivatives assembled.
144    */
side_qoi_derivative(DiffContext &,const QoISet &)145   virtual void side_qoi_derivative (DiffContext &,
146                                     const QoISet &)
147   {}
148 
149   /**
150    * Prepares the result of a build_context() call for use.
151    *
152    * FEMSystem-based problems will need to reimplement this in order to
153    * call FE::get_*() as their particular QoI requires.  Trying to
154    * evaluate a QoI without overriding init_context is both
155    * inefficient and deprecated.
156    */
init_context(DiffContext &)157   virtual void init_context(DiffContext &) { libmesh_deprecated(); }
158 
159   /**
160    * Copy of this object. User should override to copy any needed state.
161    */
162   virtual std::unique_ptr<DifferentiableQoI> clone() =0;
163 
164   /**
165    * Method to combine thread-local qois. By default, simply sums thread qois.
166    */
167   virtual void thread_join(std::vector<Number> & qoi,
168                            const std::vector<Number> & other_qoi,
169                            const QoISet & qoi_indices);
170 
171   /**
172    * Method to populate system qoi data structure with process-local qoi. By default, simply
173    * sums process qois into system qoi.
174    */
175   virtual void parallel_op(const Parallel::Communicator & communicator,
176                            std::vector<Number> & sys_qoi,
177                            std::vector<Number> & local_qoi,
178                            const QoISet & qoi_indices);
179 
180   /**
181    * Method to finalize qoi derivatives which require more than just a simple
182    * sum of element contributions.
183    */
184   virtual void finalize_derivative(NumericVector<Number> & derivatives, std::size_t qoi_index);
185 };
186 
187 } // namespace libMesh
188 
189 
190 #endif // LIBMESH_DIFF_QOI_H
191