Epetra Package Browser (Single Doxygen Collection)  Development
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 */
Epetra_SerialDenseSolver::LWORK_
int LWORK_
Definition: Epetra_SerialDenseSolver.h:400
Epetra_SerialDenseSolver::SolutionErrorsEstimated_
bool SolutionErrorsEstimated_
Definition: Epetra_SerialDenseSolver.h:382
Epetra_Object
Epetra_Object: The base Epetra class.
Definition: Epetra_Object.h:57
Epetra_SerialDenseSolver::LDB_
int LDB_
Definition: Epetra_SerialDenseSolver.h:397
Epetra_SerialDenseSolver::A_Equilibrated_
bool A_Equilibrated_
Definition: Epetra_SerialDenseSolver.h:377
Epetra_SerialDenseSolver::Min_MN_
int Min_MN_
Definition: Epetra_SerialDenseSolver.h:393
Epetra_SerialDenseSolver::A_Equilibrated
bool A_Equilibrated()
Returns true if factor is equilibrated (factor available via AF() and LDAF()).
Definition: Epetra_SerialDenseSolver.h:253
Epetra_SerialDenseSolver::Inverted_
bool Inverted_
Definition: Epetra_SerialDenseSolver.h:384
Epetra_SerialDenseSolver::WORK_
double * WORK_
Definition: Epetra_SerialDenseSolver.h:420
Epetra_SerialDenseSolver::B
double * B() const
Returns pointer to current RHS.
Definition: Epetra_SerialDenseSolver.h:305
Epetra_SerialDenseSolver::C_
double * C_
Definition: Epetra_SerialDenseSolver.h:422
Epetra_SerialDenseSolver::COLCND
double COLCND() const
Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed...
Definition: Epetra_SerialDenseSolver.h:342
Epetra_SerialDenseSolver::Matrix_
Epetra_SerialDenseMatrix * Matrix_
Definition: Epetra_SerialDenseSolver.h:411
Epetra_SerialDenseSolver::RHS_
Epetra_SerialDenseMatrix * RHS_
Definition: Epetra_SerialDenseSolver.h:413
Epetra_LAPACK
Epetra_LAPACK: The Epetra LAPACK Wrapper Class.
Definition: Epetra_LAPACK.h:67
Epetra_SerialDenseSolver::Matrix
Epetra_SerialDenseMatrix * Matrix() const
Returns pointer to current matrix.
Definition: Epetra_SerialDenseSolver.h:281
Epetra_SerialDenseSolver::ReciprocalConditionEstimated_
bool ReciprocalConditionEstimated_
Definition: Epetra_SerialDenseSolver.h:385
Epetra_SerialDenseSolver::ShouldEquilibrate
virtual bool ShouldEquilibrate()
Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system.
Definition: Epetra_SerialDenseSolver.h:259
Epetra_SerialDenseSolver::A
double * A() const
Returns pointer to the this matrix.
Definition: Epetra_SerialDenseSolver.h:299
Epetra_SerialDenseSolver::COLCND_
double COLCND_
Definition: Epetra_SerialDenseSolver.h:408
Epetra_SerialDenseSolver::AMAX
double AMAX() const
Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed).
Definition: Epetra_SerialDenseSolver.h:345
Epetra_SerialDenseSolver::SolveToRefinedSolution
void SolveToRefinedSolution(bool Flag)
Causes all solves to compute solution to best ability using iterative refinement.
Definition: Epetra_SerialDenseSolver.h:172
Epetra_CompObject.h
Epetra_SerialDenseSolver::B_
double * B_
Definition: Epetra_SerialDenseSolver.h:424
Epetra_SerialDenseSolver::Factor_
Epetra_SerialDenseMatrix * Factor_
Definition: Epetra_SerialDenseSolver.h:414
Epetra_SerialDenseSolver::AF_
double * AF_
Definition: Epetra_SerialDenseSolver.h:419
Epetra_BLAS.h
Epetra_SerialDenseSolver::RCOND_
double RCOND_
Definition: Epetra_SerialDenseSolver.h:406
Epetra_SerialDenseSolver::RefineSolution_
bool RefineSolution_
Definition: Epetra_SerialDenseSolver.h:386
Epetra_SerialDenseSolver::LDX
int LDX() const
Returns the leading dimension of the solution.
Definition: Epetra_SerialDenseSolver.h:317
Epetra_SerialDenseSolver::B_Equilibrated
bool B_Equilibrated()
Returns true if RHS is equilibrated (RHS available via B() and LDB()).
Definition: Epetra_SerialDenseSolver.h:256
Epetra_SerialDenseSolver::LDB
int LDB() const
Returns the leading dimension of the RHS.
Definition: Epetra_SerialDenseSolver.h:308
Epetra_SerialDenseSolver::SolutionRefined_
bool SolutionRefined_
Definition: Epetra_SerialDenseSolver.h:387
Epetra_SerialDenseSolver::AMAX_
double AMAX_
Definition: Epetra_SerialDenseSolver.h:409
Epetra_SerialDenseSolver::FactoredMatrix
Epetra_SerialDenseMatrix * FactoredMatrix() const
Returns pointer to factored matrix (assuming factorization has been performed).
Definition: Epetra_SerialDenseSolver.h:284
Epetra_SerialDenseSolver::AllocateWORK
void AllocateWORK()
Definition: Epetra_SerialDenseSolver.h:367
Epetra_SerialDenseSolver::LDAF_
int LDAF_
Definition: Epetra_SerialDenseSolver.h:396
Epetra_SerialDenseSolver::ROWCND_
double ROWCND_
Definition: Epetra_SerialDenseSolver.h:407
Epetra_SerialDenseSolver::N_
int N_
Definition: Epetra_SerialDenseSolver.h:392
Epetra_SerialDenseSolver::M
int M() const
Returns row dimension of system.
Definition: Epetra_SerialDenseSolver.h:293
Epetra_SerialDenseSolver::TRANS_
char TRANS_
Definition: Epetra_SerialDenseSolver.h:389
Epetra_SerialDenseSolver::NRHS_
int NRHS_
Definition: Epetra_SerialDenseSolver.h:394
Epetra_CompObject::operator=
Epetra_CompObject & operator=(const Epetra_CompObject &src)
Definition: Epetra_CompObject.h:114
Epetra_SerialDenseSolver::FERR
double * FERR() const
Returns a pointer to the forward error estimates computed by LAPACK.
Definition: Epetra_SerialDenseSolver.h:348
Epetra_SerialDenseSolver::X
double * X() const
Returns pointer to current solution.
Definition: Epetra_SerialDenseSolver.h:314
Epetra_LAPACK.h
Epetra_SerialDenseSolver::LDX_
int LDX_
Definition: Epetra_SerialDenseSolver.h:398
Epetra_Object::Print
virtual void Print(std::ostream &os) const
Print object to an output stream Print method.
Definition: Epetra_Object.cpp:97
Epetra_SerialDenseSolver::LHS_
Epetra_SerialDenseMatrix * LHS_
Definition: Epetra_SerialDenseSolver.h:412
Solve
int Solve(int, TYPE *, TYPE *, TYPE *)
Epetra_SerialDenseSolver::AF
double * AF() const
Returns pointer to the factored matrix (may be the same as A() if factorization done in place).
Definition: Epetra_SerialDenseSolver.h:320
Epetra_SerialDenseSolver::ROWCND
double ROWCND() const
Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed).
Definition: Epetra_SerialDenseSolver.h:337
Epetra_SerialDenseSolver::AllocateIWORK
void AllocateIWORK()
Definition: Epetra_SerialDenseSolver.h:368
Epetra_SerialDenseSolver::NRHS
int NRHS() const
Returns the number of current right hand sides and solution vectors.
Definition: Epetra_SerialDenseSolver.h:311
Epetra_SerialDenseSolver::SolveWithTranspose
void SolveWithTranspose(bool Flag)
If Flag is true, causes all subsequent function calls to work with the transpose of this matrix,...
Definition: Epetra_SerialDenseSolver.h:169
Epetra_SerialDenseSolver::RCOND
double RCOND() const
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed).
Definition: Epetra_SerialDenseSolver.h:332
Epetra_SerialDenseSolver::EstimateSolutionErrors_
bool EstimateSolutionErrors_
Definition: Epetra_SerialDenseSolver.h:381
Epetra_SerialDenseSolver::Transpose_
bool Transpose_
Definition: Epetra_SerialDenseSolver.h:379
Epetra_SerialDenseSolver::R
double * R() const
Returns a pointer to the row scaling vector used for equilibration.
Definition: Epetra_SerialDenseSolver.h:354
Epetra_BLAS
Epetra_BLAS: The Epetra BLAS Wrapper Class.
Definition: Epetra_BLAS.h:70
Epetra_SerialDenseSolver::R_
double * R_
Definition: Epetra_SerialDenseSolver.h:421
Epetra_SerialDenseSolver::Factored
bool Factored()
Returns true if matrix is factored (factor available via AF() and LDAF()).
Definition: Epetra_SerialDenseSolver.h:250
Epetra_SerialDenseSolver::ShouldEquilibrate_
bool ShouldEquilibrate_
Definition: Epetra_SerialDenseSolver.h:376
Epetra_SerialDenseSolver::Factored_
bool Factored_
Definition: Epetra_SerialDenseSolver.h:380
Epetra_SerialDenseSolver::IPIV_
int * IPIV_
Definition: Epetra_SerialDenseSolver.h:402
Epetra_SerialDenseSolver::BERR
double * BERR() const
Returns a pointer to the backward error estimates computed by LAPACK.
Definition: Epetra_SerialDenseSolver.h:351
Epetra_SerialDenseSolver::BERR_
double * BERR_
Definition: Epetra_SerialDenseSolver.h:418
Epetra_SerialDenseSolver::ANORM_
double ANORM_
Definition: Epetra_SerialDenseSolver.h:405
Epetra_SerialDenseSolver::IWORK_
int * IWORK_
Definition: Epetra_SerialDenseSolver.h:403
Epetra_SerialDenseSolver::FERR_
double * FERR_
Definition: Epetra_SerialDenseSolver.h:417
Epetra_SerialDenseSolver::C
double * C() const
Returns a pointer to the column scale vector used for equilibration.
Definition: Epetra_SerialDenseSolver.h:357
Epetra_SerialDenseSolver::FactorWithEquilibration
void FactorWithEquilibration(bool Flag)
Causes equilibration to be called just before the matrix factorization as part of the call to Factor.
Definition: Epetra_SerialDenseSolver.h:166
Epetra_SerialDenseSolver::INFO_
int INFO_
Definition: Epetra_SerialDenseSolver.h:399
Epetra_SerialDenseSolver::Transpose
bool Transpose()
Returns true if transpose of this matrix has and will be used.
Definition: Epetra_SerialDenseSolver.h:247
A
Epetra_SerialDenseSolver::A_
double * A_
Definition: Epetra_SerialDenseSolver.h:416
Epetra_SerialDenseSolver::ReciprocalConditionEstimated
bool ReciprocalConditionEstimated()
Returns true if the condition number of the this matrix has been computed (value available via Recipr...
Definition: Epetra_SerialDenseSolver.h:268
Epetra_SerialDenseSolver::N
int N() const
Returns column dimension of system.
Definition: Epetra_SerialDenseSolver.h:296
Epetra_SerialDenseSolver::LDA
int LDA() const
Returns the leading dimension of the this matrix.
Definition: Epetra_SerialDenseSolver.h:302
Epetra_CompObject
Epetra_CompObject: Functionality and data that is common to all computational classes.
Definition: Epetra_CompObject.h:57
Epetra_SerialDenseSolver::B_Equilibrated_
bool B_Equilibrated_
Definition: Epetra_SerialDenseSolver.h:378
Epetra_SerialDenseSolver::SolutionRefined
bool SolutionRefined()
Returns true if the current set of vectors has been refined.
Definition: Epetra_SerialDenseSolver.h:274
Epetra_SerialDenseSolver::RHS
Epetra_SerialDenseMatrix * RHS() const
Returns pointer to current RHS.
Definition: Epetra_SerialDenseSolver.h:290
Epetra_SerialDenseSolver::IPIV
int * IPIV() const
Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
Definition: Epetra_SerialDenseSolver.h:326
Epetra_SerialDenseSolver::Inverted
bool Inverted()
Returns true if matrix inverse has been computed (inverse available via AF() and LDAF()).
Definition: Epetra_SerialDenseSolver.h:265
Epetra_SerialDenseSolver::SolutionErrorsEstimated
bool SolutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
Definition: Epetra_SerialDenseSolver.h:262
Epetra_SerialDenseSolver::Solved
bool Solved()
Returns true if the current set of vectors has been solved.
Definition: Epetra_SerialDenseSolver.h:271
Epetra_SerialDenseSolver::LHS
Epetra_SerialDenseMatrix * LHS() const
Returns pointer to current LHS.
Definition: Epetra_SerialDenseSolver.h:287
Epetra_Object.h
Epetra_SerialDenseSolver::LDAF
int LDAF() const
Returns the leading dimension of the factored matrix.
Definition: Epetra_SerialDenseSolver.h:323
Epetra_SerialDenseMatrix
Epetra_SerialDenseMatrix: A class for constructing and using real double precision general dense matr...
Definition: Epetra_SerialDenseMatrix.h:107
Epetra_SerialDenseSolver::ANORM
double ANORM() const
Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
Definition: Epetra_SerialDenseSolver.h:329
Epetra_SerialDenseSolver::Equilibrate_
bool Equilibrate_
Definition: Epetra_SerialDenseSolver.h:375
Epetra_SerialDenseSolver::M_
int M_
Definition: Epetra_SerialDenseSolver.h:391
Epetra_SerialDenseSolver::Solved_
bool Solved_
Definition: Epetra_SerialDenseSolver.h:383
Epetra_SerialDenseSolver::LDA_
int LDA_
Definition: Epetra_SerialDenseSolver.h:395
B
Epetra_SerialDenseSolver
Epetra_SerialDenseSolver: A class for solving dense linear problems.
Definition: Epetra_SerialDenseSolver.h:132
Epetra_SerialDenseSolver::X_
double * X_
Definition: Epetra_SerialDenseSolver.h:425