Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Epetra_SerialDenseSolver.h
Go to the documentation of this file.
1 /*
2 //@HEADER
3 // ************************************************************************
4 //
5 // Epetra: Linear Algebra Services Package
6 // Copyright 2011 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39 //
40 // ************************************************************************
41 //@HEADER
42 */
43 
44 #ifndef EPETRA_SERIALDENSESOLVER_H
45 #define EPETRA_SERIALDENSESOLVER_H
47 #include "Epetra_Object.h"
48 #include "Epetra_CompObject.h"
49 #include "Epetra_BLAS.h"
50 #include "Epetra_LAPACK.h"
51 
52 
54 
131 //=========================================================================
132 class EPETRA_LIB_DLL_EXPORT Epetra_SerialDenseSolver :
133  public Epetra_CompObject, public Epetra_BLAS,
134  public Epetra_LAPACK, public Epetra_Object {
135  public:
136 
138 
141 
142 
144  virtual ~Epetra_SerialDenseSolver();
146 
148 
149 
151  int SetMatrix(Epetra_SerialDenseMatrix & A);
152 
154 
157  int SetVectors(Epetra_SerialDenseMatrix & X, Epetra_SerialDenseMatrix & B);
159 
161 
162 
164 
166  void FactorWithEquilibration(bool Flag) {Equilibrate_ = Flag; return;};
167 
169  void SolveWithTranspose(bool Flag) {Transpose_ = Flag; if (Flag) TRANS_ = 'T'; else TRANS_ = 'N'; return;};
170 
172  void SolveToRefinedSolution(bool Flag) {RefineSolution_ = Flag; return;};
173 
175 
178  void EstimateSolutionErrors(bool Flag) ;
180 
182 
183 
185 
188  virtual int Factor(void);
189 
191 
194  virtual int Solve(void);
195 
197 
200  virtual int Invert(void);
201 
203 
206  virtual int ComputeEquilibrateScaling(void);
207 
209 
212  virtual int EquilibrateMatrix(void);
213 
215 
218  int EquilibrateRHS(void);
219 
220 
222 
225  virtual int ApplyRefinement(void);
226 
228 
231  int UnequilibrateLHS(void);
232 
234 
240  virtual int ReciprocalConditionEstimate(double & Value);
242 
244 
245 
247  bool Transpose() {return(Transpose_);};
248 
250  bool Factored() {return(Factored_);};
251 
253  bool A_Equilibrated() {return(A_Equilibrated_);};
254 
256  bool B_Equilibrated() {return(B_Equilibrated_);};
257 
259  virtual bool ShouldEquilibrate() {ComputeEquilibrateScaling(); return(ShouldEquilibrate_);};
260 
262  bool SolutionErrorsEstimated() {return(SolutionErrorsEstimated_);};
263 
265  bool Inverted() {return(Inverted_);};
266 
268  bool ReciprocalConditionEstimated() {return(ReciprocalConditionEstimated_);};
269 
271  bool Solved() {return(Solved_);};
272 
274  bool SolutionRefined() {return(SolutionRefined_);};
276 
278 
279 
281  Epetra_SerialDenseMatrix * Matrix() const {return(Matrix_);};
282 
284  Epetra_SerialDenseMatrix * FactoredMatrix() const {return(Factor_);};
285 
287  Epetra_SerialDenseMatrix * LHS() const {return(LHS_);};
288 
290  Epetra_SerialDenseMatrix * RHS() const {return(RHS_);};
291 
293  int M() const {return(M_);};
294 
296  int N() const {return(N_);};
297 
299  double * A() const {return(A_);};
300 
302  int LDA() const {return(LDA_);};
303 
305  double * B() const {return(B_);};
306 
308  int LDB() const {return(LDB_);};
309 
311  int NRHS() const {return(NRHS_);};
312 
314  double * X() const {return(X_);};
315 
317  int LDX() const {return(LDX_);};
318 
320  double * AF() const {return(AF_);};
321 
323  int LDAF() const {return(LDAF_);};
324 
326  int * IPIV() const {return(IPIV_);};
327 
329  double ANORM() const {return(ANORM_);};
330 
332  double RCOND() const {return(RCOND_);};
333 
335 
337  double ROWCND() const {return(ROWCND_);};
338 
340 
342  double COLCND() const {return(COLCND_);};
343 
345  double AMAX() const {return(AMAX_);};
346 
348  double * FERR() const {return(FERR_);};
349 
351  double * BERR() const {return(BERR_);};
352 
354  double * R() const {return(R_);};
355 
357  double * C() const {return(C_);};
359 
361 
362  virtual void Print(std::ostream& os) const;
365  protected:
366 
367  void AllocateWORK() {if (WORK_==0) {LWORK_ = 4*N_; WORK_ = new double[LWORK_];} return;};
368  void AllocateIWORK() {if (IWORK_==0) IWORK_ = new int[N_]; return;};
369  void InitPointers();
370  void DeleteArrays();
371  void ResetMatrix();
372  void ResetVectors();
373 
374 
380  bool Factored_;
383  bool Solved_;
384  bool Inverted_;
388 
389  char TRANS_;
390 
391  int M_;
392  int N_;
393  int Min_MN_;
394  int NRHS_;
395  int LDA_;
396  int LDAF_;
397  int LDB_;
398  int LDX_;
399  int INFO_;
400  int LWORK_;
401 
402  int * IPIV_;
403  int * IWORK_;
404 
405  double ANORM_;
406  double RCOND_;
407  double ROWCND_;
408  double COLCND_;
409  double AMAX_;
410 
415 
416  double * A_;
417  double * FERR_;
418  double * BERR_;
419  double * AF_;
420  double * WORK_;
421  double * R_;
422  double * C_;
423 
424  double * B_;
425  double * X_;
426 
427 
428  private:
429  // Epetra_SerialDenseSolver copy constructor (put here because we don't want user access)
430 
433 };
434 
435 #endif /* EPETRA_SERIALDENSESOLVER_H */
bool Transpose()
Returns true if transpose of this matrix has and will be used.
double * B() const
Returns pointer to current RHS.
virtual void Print(std::ostream &os) const
Print object to an output stream Print method.
double ROWCND() const
Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed)...
bool B_Equilibrated()
Returns true if RHS is equilibrated (RHS available via B() and LDB()).
bool SolutionRefined()
Returns true if the current set of vectors has been refined.
double COLCND() const
Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed...
double * C() const
Returns a pointer to the column scale vector used for equilibration.
double * BERR() const
Returns a pointer to the backward error estimates computed by LAPACK.
int M() const
Returns row dimension of system.
void SolveWithTranspose(bool Flag)
If Flag is true, causes all subsequent function calls to work with the transpose of this matrix...
double * R() const
Returns a pointer to the row scaling vector used for equilibration.
bool SolutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
Epetra_CompObject & operator=(const Epetra_CompObject &src)
Epetra_SerialDenseMatrix: A class for constructing and using real double precision general dense matr...
Epetra_SerialDenseMatrix * LHS() const
Returns pointer to current LHS.
double ANORM() const
Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
int * IPIV() const
Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
Epetra_SerialDenseMatrix * RHS_
double RCOND() const
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed)...
Epetra_BLAS: The Epetra BLAS Wrapper Class.
Definition: Epetra_BLAS.h:70
int N() const
Returns column dimension of system.
Epetra_SerialDenseMatrix * RHS() const
Returns pointer to current RHS.
Epetra_Object: The base Epetra class.
Definition: Epetra_Object.h:57
virtual bool ShouldEquilibrate()
Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system...
Epetra_SerialDenseMatrix * Matrix_
Epetra_CompObject: Functionality and data that is common to all computational classes.
int LDB() const
Returns the leading dimension of the RHS.
bool Inverted()
Returns true if matrix inverse has been computed (inverse available via AF() and LDAF()).
double * X() const
Returns pointer to current solution.
int LDA() const
Returns the leading dimension of the this matrix.
bool ReciprocalConditionEstimated()
Returns true if the condition number of the this matrix has been computed (value available via Recipr...
Epetra_LAPACK: The Epetra LAPACK Wrapper Class.
Definition: Epetra_LAPACK.h:67
int LDAF() const
Returns the leading dimension of the factored matrix.
void SolveToRefinedSolution(bool Flag)
Causes all solves to compute solution to best ability using iterative refinement. ...
int LDX() const
Returns the leading dimension of the solution.
Epetra_SerialDenseMatrix * FactoredMatrix() const
Returns pointer to factored matrix (assuming factorization has been performed).
bool A_Equilibrated()
Returns true if factor is equilibrated (factor available via AF() and LDAF()).
Epetra_SerialDenseMatrix * Factor_
void FactorWithEquilibration(bool Flag)
Causes equilibration to be called just before the matrix factorization as part of the call to Factor...
bool Factored()
Returns true if matrix is factored (factor available via AF() and LDAF()).
Epetra_SerialDenseMatrix * Matrix() const
Returns pointer to current matrix.
double AMAX() const
Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed)...
double * FERR() const
Returns a pointer to the forward error estimates computed by LAPACK.
bool Solved()
Returns true if the current set of vectors has been solved.
int Solve(int, TYPE *, TYPE *, TYPE *)
double * AF() const
Returns pointer to the factored matrix (may be the same as A() if factorization done in place)...
int NRHS() const
Returns the number of current right hand sides and solution vectors.
Epetra_SerialDenseMatrix * LHS_
Epetra_SerialDenseSolver: A class for solving dense linear problems.
double * A() const
Returns pointer to the this matrix.