Line data Source code
1 : // The libMesh Finite Element Library. 2 : // Copyright (C) 2002-2025 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 : */ 65 828 : virtual ~DifferentiableQoI () = default; 66 : 67 : #ifdef LIBMESH_ENABLE_DEPRECATED 68 : /** 69 : * Initialize system qoi. This version of the function required 70 : * direct vector access, and is now deprecated. 71 : */ 72 137 : virtual void init_qoi( std::vector<Number> & /*sys_qoi*/){} 73 : #else 74 : /** 75 : * Non-virtual, to try to help deprecated user code catch this 76 : * change at compile time (if they specified override) 77 : */ 78 : void init_qoi( std::vector<Number> & /*sys_qoi*/){} 79 : #endif 80 : 81 : /** 82 : * Initialize system qoi. Often this will just call 83 : * sys.init_qois(some_desired_number_of_qois) 84 : */ 85 0 : virtual void init_qoi_count( System & /*sys*/){} 86 : 87 : /** 88 : * Clear all the data structures associated with 89 : * the QoI. 90 : */ 91 0 : virtual void clear_qoi () {} 92 : 93 : /** 94 : * If \p assemble_qoi_sides is true (it is false by default), the 95 : * assembly loop for a quantity of interest or its derivatives will 96 : * loop over domain boundary sides. To add domain interior sides, 97 : * also set assemble_qoi_internal_sides to true. 98 : */ 99 : bool assemble_qoi_sides; 100 : 101 : /** 102 : * If \p assemble_qoi_internal_sides is true (it is false by 103 : * default), the assembly loop for a quantity of interest or its 104 : * derivatives will loop over element sides which do not fall on 105 : * domain boundaries. 106 : */ 107 : bool assemble_qoi_internal_sides; 108 : 109 : /** 110 : * If \p assemble_qoi_elements is false (it is true by default), the 111 : * assembly loop for a quantity of interest or its derivatives will 112 : * skip computing on mesh elements, and will only compute on mesh 113 : * sides. 114 : */ 115 : bool assemble_qoi_elements; 116 : 117 : /** 118 : * Does any work that needs to be done on \p elem in a quantity of 119 : * interest assembly loop, outputting to elem_qoi. 120 : * 121 : * Only qois included in the supplied \p QoISet need to be 122 : * assembled. 123 : */ 124 11799030 : virtual void element_qoi (DiffContext &, 125 : const QoISet &) 126 11799030 : {} 127 : 128 : /** 129 : * Does any work that needs to be done on \p elem in a quantity of 130 : * interest derivative assembly loop, outputting to 131 : * elem_qoi_derivative 132 : * 133 : * Only qois included in the supplied \p QoISet need their 134 : * derivatives assembled. 135 : */ 136 64583 : virtual void element_qoi_derivative (DiffContext &, 137 : const QoISet &) 138 64583 : {} 139 : 140 : /** 141 : * Does any work that needs to be done on \p side of \p elem in a 142 : * quantity of interest assembly loop, outputting to elem_qoi. 143 : * 144 : * Only qois included in the supplied \p QoISet need to be 145 : * assembled. 146 : */ 147 0 : virtual void side_qoi (DiffContext &, 148 : const QoISet &) 149 0 : {} 150 : 151 : /** 152 : * Does any work that needs to be done on \p side of \p elem in a 153 : * quantity of interest derivative assembly loop, outputting to 154 : * elem_qoi_derivative. 155 : * 156 : * Only qois included in the supplied \p QoISet need their 157 : * derivatives assembled. 158 : */ 159 6160 : virtual void side_qoi_derivative (DiffContext &, 160 : const QoISet &) 161 6160 : {} 162 : 163 : /** 164 : * Prepares the result of a build_context() call for use. 165 : * 166 : * FEMSystem-based problems will need to reimplement this in order to 167 : * call FE::get_*() as their particular QoI requires. Trying to 168 : * evaluate a QoI without overriding init_context is both 169 : * inefficient and deprecated. 170 : */ 171 0 : virtual void init_context(DiffContext &) { libmesh_deprecated(); } 172 : 173 : /** 174 : * Copy of this object. User should override to copy any needed state. 175 : */ 176 : virtual std::unique_ptr<DifferentiableQoI> clone() =0; 177 : 178 : /** 179 : * Method to combine thread-local qois. By default, simply sums thread qois. 180 : */ 181 : virtual void thread_join(std::vector<Number> & qoi, 182 : const std::vector<Number> & other_qoi, 183 : const QoISet & qoi_indices); 184 : 185 : /** 186 : * Method to populate system qoi data structure with process-local qoi. By default, simply 187 : * sums process qois into system qoi. 188 : */ 189 : virtual void parallel_op(const Parallel::Communicator & communicator, 190 : std::vector<Number> & sys_qoi, 191 : std::vector<Number> & local_qoi, 192 : const QoISet & qoi_indices); 193 : 194 : /** 195 : * Method to finalize qoi derivatives which require more than just a simple 196 : * sum of element contributions. 197 : */ 198 : virtual void finalize_derivative(NumericVector<Number> & derivatives, std::size_t qoi_index); 199 : }; 200 : 201 : } // namespace libMesh 202 : 203 : 204 : #endif // LIBMESH_DIFF_QOI_H