53 #include "../epetra_test_err.h"
59 int N1,
int NRHS1,
double OneNorm1,
60 double *
B1,
int LDB1,
61 double * X1,
int LDX1,
62 bool Transpose,
bool verbose);
66 bool Residual(
int N,
int NRHS,
double *
A,
int LDA,
bool Transpose,
67 double * X,
int LDX,
double *
B,
int LDB,
double * resid);
79 int main(
int argc,
char *argv[])
81 int ierr = 0, i, j, k;
85 MPI_Init(&argc,&argv);
94 if (argc>1)
if (argv[1][0]==
'-' && argv[1][1]==
'v') verbose =
true;
96 if (verbose && Comm.
MyPID()==0)
99 int rank = Comm.
MyPID();
106 if (verbose) cout << Comm <<endl;
111 if (verbose && rank!=0) verbose =
false;
115 double *
A =
new double[N*N];
116 double * A1 =
new double[N*N];
117 double * X =
new double[(N+1)*NRHS];
118 double * X1 =
new double[(N+1)*NRHS];
121 double *
B =
new double[N*NRHS];
122 double *
B1 =
new double[N*NRHS];
129 bool Transpose =
false;
133 for (
int kk=0; kk<2; kk++) {
134 for (i=1; i<=N; i++) {
137 for (j=1; j<=i; j++) OneNorm1 += 1.0/((
double) j);
157 for (k=0; k<NRHS; k++)
158 for (j=0; j<i; j++) {
159 B[j+k*LDB] = 1.0/((double) (k+3)*(j+3));
160 B1[j+k*LDB1] = B[j+k*LDB1];
168 ierr =
check(solver, A1, LDA1, i, NRHS, OneNorm1, B1, LDB1, X1, LDX1, Transpose, verbose);
172 if (verbose) cout <<
"Factorization failed due to bad conditioning. This is normal if RCOND is small."
191 double ScalarA = 2.0;
197 for (i=0; i<DM; i++) D[j][i] = (
double) (1+i+j*DM) ;
201 double NormInfD_ref = (double)(DM*(DN*(DN+1))/2);
202 double NormOneD_ref = (double)((DM*DN*(DM*DN+1))/2 - (DM*(DN-1)*(DM*(DN-1)+1))/2 );
208 cout <<
" *** Before scaling *** " << endl
209 <<
" Computed one-norm of test matrix = " << NormOneD << endl
210 <<
" Expected one-norm = " << NormOneD_ref << endl
211 <<
" Computed inf-norm of test matrix = " << NormInfD << endl
212 <<
" Expected inf-norm = " << NormInfD_ref << endl;
218 cout <<
" *** After scaling *** " << endl
219 <<
" Computed one-norm of test matrix = " << NormOneD << endl
220 <<
" Expected one-norm = " << NormOneD_ref*ScalarA << endl
221 <<
" Computed inf-norm of test matrix = " << NormInfD << endl
222 <<
" Expected inf-norm = " << NormInfD_ref*ScalarA << endl;
241 smallA(i,j) = 1.0*i+2.0*j+1.0;
250 if (verbose) cout <<
"err in Epetra_SerialDenseMatrix::operator==, "
251 <<
"erroneously returned true." << std::endl;
257 if (verbose) cout <<
"err in Epetra_SerialDenseMatrix::operator==, "
258 <<
"erroneously returned true." << std::endl;
262 int err1 = smallA.
Multiply(
false, x, y1);
263 int err2 = y2.
Multiply(
'N',
'N', 1.0, smallA, x, 0.0);
264 if (err1 != 0 || err2 != 0) {
265 if (verbose) cout <<
"err in Epetra_SerialDenseMatrix::Multiply"<<endl;
270 if (y1(i,0) != y2(i,0)) {
271 if (verbose) cout <<
"different versions of Multiply don't match."<<endl;
287 if (verbose) cout <<
"\n\nComputing factor of an " << N <<
" x " << N <<
" general matrix...Please wait.\n\n" << endl;
291 A =
new double[LDA*N];
292 X =
new double[LDB*NRHS];
294 for (j=0; j<N; j++) {
295 for (k=0; k<NRHS; k++) X[j+k*LDX] = 1.0/((
double) (j+5+k));
296 for (i=0; i<N; i++) {
297 if (i==((j+2)%N)) A[i+j*LDA] = 100.0 + i;
298 else A[i+j*LDA] = -11.0/((double) (i+5)*(j+2));
317 ierr = BigSolver.
Factor();
318 if (ierr!=0 && verbose) cout <<
"Error in factorization = "<<ierr<< endl;
322 double FLOPS = counter.
Flops();
323 double MFLOPS = FLOPS/time/1000000.0;
324 if (verbose) cout <<
"MFLOPS for Factorization = " << MFLOPS << endl;
336 RHS.
Multiply(
'N',
'N', 1.0, OrigBigMatrix, LHS, 0.0);
341 FLOPS = RHS_counter.
Flops();
342 MFLOPS = FLOPS/time/1000000.0;
343 if (verbose) cout <<
"MFLOPS to build RHS (NRHS = " << NRHS <<
") = " << MFLOPS << endl;
349 ierr = BigSolver.
Solve();
350 if (ierr==1 && verbose) cout <<
"LAPACK guidelines suggest this matrix might benefit from equilibration." << endl;
351 else if (ierr!=0 && verbose) cout <<
"Error in solve = "<<ierr<< endl;
355 FLOPS = BigSolver.
Flops();
356 MFLOPS = FLOPS/time/1000000.0;
357 if (verbose) cout <<
"MFLOPS for Solve (NRHS = " << NRHS <<
") = " << MFLOPS << endl;
359 double * resid =
new double[NRHS];
361 OrigRHS.
A(), OrigRHS.
LDA(), resid);
364 if (!OK) cout <<
"************* Residual do not meet tolerance *************" << endl;
365 for (i=0; i<NRHS; i++)
366 cout <<
"Residual[" << i <<
"] = "<< resid[i] << endl;
374 X2.
Size(BigMatrix.
N());
375 B2.
Size(BigMatrix.
M());
376 int length = BigMatrix.
N();
377 {
for (
int kk=0; kk<length; kk++) X2[kk] = ((
double ) kk)/ ((double) length);}
382 B2.
Multiply(
'N',
'N', 1.0, OrigBigMatrix, X2, 0.0);
387 FLOPS = RHS_counter.
Flops();
388 MFLOPS = FLOPS/time/1000000.0;
389 if (verbose) cout <<
"MFLOPS to build single RHS = " << MFLOPS << endl;
395 ierr = BigSolver.
Solve();
397 if (ierr==1 && verbose) cout <<
"LAPACK guidelines suggest this matrix might benefit from equilibration." << endl;
398 else if (ierr!=0 && verbose) cout <<
"Error in solve = "<<ierr<< endl;
401 FLOPS = counter.
Flops();
402 MFLOPS = FLOPS/time/1000000.0;
403 if (verbose) cout <<
"MFLOPS to solve single RHS = " << MFLOPS << endl;
406 OrigB2.
LDA(), resid);
409 if (!OK) cout <<
"************* Residual do not meet tolerance *************" << endl;
410 cout <<
"Residual = "<< resid[0] << endl;
423 double * C1 =
new double[N*N];
435 for (j=0; j<N; j++) {
436 assert(C(i,j) == C1[i+j*N]);
437 assert(C(i,j) == C[j][i]);
441 cout <<
"Default constructor and index operator check OK. Values of Hilbert matrix = "
443 <<
"Values should be 1/(i+j+1), except value (1,2) should be 1000" << endl;
449 assert(shapedMatrix.
M() == 10);
450 assert(shapedMatrix.
N() == 12);
451 for(i = 0; i < 10; i++)
452 for(j = 0; j < 12; j++)
453 assert(shapedMatrix(i, j) == 0.0);
455 assert(sizedVector.
Length() == 20);
456 for(i = 0; i < 20; i++)
457 assert(sizedVector(i) == 0.0);
459 cout <<
"Shaped/sized constructors check OK." << endl;
464 if(verbose && temperr == 0)
465 cout <<
"Operator = checked OK." << endl;
468 if(verbose && temperr == 0)
469 cout <<
"Copy ctr checked OK." << endl;
485 if (v1.
Norm1()!=6.0) temperr++;
486 if (fabs(sqrt(14.0)-v1.
Norm2())>1.0e-6) temperr++;
487 if (v1.
NormInf()!=3.0) temperr++;
488 if(verbose && temperr == 0)
489 cout <<
"Vector Norms checked OK." << endl;
491 if (v1.
Dot(v2)!=1.0) temperr++;
492 if(verbose && temperr == 0)
493 cout <<
"Vector Dot product checked OK." << endl;
505 int N1,
int NRHS1,
double OneNorm1,
506 double *
B1,
int LDB1,
507 double * X1,
int LDX1,
508 bool Transpose,
bool verbose) {
515 if (verbose) cout <<
"\n\nNumber of Rows = " << M << endl<< endl;
519 if (verbose) cout <<
"\n\nNumber of Equations = " << N << endl<< endl;
522 int LDA = solver.
LDA();
523 if (verbose) cout <<
"\n\nLDA = " << LDA << endl<< endl;
526 int LDB = solver.
LDB();
527 if (verbose) cout <<
"\n\nLDB = " << LDB << endl<< endl;
530 int LDX = solver.
LDX();
531 if (verbose) cout <<
"\n\nLDX = " << LDX << endl<< endl;
534 int NRHS = solver.
NRHS();
535 if (verbose) cout <<
"\n\nNRHS = " << NRHS << endl<< endl;
538 assert(solver.
ANORM()==-1.0);
539 assert(solver.
RCOND()==-1.0);
541 assert(solver.
ROWCND()==-1.0);
542 assert(solver.
COLCND()==-1.0);
543 assert(solver.
AMAX()==-1.0);
559 int ierr = solver.
Factor();
561 if (ierr!=0)
return(ierr);
567 double rcond1 = 1.0/exp(3.5*((
double)N));
568 if (N==1) rcond1 = 1.0;
569 cout <<
"\n\nRCOND = "<< rcond <<
" should be approx = "
570 << rcond1 << endl << endl;
573 ierr = solver.
Solve();
575 if (ierr!=0 && verbose) cout <<
"LAPACK rules suggest system should be equilibrated." << endl;
584 cout <<
"\n\nFERR[0] = "<< solver.
FERR()[0] << endl;
585 cout <<
"\n\nBERR[0] = "<< solver.
BERR()[0] << endl<< endl;
589 double * resid =
new double[NRHS];
592 if (!OK) cout <<
"************* Residual do not meet tolerance *************" << endl;
603 cout <<
"\n\nResiduals using factorization to solve" << endl;
604 for (i=0; i<NRHS; i++)
605 cout <<
"Residual[" << i <<
"] = "<< resid[i] << endl;
623 assert(solver.
Solve()>-1);
630 if (!OK) cout <<
"************* Residual do not meet tolerance *************" << endl;
631 cout <<
"Residuals using inverse to solve" << endl;
632 for (i=0; i<NRHS; i++)
633 cout <<
"Residual[" << i <<
"] = "<< resid[i] << endl;
643 for (
int j=0; j<N; j++)
644 for (
int i=0; i<N; i++)
645 A[i+j*LDA] = 1.0/((
double)(i+j+1));
649 bool Residual(
int N,
int NRHS,
double *
A,
int LDA,
bool Transpose,
650 double * X,
int LDX,
double *
B,
int LDB,
double * resid) {
654 if (Transpose) Transa =
'T';
655 Blas.
GEMM(Transa,
'N', N, NRHS, N, -1.0, A, LDA,
656 X, LDX, 1.0, B, LDB);
658 for (
int i=0; i<NRHS; i++) {
659 resid[i] = Blas.
NRM2(N, B+i*LDB);
660 if (resid[i]>1.0
E-7) OK =
false;
681 if(verbose) cout <<
"Checking copy->copy (new alloc)" << endl;
686 cout <<
"before assignment:" << endl;
692 cout <<
"after assignment:" << endl;
702 if(verbose) cout <<
"Checked OK." << endl;
708 if(verbose) cout <<
"\nChecking copy->copy (no alloc)" << endl;
713 double* origA = lhs.A();
714 int origLDA = lhs.LDA();
716 cout <<
"before assignment:" << endl;
722 cout <<
"after assignment:" << endl;
738 if(verbose) cout <<
"Checked OK." << endl;
744 if(verbose) cout <<
"\nChecking view->copy" << endl;
750 cout <<
"before assignment:" << endl;
756 cout <<
"after assignment:" << endl;
767 if(verbose) cout <<
"Checked OK." << endl;
773 if(verbose) cout <<
"\nChecking copy->view" << endl;
778 cout <<
"before assignment:" << endl;
784 cout <<
"after assignment:" << endl;
794 if(verbose) cout <<
"Checked OK." << endl;
800 if(verbose) cout <<
"\nChecking view->view" << endl;
806 cout <<
"before assignment:" << endl;
812 cout <<
"after assignment:" << endl;
823 if(verbose) cout <<
"Checked OK." << endl;
832 const int m1rows = 5;
833 const int m1cols = 4;
834 const int m2rows = 2;
835 const int m2cols = 6;
839 if(verbose)
printHeading(
"Testing matrix copy constructors");
841 if(verbose) cout <<
"checking copy constructor (view)" << endl;
843 if(debug)
printArray(m1rand, m1rows * m1cols);
846 cout <<
"original matrix:" << endl;
851 cout <<
"clone matrix:" << endl;
854 if(verbose) cout <<
"making sure signatures match" << endl;
859 if(verbose) cout <<
"Checked OK." << endl;
862 if(verbose) cout <<
"\nchecking copy constructor (copy)" << endl;
864 if(debug)
printArray(m2rand, m2rows * m2cols);
867 cout <<
"original matrix:" << endl;
872 cout <<
"clone matrix:" << endl;
875 if(verbose) cout <<
"checking that signatures match" << endl;
879 if(verbose) cout <<
"Checked OK." << endl;
882 if(verbose) cout <<
"\nmodifying entry in m2, m2clone should be unchanged" << endl;
886 cout <<
"orig:" << endl;
888 cout <<
"clone:" << endl;
894 if(verbose) cout <<
"Checked OK." << endl;
903 cout <<
"\n==================================================================\n";
904 cout << heading << endl;
905 cout <<
"==================================================================\n";
912 cout <<
"*** " << name <<
" ***" << endl;
920 cout <<
"user array (size " << length <<
"): ";
921 for(
int i = 0; i < length; i++)
922 cout << array[i] <<
" ";
930 const double a = 16807.0;
931 const double BigInt = 2147483647.0;
932 const double DbleOne = 1.0;
933 const double DbleTwo = 2.0;
934 double seed = rand();
936 double* array =
new double[length];
938 for(
int i = 0; i < length; i++) {
939 seed = fmod(a * seed, BigInt);
940 array[i] = DbleTwo * (seed / BigInt) - DbleOne;
950 if((a.
M() != b.
M() )||
966 for(
int i = 0; i < m; i++)
967 for(
int j = 0; j < n; j++) {
984 double orig_a = a(r,c);
985 double new_value = a(r,c) + 1;
986 if(b(r,c) == new_value)
990 if(b(r,c) == new_value)
void SetFlopCounter(const Epetra_Flops &FlopCounter_in)
Set the internal Epetra_Flops() pointer.
virtual double NormInf() const
Computes the Infinity-Norm of the this matrix.
bool Transpose()
Returns true if transpose of this matrix has and will be used.
void GEMM(const char TRANSA, const char TRANSB, const int M, const int N, const int K, const float ALPHA, const float *A, const int LDA, const float *B, const int LDB, const float BETA, float *C, const int LDC) const
Epetra_BLAS matrix-matrix multiply function (SGEMM)
int Size(int Length_in)
Set length of a Epetra_SerialDenseVector object; init values to zero.
double ROWCND() const
Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed)...
void printMat(const char *name, Epetra_IntSerialDenseMatrix &matrix)
int matrixAssignment(bool verbose, bool debug)
double Norm2() const
Compute 2-norm of each vector in multi-vector.
float NRM2(const int N, const float *X, const int INCX=1) const
Epetra_BLAS norm function (SNRM2).
#define EPETRA_TEST_ERR(a, b)
bool B_Equilibrated()
Returns true if RHS is equilibrated (RHS available via B() and LDB()).
double ElapsedTime(void) const
Epetra_Time elapsed time function.
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...
static void SetTracebackMode(int TracebackModeValue)
Set the value of the Epetra_Object error traceback report mode.
double * BERR() const
Returns a pointer to the backward error estimates computed by LAPACK.
Epetra_DataAccess CV() const
Returns the data access mode of the this matrix.
bool seperateData(Epetra_IntSerialDenseMatrix &a, Epetra_IntSerialDenseMatrix &b)
virtual double NormOne() const
Computes the 1-Norm of the this matrix.
double NormInf() const
Compute Inf-norm of each vector in multi-vector.
int SetVectors(Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &B)
Sets the pointers for left and right hand side vector(s).
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...
bool SolutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
int MyPID() const
Return my process ID.
Epetra_MpiComm: The Epetra MPI Communication Class.
Epetra_SerialDenseMatrix: A class for constructing and using real double precision general dense matr...
std::string Epetra_Version()
int check(Epetra_CrsGraph &A, int NumMyRows1, int NumGlobalRows1, int NumMyNonzeros1, int NumGlobalNonzeros1, int *MyGlobalElements, bool verbose)
double ANORM() const
Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
Epetra_Time: The Epetra Timing Class.
void printHeading(const char *heading)
double RCOND() const
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed)...
double * A() const
Returns pointer to the this matrix.
Epetra_BLAS: The Epetra BLAS Wrapper Class.
Epetra_SerialDenseVector: A class for constructing and using dense vectors.
bool Residual(int N, int NRHS, double *A, int LDA, bool Transpose, double *X, int LDX, double *B, int LDB, double *resid)
int N() const
Returns column dimension of system.
int matrixCpyCtr(bool verbose, bool debug)
int Length() const
Returns length of vector.
void printArray(int *array, int length)
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_SerialComm: The Epetra Serial Communication Class.
void ResetFlops()
Resets the number of floating point operations to zero for this multi-vector.
double Norm1() const
Compute 1-norm of each vector in multi-vector.
virtual int Solve(void)
Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()..
Epetra_Flops: The Epetra Floating Point Operations Class.
virtual int Invert(void)
Inverts the this matrix.
void SolveToRefinedSolution(bool Flag)
Causes all solves to compute solution to best ability using iterative refinement. ...
int Scale(double ScalarA)
Inplace scalar-matrix product A = a A.
int SetMatrix(Epetra_SerialDenseMatrix &A)
Sets the pointers for coefficient matrix.
int LDA() const
Returns the leading dimension of the this matrix.
int LDX() const
Returns the leading dimension of the solution.
double Flops() const
Returns the number of floating point operations with this object and resets the count.
int Multiply(char TransA, char TransB, double ScalarAB, const Epetra_SerialDenseMatrix &A, const Epetra_SerialDenseMatrix &B, double ScalarThis)
Matrix-Matrix multiplication, this = ScalarThis*this + ScalarAB*A*B.
bool A_Equilibrated()
Returns true if factor is equilibrated (factor available via AF() and LDAF()).
double Flops() const
Returns the number of floating point operations with this multi-vector.
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()).
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 main(int argc, char *argv[])
double Dot(const Epetra_SerialDenseVector &x) const
Compute 1-norm of each vector in multi-vector.
int N() const
Returns column dimension of system.
int NRHS() const
Returns the number of current right hand sides and solution vectors.
int * getRandArray(int length)
int Shape(int NumRows, int NumCols)
Set dimensions of a Epetra_SerialDenseMatrix object; init values to zero.
bool identicalSignatures(Epetra_IntSerialDenseMatrix &a, Epetra_IntSerialDenseMatrix &b, bool testLDA=true)
virtual int ReciprocalConditionEstimate(double &Value)
Returns the reciprocal of the 1-norm condition number of the this matrix.
virtual int Factor(void)
Computes the in-place LU factorization of the matrix using the LAPACK routine DGETRF.
int M() const
Returns row dimension of system.
void GenerateHilbert(double *A, int LDA, int N)
Epetra_SerialDenseSolver: A class for solving dense linear problems.