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_RB_TEMPORAL_DISCRETIZATION_H 21 : #define LIBMESH_RB_TEMPORAL_DISCRETIZATION_H 22 : 23 : // libMesh includes 24 : #include "libmesh/libmesh_common.h" 25 : 26 : #include <vector> 27 : 28 : namespace libMesh 29 : { 30 : 31 : /** 32 : * Define a class that encapsulates the details of a 33 : * "generalized Euler" temporal discretization to be 34 : * used in the rbOOmit framework. 35 : * 36 : * \author David Knezevic 37 : * \date 2012 38 : * \brief Encapsulates the details of the generalized Euler discretization. 39 : */ 40 : class RBTemporalDiscretization 41 : { 42 : public: 43 : 44 : /** 45 : * Constructor. 46 : */ 47 : RBTemporalDiscretization(); 48 : 49 : /** 50 : * Special functions can all be defaulted for this simple class. 51 : */ 52 : RBTemporalDiscretization (RBTemporalDiscretization &&) = default; 53 : RBTemporalDiscretization (const RBTemporalDiscretization &) = default; 54 : RBTemporalDiscretization & operator= (const RBTemporalDiscretization &) = default; 55 : RBTemporalDiscretization & operator= (RBTemporalDiscretization &&) = default; 56 288 : virtual ~RBTemporalDiscretization () = default; 57 : 58 : /** 59 : * Get/set delta_t, the time-step size. 60 : */ 61 : Real get_delta_t() const; 62 : void set_delta_t(const Real delta_t_in); 63 : 64 : /** 65 : * Get/set euler_theta, parameter that determines 66 : * the temporal discretization. 67 : */ 68 : Real get_euler_theta() const; 69 : void set_euler_theta(const Real euler_theta_in); 70 : 71 : /** 72 : * Get/set the current time-step. 73 : */ 74 : unsigned int get_time_step() const; 75 : void set_time_step(const unsigned int k); 76 : 77 : /** 78 : * Get/set the total number of time-steps. 79 : */ 80 : unsigned int get_n_time_steps() const; 81 : void set_n_time_steps(const unsigned int K); 82 : 83 : /** 84 : * Get/set the RHS control. 85 : */ 86 : Real get_control(const unsigned int k) const; 87 : void set_control(const std::vector<Real> & control); 88 : 89 : /** 90 : * Read in and initialize parameters from \p parameters_filename. 91 : */ 92 : void process_temporal_parameters_file (const std::string & parameters_filename); 93 : 94 : /** 95 : * Pull the temporal discretization data from \p other. 96 : */ 97 : void pull_temporal_discretization_data(RBTemporalDiscretization & other); 98 : 99 : private: 100 : 101 : /** 102 : * The time-step size. 103 : */ 104 : Real _delta_t; 105 : 106 : /** 107 : * The parameter that determines the generalized Euler scheme 108 : * discretization that we employ. 109 : * euler_theta = 0 ---> Forward Euler 110 : * euler_theta = 0.5 ---> Crank-Nicolson 111 : * euler_theta = 1 ---> Backward Euler 112 : */ 113 : Real _euler_theta; 114 : 115 : /** 116 : * The current time-step. 117 : */ 118 : unsigned int _current_time_step; 119 : 120 : /** 121 : * The number of time-steps. 122 : */ 123 : unsigned int _n_time_steps; 124 : 125 : /** 126 : * The RHS control (scalar function of time). 127 : * A function h(t) that is used in the RHS as h(t)*f(x,\f$ \mu \f$). 128 : * See Martin Grepl's thesis 129 : */ 130 : std::vector<Real> _control; 131 : }; 132 : 133 : } 134 : 135 : #endif // LIBMESH_RB_TEMPORAL_DISCRETIZATION_H