00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #ifndef EPETRA_EXT_DIAGONAL_TRANSIENT_MODEL_HPP
00030 #define EPETRA_EXT_DIAGONAL_TRANSIENT_MODEL_HPP
00031
00032
00033 #include "EpetraExt_ModelEvaluator.h"
00034 #include "Teuchos_VerboseObject.hpp"
00035 #include "Teuchos_ParameterListAcceptor.hpp"
00036 #include "Teuchos_Array.hpp"
00037
00038
00039 class Epetra_Comm;
00040 class Epetra_CrsGraph;
00041
00042
00043 namespace EpetraExt {
00044
00045
00086 class DiagonalTransientModel
00087 : public ::EpetraExt::ModelEvaluator,
00088 public Teuchos::VerboseObject<DiagonalTransientModel>,
00089 public Teuchos::ParameterListAcceptor
00090 {
00091 public:
00092
00095
00097 DiagonalTransientModel(
00098 Teuchos::RCP<Epetra_Comm> const& epetra_comm
00099 );
00100
00102 Teuchos::RCP<const Epetra_Vector> get_gamma() const;
00103
00105 Teuchos::RCP<const Epetra_Vector>
00106 getExactSolution(
00107 const double t, const Epetra_Vector *coeff_s_p = 0
00108 ) const;
00109
00111 Teuchos::RCP<const Epetra_MultiVector>
00112 getExactSensSolution(
00113 const double t, const Epetra_Vector *coeff_s_p = 0
00114 ) const;
00115
00117
00120
00122 void setParameterList(Teuchos::RCP<Teuchos::ParameterList> const& paramList);
00124 Teuchos::RCP<Teuchos::ParameterList> getNonconstParameterList();
00126 Teuchos::RCP<Teuchos::ParameterList> unsetParameterList();
00128 Teuchos::RCP<const Teuchos::ParameterList> getParameterList() const;
00130 Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
00131
00133
00136
00138 Teuchos::RCP<const Epetra_Map> get_x_map() const;
00140 Teuchos::RCP<const Epetra_Map> get_f_map() const;
00142 Teuchos::RCP<const Epetra_Map> get_p_map(int l) const;
00144 Teuchos::RCP<const Teuchos::Array<std::string> > get_p_names(int l) const;
00146 Teuchos::RCP<const Epetra_Map> get_g_map(int j) const;
00148 Teuchos::RCP<const Epetra_Vector> get_x_init() const;
00150 Teuchos::RCP<const Epetra_Vector> get_x_dot_init() const;
00152 Teuchos::RCP<const Epetra_Vector> get_p_init(int l) const;
00154 Teuchos::RCP<Epetra_Operator> create_W() const;
00156 InArgs createInArgs() const;
00158 OutArgs createOutArgs() const;
00160 void evalModel( const InArgs& inArgs, const OutArgs& outArgs ) const;
00161
00163
00164 public:
00165
00166 enum EGammaFit { GAMMA_FIT_LINEAR, GAMMA_FIT_RANDOM };
00167
00168 private:
00169
00170
00171
00172
00173 typedef Teuchos::Array<double> coeff_s_t;
00174 typedef Teuchos::Array<int> coeff_s_idx_t;
00175 typedef Teuchos::Array<Teuchos::RCP<const Epetra_Map> > RCP_Eptra_Map_Array_t;
00176 typedef Teuchos::Array<Teuchos::RCP<Epetra_Vector> > RCP_Eptra_Vector_Array_t;
00177 typedef Teuchos::Array<Teuchos::RCP<Teuchos::Array<std::string> > > RCP_Array_String_Array_t;
00178
00179
00180
00181
00182
00183 Teuchos::RCP<Teuchos::ParameterList> paramList_;
00184 Teuchos::RCP<Epetra_Comm> epetra_comm_;
00185 Teuchos::RCP<Epetra_Map> epetra_map_;
00186 bool implicit_;
00187 int numElements_;
00188 double gamma_min_;
00189 double gamma_max_;
00190 coeff_s_t coeff_s_;
00191 coeff_s_idx_t coeff_s_idx_;
00192 EGammaFit gamma_fit_;
00193 double x0_;
00194 bool exactSolutionAsResponse_;
00195 Teuchos::RCP<Epetra_Vector> gamma_;
00196 Teuchos::RCP<Epetra_CrsGraph> W_graph_;
00197 int Np_;
00198 int np_;
00199 int Ng_;
00200 RCP_Eptra_Map_Array_t map_p_;
00201 RCP_Array_String_Array_t names_p_;
00202 RCP_Eptra_Map_Array_t map_g_;
00203 RCP_Eptra_Vector_Array_t p_init_;
00204 Teuchos::RCP<Epetra_Vector> x_init_;
00205 Teuchos::RCP<Epetra_Vector> x_dot_init_;
00206
00207 mutable Teuchos::RCP<const Epetra_Vector> coeff_s_p_;
00208
00209 bool isIntialized_;
00210
00211
00212
00213
00214 void initialize();
00215
00216 void set_coeff_s_p(
00217 const Teuchos::RCP<const Epetra_Vector> &coeff_s_p
00218 ) const;
00219
00220 void unset_coeff_s_p() const;
00221
00222 int coeff_s_idx(int i) const
00223 {
00224 return coeff_s_idx_[i];
00225 }
00226
00227 double coeff_s(int i) const
00228 {
00229 return (*coeff_s_p_)[coeff_s_idx(i)];
00230 }
00231
00232 };
00233
00234
00239 Teuchos::RCP<DiagonalTransientModel>
00240 diagonalTransientModel(
00241 Teuchos::RCP<Epetra_Comm> const& epetra_comm,
00242 Teuchos::RCP<Teuchos::ParameterList> const& paramList = Teuchos::null
00243 );
00244
00245
00246 }
00247
00248
00249
00250
00251
00252
00253
00254 #endif // EPETRA_EXT_DIAGONAL_TRANSIENT_MODEL_HPP