Line data Source code
1 : // rbOOmit: An implementation of the Certified Reduced Basis method. 2 : // Copyright (C) 2009, 2010 David J. Knezevic 3 : 4 : // This file is part of rbOOmit. 5 : 6 : // rbOOmit is free software; you can redistribute it and/or 7 : // modify it under the terms of the GNU Lesser General Public 8 : // License as published by the Free Software Foundation; either 9 : // version 2.1 of the License, or (at your option) any later version. 10 : 11 : // rbOOmit is distributed in the hope that it will be useful, 12 : // but WITHOUT ANY WARRANTY; without even the implied warranty of 13 : // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 14 : // Lesser General Public License for more details. 15 : 16 : // You should have received a copy of the GNU Lesser General Public 17 : // License along with this library; if not, write to the Free Software 18 : // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 19 : 20 : #ifndef LIBMESH_TRANSIENT_RB_EVALUATION_H 21 : #define LIBMESH_TRANSIENT_RB_EVALUATION_H 22 : 23 : // rbOOmit includes 24 : #include "libmesh/rb_evaluation.h" 25 : #include "libmesh/rb_temporal_discretization.h" 26 : 27 : // libMesh includes 28 : 29 : // C++ includes 30 : 31 : namespace libMesh 32 : { 33 : 34 : class TransientRBThetaExpansion; 35 : 36 : /** 37 : * This class is part of the rbOOmit framework. 38 : * 39 : * TransientRBEvaluation extends RBEvaluation to 40 : * encapsulate the code and data required 41 : * to perform "online" RB evaluations for 42 : * Linear Time Invariant (LTI) transient problems. 43 : * 44 : * We can handle time controls on the RHS as h(t)*f(x,\f$ \mu \f$). 45 : * See Martin Grepl's thesis for more details. 46 : * 47 : * \author David J. Knezevic 48 : * \date 2011 49 : */ 50 32 : class TransientRBEvaluation : public RBEvaluation, public RBTemporalDiscretization 51 : { 52 : public: 53 : 54 : /** 55 : * Constructor. 56 : */ 57 : TransientRBEvaluation (const Parallel::Communicator & comm_in); 58 : 59 : /** 60 : * Special functions. 61 : * - This class contains unique_ptrs, so it can't be default copy 62 : * constructed/assigned. 63 : * - The destructor is defaulted out of line. 64 : */ 65 : TransientRBEvaluation (TransientRBEvaluation &&) = default; 66 : TransientRBEvaluation (const TransientRBEvaluation &) = delete; 67 : TransientRBEvaluation & operator= (const TransientRBEvaluation &) = delete; 68 : TransientRBEvaluation & operator= (TransientRBEvaluation &&) = default; 69 : virtual ~TransientRBEvaluation (); 70 : 71 : /** 72 : * The type of the parent. 73 : */ 74 : typedef RBEvaluation Parent; 75 : 76 : /** 77 : * Clear this TransientRBEvaluation object. 78 : * Override to also clear the M_q representors 79 : */ 80 : virtual void clear() override; 81 : 82 : /** 83 : * Resize and clear the data vectors corresponding to the 84 : * value of \p Nmax. Optionally resize the data structures 85 : * required for the error bound. 86 : * Overridden to resize data relevant in the time-dependent 87 : * case. 88 : */ 89 : virtual void resize_data_structures(const unsigned int Nmax, 90 : bool resize_error_bound_data=true) override; 91 : 92 : /** 93 : * Perform online solve for current_params 94 : * with the N basis functions. Overridden 95 : * to perform a time-dependent solve. 96 : */ 97 : using RBEvaluation::rb_solve; 98 : virtual Real rb_solve(unsigned int N) override; 99 : 100 : /** 101 : * If a solve has already been performed, then we cached some data 102 : * and we can perform a new solve much more rapidly 103 : * (with the same parameters but a possibly different initial condition/rhs control). 104 : */ 105 : virtual Real rb_solve_again(); 106 : 107 : /** 108 : * \returns A scaling factor that we can use to provide a consistent 109 : * scaling of the RB error bound across different parameter values. 110 : */ 111 : virtual Real get_error_bound_normalization() override; 112 : 113 : /** 114 : * Specifies the residual scaling on the numerator to 115 : * be used in the a posteriori error bound. Override 116 : * in subclass in order to obtain the desired error bound. 117 : */ 118 : virtual Real residual_scaling_numer(Real alpha_LB); 119 : 120 : /** 121 : * Compute the dual norm of the residual for the solution 122 : * saved in RB_solution. This function uses the cached time-independent 123 : * data. 124 : */ 125 : using RBEvaluation::compute_residual_dual_norm; 126 : virtual Real compute_residual_dual_norm(const unsigned int N) override; 127 : 128 : /** 129 : * Compute the dual norm of the residual for the solution 130 : * saved in RB_solution. This function does not used the cached 131 : * data and therefore also works when the parameter changes as 132 : * a function of time. 133 : */ 134 : virtual Real uncached_compute_residual_dual_norm(const unsigned int N); 135 : 136 : /** 137 : * Helper function for caching the terms in the 138 : * online residual assembly that do not change in time. 139 : * (This is only useful when the parameter is fixed in time.) 140 : */ 141 : void cache_online_residual_terms(const unsigned int N); 142 : 143 : /** 144 : * Clear all the Riesz representors that are used to compute the RB residual 145 : * (and hence error bound). This is useful since once we complete the Greedy 146 : * we may not need the representors any more. 147 : * Override to clear the M_q representors. 148 : */ 149 : virtual void clear_riesz_representors() override; 150 : 151 : /** 152 : * Write out all the data to text files in order to segregate the 153 : * Offline stage from the Online stage. 154 : * 155 : * \note This is a legacy method, use RBDataSerialization instead. 156 : */ 157 : virtual void legacy_write_offline_data_to_files(const std::string & directory_name = "offline_data", 158 : const bool write_binary_data=true) override; 159 : 160 : /** 161 : * Read in the saved Offline reduced basis data 162 : * to initialize the system for Online solves. 163 : * 164 : * \note This is a legacy method, use RBDataSerialization instead. 165 : */ 166 : virtual void legacy_read_offline_data_from_files(const std::string & directory_name = "offline_data", 167 : bool read_error_bound_data=true, 168 : const bool read_binary_data=true) override; 169 : 170 : //----------- PUBLIC DATA MEMBERS -----------// 171 : 172 : /** 173 : * Dense RB L2 matrix. 174 : */ 175 : DenseMatrix<Number> RB_L2_matrix; 176 : 177 : /** 178 : * Cached data for subsequent solves. 179 : */ 180 : DenseMatrix<Number> RB_LHS_matrix; 181 : DenseMatrix<Number> RB_RHS_matrix; 182 : DenseVector<Number> RB_RHS_save; 183 : 184 : /** 185 : * Dense matrices for the RB mass matrices. 186 : */ 187 : std::vector<DenseMatrix<Number>> RB_M_q_vector; 188 : 189 : /** 190 : * The RB outputs for all time-levels from the 191 : * most recent rb_solve. 192 : */ 193 : std::vector<std::vector<Number>> RB_outputs_all_k; 194 : 195 : /** 196 : * The error bounds for each RB output for all 197 : * time-levels from the most recent rb_solve. 198 : */ 199 : std::vector<std::vector<Real>> RB_output_error_bounds_all_k; 200 : 201 : /** 202 : * The RB solution at the previous time-level. 203 : */ 204 : DenseVector<Number> old_RB_solution; 205 : 206 : /** 207 : * Array storing the solution data at each time level from the most recent solve. 208 : */ 209 : std::vector<DenseVector<Number>> RB_temporal_solution_data; 210 : 211 : /** 212 : * The error bound data for all time-levels from the 213 : * most recent rb_solve. 214 : */ 215 : std::vector<Real > error_bound_all_k; 216 : 217 : /** 218 : * Vector storing initial L2 error for all 219 : * 1 <= N <= RB_size. 220 : */ 221 : std::vector<Real> initial_L2_error_all_N; 222 : 223 : /** 224 : * The RB initial conditions (i.e. L2 projection of the truth 225 : * initial condition) for each N. 226 : */ 227 : std::vector<DenseVector<Number>> RB_initial_condition_all_N; 228 : 229 : /** 230 : * Vectors storing the residual representor inner products 231 : * to be used in computing the residuals online. 232 : */ 233 : std::vector<std::vector<std::vector<Number>>> Fq_Mq_representor_innerprods; 234 : std::vector<std::vector<std::vector<Number>>> Mq_Mq_representor_innerprods; 235 : std::vector<std::vector<std::vector<std::vector<Number>>>> Aq_Mq_representor_innerprods; 236 : 237 : 238 : /** 239 : * Cached residual terms. These can be used to accelerate residual calculations 240 : * when we have an LTI system. 241 : */ 242 : Number cached_Fq_term; 243 : DenseVector<Number> cached_Fq_Aq_vector; 244 : DenseMatrix<Number> cached_Aq_Aq_matrix; 245 : DenseVector<Number> cached_Fq_Mq_vector; 246 : DenseMatrix<Number> cached_Aq_Mq_matrix; 247 : DenseMatrix<Number> cached_Mq_Mq_matrix; 248 : 249 : /** 250 : * Vector storing the mass matrix representors. 251 : * These are basis dependent and hence stored here. 252 : */ 253 : std::vector<std::vector<std::unique_ptr<NumericVector<Number>>>> M_q_representor; 254 : 255 : /** 256 : * Check that the data has been cached in case of using rb_solve_again 257 : */ 258 : bool _rb_solve_data_cached; 259 : }; 260 : 261 : } 262 : 263 : #endif // LIBMESH_TRANSIENT_RB_EVALUATION_H