IFPACK Development
Loading...
Searching...
No Matches
Ifpack_SerialTriDiSolver.h
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 IFPACK_SERIALTRIDISOLVER_H
45#define IFPACK_SERIALTRIDISOLVER_H
46
47#if defined(Ifpack_SHOW_DEPRECATED_WARNINGS)
48#ifdef __GNUC__
49#warning "The Ifpack package is deprecated"
50#endif
51#endif
53//class Epetra_SerialDenseVector;
55
56#include "Epetra_Object.h"
57#include "Epetra_CompObject.h"
58#include "Epetra_BLAS.h"
59#include "Teuchos_LAPACK.hpp"
60
61
63
135//=========================================================================
137 public Epetra_CompObject, public Epetra_BLAS,
138 public Epetra_Object {
139 public:
140
142
143
145
146
150
152
153
156
158
163
164
166
167
169
171 // void FactorWithEquilibration(bool Flag) {Equilibrate_ = Flag; return;};
172
174 void SolveWithTranspose(bool Flag) {Transpose_ = Flag; if (Flag) TRANS_ = 'T'; else TRANS_ = 'N'; return;};
175
177 void SolveToRefinedSolution(bool Flag) {RefineSolution_ = Flag; return;};
178
180
183 void EstimateSolutionErrors(bool Flag) ;
185
187
188
190
193 virtual int Factor(void);
194
196
199 virtual int Solve(void);
200
202
205 virtual int Invert(void);
206
208
211 virtual int ApplyRefinement(void);
212
214
217 // int UnequilibrateLHS(void);
218
220
226 virtual int ReciprocalConditionEstimate(double & Value);
228
230
231
233 bool Transpose() {return(Transpose_);};
234
236 bool Factored() {return(Factored_);};
237
239 bool SolutionErrorsEstimated() {return(SolutionErrorsEstimated_);};
240
242 bool Inverted() {return(Inverted_);};
243
245 bool ReciprocalConditionEstimated() {return(ReciprocalConditionEstimated_);};
246
248 bool Solved() {return(Solved_);};
249
251 bool SolutionRefined() {return(SolutionRefined_);};
253
255
256
258 Ifpack_SerialTriDiMatrix * Matrix() const {return(Matrix_);};
259
261 Ifpack_SerialTriDiMatrix * FactoredMatrix() const {return(Factor_);};
262
264 Epetra_SerialDenseMatrix * LHS() const {return(LHS_);};
265
267 Epetra_SerialDenseMatrix * RHS() const {return(RHS_);};
268
270 int N() const {return(N_);};
271
273 double * A() const {return(A_);};
274
276 int LDA() const {return(LDA_);};
277
279 double * B() const {return(B_);};
280
282 int LDB() const {return(LDB_);};
283
285 int NRHS() const {return(NRHS_);};
286
288 double * X() const {return(X_);};
289
291 int LDX() const {return(LDX_);};
292
294 double * AF() const {return(AF_);};
295
297 int LDAF() const {return(LDAF_);};
298
300 int * IPIV() const {return(IPIV_);};
301
303 double ANORM() const {return(ANORM_);};
304
306 double RCOND() const {return(RCOND_);};
307
309
311 double ROWCND() const {return(ROWCND_);};
312
314
316 double COLCND() const {return(COLCND_);};
317
319 double AMAX() const {return(AMAX_);};
320
322 double * FERR() const {return(FERR_);};
323
325 double * BERR() const {return(BERR_);};
326
328
330
331
332 virtual void Print(std::ostream& os) const;
334 protected:
335
336 void AllocateWORK() {if (WORK_==0) {LWORK_ = 4*N_; WORK_ = new double[LWORK_];} return;};
337 void AllocateIWORK() {if (IWORK_==0) IWORK_ = new int[N_]; return;};
338 void InitPointers();
339 void DeleteArrays();
340 void ResetMatrix();
341 void ResetVectors();
342
343 bool Transpose_;
344 bool Factored_;
345 bool EstimateSolutionErrors_;
346 bool SolutionErrorsEstimated_;
347 bool Solved_;
348 bool Inverted_;
349 bool ReciprocalConditionEstimated_;
350 bool RefineSolution_;
351 bool SolutionRefined_;
352
353 char TRANS_;
354
355 int N_;
356 int Min_MN_;
357 int NRHS_;
358 int LDA_;
359 int LDAF_;
360 int LDB_;
361 int LDX_;
362 int INFO_;
363 int LWORK_;
364
365 int * IPIV_;
366 int * IWORK_;
367
368 double ANORM_;
369 double RCOND_;
370 double ROWCND_;
371 double COLCND_;
372 double AMAX_;
373
374 Ifpack_SerialTriDiMatrix * Matrix_;
377 Ifpack_SerialTriDiMatrix * Factor_;
378
379 double * A_;
380 double * FERR_;
381 double * BERR_;
382 double * AF_;
383 double * WORK_;
384
385 double * B_;
386 double * X_;
387
388
389 private:
390 Teuchos::LAPACK<int,double> lapack;
391 // Ifpack_SerialTriDiSolver copy constructor (put here because we don't want user access)
392
394 Ifpack_SerialTriDiSolver & operator=(const Ifpack_SerialTriDiSolver& Source);
395};
396
397#endif /* EPETRA_SERIALTRIDISOLVER_H */
Ifpack_SerialTriDiMatrix: A class for constructing and using real double precision general TriDi matr...
Ifpack_SerialTriDiSolver: A class for solving TriDi linear problems.
bool ReciprocalConditionEstimated()
Returns true if the condition number of the this matrix has been computed (value available via Recipr...
double * BERR() const
Returns a pointer to the backward error estimates computed by LAPACK.
bool Inverted()
Returns true if matrix inverse has been computed (inverse available via AF() and LDAF()).
bool Solved()
Returns true if the current set of vectors has been solved.
void EstimateSolutionErrors(bool Flag)
Causes all solves to estimate the forward and backward solution error.
int N() const
Returns column dimension of system.
int LDX() const
Returns the leading dimension of the solution.
virtual int ApplyRefinement(void)
Apply Iterative Refinement.
int NRHS() const
Returns the number of current right hand sides and solution vectors.
double * X() const
Returns pointer to current solution.
double ROWCND() const
Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed).
virtual void Print(std::ostream &os) const
Print service methods; defines behavior of ostream << operator.
double * A() const
Returns pointer to the this matrix.
bool Factored()
Returns true if matrix is factored (factor available via AF() and LDAF()).
virtual int Invert(void)
Inverts the this matrix.
int SetMatrix(Ifpack_SerialTriDiMatrix &A)
Sets the pointers for coefficient matrix.
double AMAX() const
Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed).
virtual int Solve(void)
Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()....
bool Transpose()
Returns true if transpose of this matrix has and will be used.
Epetra_SerialDenseMatrix * LHS() const
Returns pointer to current LHS.
Ifpack_SerialTriDiMatrix * FactoredMatrix() const
Returns pointer to factored matrix (assuming factorization has been performed).
virtual int Factor(void)
Computes the in-place LU factorization of the matrix using the LAPACK routine DGETRF.
int LDAF() const
Returns the leading dimension of the factored matrix.
double * FERR() const
Returns a pointer to the forward error estimates computed by LAPACK.
virtual int ReciprocalConditionEstimate(double &Value)
Unscales the solution vectors if equilibration was used to solve the system.
int * IPIV() const
Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
int LDB() const
Returns the leading dimension of the RHS.
Epetra_SerialDenseMatrix * RHS() const
Returns pointer to current RHS.
double COLCND() const
Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed...
int SetVectors(Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &B)
Sets the pointers for left and right hand side vector(s).
Ifpack_SerialTriDiSolver()
Default constructor; matrix should be set using SetMatrix(), LHS and RHS set with SetVectors().
virtual ~Ifpack_SerialTriDiSolver()
Ifpack_SerialTriDiSolver destructor.
void SolveToRefinedSolution(bool Flag)
Causes all solves to compute solution to best ability using iterative refinement.
double * AF() const
Returns pointer to the factored matrix (may be the same as A() if factorization done in place).
void SolveWithTranspose(bool Flag)
Causes equilibration to be called just before the matrix factorization as part of the call to Factor.
double RCOND() const
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed).
bool SolutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
double * B() const
Returns pointer to current RHS.
int LDA() const
Returns the leading dimension of the this matrix.
bool SolutionRefined()
Returns true if the current set of vectors has been refined.
double ANORM() const
Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
Ifpack_SerialTriDiMatrix * Matrix() const
Returns pointer to current matrix.