RBGen  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
RBGen_IncSVDPOD.h
1 #ifndef RBGEN_INCSVD_POD_H
2 #define RBGEN_INCSVD_POD_H
3 
4 #include "RBGen_PODMethod.hpp"
5 #include "RBGen_Method.hpp"
6 #include "RBGen_Filter.hpp"
7 #include "AnasaziEpetraAdapter.hpp"
8 #include "AnasaziOrthoManager.hpp"
9 #include "Epetra_MultiVector.h"
10 #include "Epetra_Operator.h"
11 #include "Epetra_SerialDenseMatrix.h"
13 #include "Teuchos_Time.hpp"
14 
15 
16 //
17 // computeBasis()
18 // |
19 // makePass() ___
20 // while () { / expand()
21 // incStep() ---- SVD()
22 // } \___ shrink()
23 //
24 // makePass(), expand() and shrink() are pure virtual
25 //
26 // makePass() is implemented in a base class that decides
27 // if a method is Multipass, and if so, in what manner
28 // makePass() puts the data for the next pass into the proper columns
29 // of U_, then calls incstep (telling it how many new columns)
30 //
31 // incstep calls expand to construct expanded U_ and V_, and B_
32 // it computes the SVD of B_
33 // it passes this data to shrink(), which shrink according to the
34 // base class.
35 //
36 // updateBasis() allows the current factorization to be updated
37 // by appending new snapshots. This allows a truly incremental
38 // factorization. It is also pure virtual, and is implemented
39 // along side makePass(). It may not be implemented in all
40 // IncSVDPOD methods; in fact, it is currently implemented
41 // only by derived classes of ISVDSingle
42 //
43 // expand(),shrink() are implemented in a base class that
44 // decides the representation: UDV, QRW, QBW
45 //
46 // IncSVDPOD
47 // |
48 // -------------------------------------------------------------------
49 // | | | | | | |
50 // ISVDUDV ISVDQRW ISVDQBW ISVDMultiCD ISVDMultiSDA ISVDMultiSDB ISVDSingle
51 // | | | | | | |
52 // ------------------ --------------------------------------
53 // \ /
54 // \ /
55 // \ /
56 // \ /
57 // \ /
58 // \--- Concrete Base Class ----/
59 //
60 // Then a concrete base class (one of 3x4==12 varieties) is formed through
61 // inheritence. This is the Template Pattern type of Design Pattern.
62 //
63 
64 namespace RBGen {
65 
67  class IncSVDPOD : public virtual Method<Epetra_MultiVector,Epetra_Operator>, public virtual PODMethod<double> {
68 
69  public:
71 
72 
74  IncSVDPOD();
75 
77  virtual ~IncSVDPOD() {};
79 
81 
82 
84  void computeBasis();
85 
87  virtual void updateBasis( const Teuchos::RCP< Epetra_MultiVector >& update_ss ) = 0;
88 
90 
92 
93 
96 
99 
101  std::vector<double> getSingularValues() const;
102 
104  double getCompTime() const { return timerComp_.totalElapsedTime(); }
105 
107  const std::vector<double> & getResNorms();
108 
110 
112 
113 
117  const Teuchos::RCP< RBGen::FileIOHandler< Epetra_Operator > >& fileio = Teuchos::null );
118 
119  void Reset( const Teuchos::RCP<Epetra_MultiVector>& new_ss );
120 
122 
124 
125 
126  bool isInitialized() { return isInitialized_; }
127 
129 
130  protected:
131 
132  // private member for performing inc steps
133  void incStep(int lup);
134  virtual void makePass() = 0;
135  virtual void expand(const int lup) = 0;
136  virtual void shrink(const int down, std::vector<double> &S, Epetra_SerialDenseMatrix &U, Epetra_SerialDenseMatrix &V) = 0;
137 
138  // Is this object initialized?
139  bool isInitialized_;
140 
141  // Singular value filter
143 
144  // Max allocation size.
145  // The maximum rank DURING any step:
146  // lup <= maxBasisSize_ - curRank_
147  int maxBasisSize_;
148 
149  // Current rank of the factorization
150  int curRank_;
151 
152  // Pointers to the snapshots and reduced basis.
155 
156  // SerialDenseMatrix holding current core matrix B
158 
159  // Vector holding singular values.
160  std::vector<double> sigma_;
161 
162  // Number of snapshots processed thus far
163  int numProc_;
164 
165  // Maximum allowable number of passes through A
166  int maxNumPasses_;
167 
168  // Current number of passes through A
169  int curNumPasses_;
170 
171  // Convergence tolerance
172  double tol_;
173 
174  // min,max number of update vectors
175  int lmin_;
176  int lmax_;
177  int startRank_;
178 
179  // ortho manager
181 
182  // cummulative timer
183  Teuchos::Time timerComp_;
184 
185  // debug flag
186  bool debug_;
187 
188  // verb level
189  int verbLevel_;
190 
191  // residual norms
192  std::vector<double> resNorms_;
193 
194  // maximum number of data set columns
195  // i.e., maximum number of rows in V
196  int maxNumCols_;
197 
198  // is V locally replicated or globally distributed?
199  bool Vlocal_;
200  };
201 
202 } // end of RBGen namespace
203 
204 #endif // RBGEN_INCSVD_POD_H
void Reset(const Teuchos::RCP< Epetra_MultiVector > &new_ss)
Reset the snapshot set used to compute the reduced basis.
std::vector< double > getSingularValues() const
Return the singular values.
double getCompTime() const
Return the cummulative wall-clock time.
Teuchos::RCP< const Epetra_MultiVector > getRightBasis() const
Return a basis for the right singular subspace.
Teuchos::RCP< const Epetra_MultiVector > getBasis() const
Return a basis for the left singular subspace.
void computeBasis()
Computes bases for the left and (optionally) right singular subspaces, along with singular vaues...
virtual void updateBasis(const Teuchos::RCP< Epetra_MultiVector > &update_ss)=0
Update the current basis by appending new snapshots.
Abstract base class for reduced basis methods.
Abstract base class for reduced basis POD methods.
void Initialize(const Teuchos::RCP< Teuchos::ParameterList > &params, const Teuchos::RCP< const Epetra_MultiVector > &init, const Teuchos::RCP< RBGen::FileIOHandler< Epetra_Operator > > &fileio=Teuchos::null)
Initialize the method with the given parameter list and snapshot set.
const std::vector< double > & getResNorms()
Return the scaled residual norms.
IncSVDPOD()
Default constructor.
Class for producing a basis using the Incremental SVD.
double totalElapsedTime(bool readCurrentTime=false) const
virtual ~IncSVDPOD()
Destructor.