EpetraExt Package Browser (Single Doxygen Collection) Development
test/Transpose/cxx_main.cpp
Go to the documentation of this file.
00001 //@HEADER
00002 // ***********************************************************************
00003 //
00004 //     EpetraExt: Epetra Extended - Linear Algebra Services Package
00005 //                 Copyright (2011) Sandia Corporation
00006 //
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00038 //
00039 // ***********************************************************************
00040 //@HEADER
00041 
00042 // Epetra_BlockMap Test routine
00043 
00044 #include "EpetraExt_Version.h"
00045 
00046 #include "Epetra_CrsMatrix.h"
00047 #include "Epetra_VbrMatrix.h"
00048 #include "Epetra_Vector.h"
00049 #include "Epetra_MultiVector.h"
00050 #include "Epetra_LocalMap.h"
00051 #include "Epetra_IntVector.h"
00052 #include "Epetra_Map.h"
00053 
00054 #ifdef EPETRA_MPI
00055 #include "Epetra_MpiComm.h"
00056 #include <mpi.h>
00057 #else
00058 #include "Epetra_SerialComm.h"
00059 #endif
00060 
00061 #include "Epetra_Time.h"
00062 #include "EpetraExt_Transpose_RowMatrix.h"
00063 #include "Trilinos_Util.h"
00064 
00065 int checkResults(Epetra_RowMatrix * A, Epetra_CrsMatrix * transA, 
00066                  Epetra_Vector * xexact, bool verbose);
00067 
00068 int main(int argc, char *argv[])
00069 {
00070   int i;
00071 
00072 #ifdef EPETRA_MPI
00073   // Initialize MPI
00074   MPI_Init(&argc,&argv);
00075   Epetra_MpiComm comm(MPI_COMM_WORLD);
00076 #else
00077   Epetra_SerialComm comm;
00078 #endif
00079 
00080   // Uncomment to debug in parallel int tmp; if (comm.MyPID()==0) cin >> tmp; comm.Barrier();
00081 
00082   bool verbose = false;
00083 
00084   // Check if we should print results to standard out
00085   if (argc>1) if (argv[1][0]=='-' && argv[1][1]=='v') verbose = true;
00086 
00087   if (!verbose) comm.SetTracebackMode(0); // This should shut down any error traceback reporting
00088 
00089   if (verbose) cout << comm << endl << flush;
00090 
00091   if (verbose) verbose = (comm.MyPID()==0);
00092 
00093   if (verbose)
00094     cout << EpetraExt::EpetraExt_Version() << endl << endl;
00095 
00096   int nx = 128;
00097   int ny = comm.NumProc()*nx; // Scale y grid with number of processors
00098 
00099   // Create funky stencil to make sure the matrix is non-symmetric (transpose non-trivial):
00100 
00101   // (i-1,j-1) (i-1,j  )
00102   // (i  ,j-1) (i  ,j  ) (i  ,j+1)
00103   // (i+1,j-1) (i+1,j  )
00104 
00105   int npoints = 7;
00106 
00107   int xoff[] = {-1,  0,  1, -1,  0,  1,  0};
00108   int yoff[] = {-1, -1, -1,  0,  0,  0,  1};
00109 
00110   Epetra_Map * map;
00111   Epetra_CrsMatrix * A;
00112   Epetra_Vector * x, * b, * xexact;
00113   
00114   Trilinos_Util_GenerateCrsProblem(nx, ny, npoints, xoff, yoff, comm, map, A, x, b, xexact);
00115 
00116   if (nx<8)
00117   {
00118     cout << *A << endl;
00119     cout << "X exact = " << endl << *xexact << endl;
00120     cout << "B       = " << endl << *b << endl;
00121   }
00122 
00123   // Construct transposer 
00124   Epetra_Time timer(comm);
00125 
00126   double start = timer.ElapsedTime();
00127 
00128   //bool IgnoreNonLocalCols = false;
00129   bool MakeDataContiguous = true;
00130   EpetraExt::RowMatrix_Transpose transposer( MakeDataContiguous );
00131 
00132   if (verbose) cout << "\nTime to construct transposer  = " << timer.ElapsedTime() - start << endl;
00133   
00134   Epetra_CrsMatrix & transA = dynamic_cast<Epetra_CrsMatrix&>(transposer(*A));
00135 
00136   start = timer.ElapsedTime();
00137   if (verbose) cout << "\nTime to create transpose matrix  = " << timer.ElapsedTime() - start << endl;
00138   
00139   // Now test output of transposer by performing matvecs
00140   int ierr = 0;
00141   ierr += checkResults(A, &transA, xexact, verbose);
00142 
00143 
00144   // Now change values in original matrix and test update facility of transposer
00145   // Add 2 to the diagonal of each row
00146   double Value = 2.0;
00147   for (i=0; i< A->NumMyRows(); i++)
00148   A->SumIntoMyValues(i, 1, &Value, &i);
00149 
00150   start = timer.ElapsedTime();
00151   transposer.fwd();
00152 
00153   if (verbose) cout << "\nTime to update transpose matrix  = " << timer.ElapsedTime() - start << endl;
00154   
00155   ierr += checkResults(A, &transA, xexact, verbose);
00156 
00157   delete A;
00158   delete b;
00159   delete x;
00160   delete xexact;
00161   delete map;
00162 
00163   if (verbose) cout << endl << "Checking transposer for VbrMatrix objects" << endl<< endl;
00164 
00165   int nsizes = 4;
00166   int sizes[] = {4, 6, 5, 3};
00167 
00168   Epetra_VbrMatrix * Avbr;
00169   Epetra_BlockMap * bmap;
00170 
00171   Trilinos_Util_GenerateVbrProblem(nx, ny, npoints, xoff, yoff, nsizes, sizes,
00172                                    comm, bmap, Avbr, x, b, xexact);
00173 
00174   if (nx<8)
00175   {
00176     cout << *Avbr << endl;
00177     cout << "X exact = " << endl << *xexact << endl;
00178     cout << "B       = " << endl << *b << endl;
00179   }
00180 
00181   start = timer.ElapsedTime();
00182   EpetraExt::RowMatrix_Transpose transposer1( MakeDataContiguous );
00183 
00184   Epetra_CrsMatrix & transA1 = dynamic_cast<Epetra_CrsMatrix&>(transposer1(*Avbr));
00185   if (verbose) cout << "\nTime to create transpose matrix  = " << timer.ElapsedTime() - start << endl;
00186   
00187   // Now test output of transposer by performing matvecs
00188 ;
00189   ierr += checkResults(Avbr, &transA1, xexact, verbose);
00190 
00191   // Now change values in original matrix and test update facility of transposer
00192   // Scale matrix on the left by rowsums
00193 
00194   Epetra_Vector invRowSums(Avbr->RowMap());
00195 
00196   Avbr->InvRowSums(invRowSums);
00197   Avbr->LeftScale(invRowSums);
00198 
00199   start = timer.ElapsedTime();
00200   transposer1.fwd();
00201   if (verbose) cout << "\nTime to update transpose matrix  = " << timer.ElapsedTime() - start << endl;
00202   
00203   ierr += checkResults(Avbr, &transA1, xexact, verbose);
00204 
00205   delete Avbr;
00206   delete b;
00207   delete x;
00208   delete xexact;
00209   delete bmap;
00210 
00211 #ifdef EPETRA_MPI
00212   MPI_Finalize();
00213 #endif
00214 
00215   return ierr;
00216 }
00217 
00218 int checkResults(Epetra_RowMatrix * A, Epetra_CrsMatrix * transA, 
00219                  Epetra_Vector * xexact, bool verbose) {
00220 
00221   int n = A->NumGlobalRows();
00222 
00223   if (n<100) cout << "A transpose = " << endl << *transA << endl;
00224 
00225   Epetra_Vector x1(View,A->OperatorDomainMap(), &((*xexact)[0]));
00226   Epetra_Vector b1(A->OperatorRangeMap());
00227 
00228   A->SetUseTranspose(true);
00229 
00230   Epetra_Time timer(A->Comm());
00231   double start = timer.ElapsedTime();
00232   A->Apply(x1, b1);
00233   if (verbose) cout << "\nTime to compute b1: matvec with original matrix using transpose flag  = " << timer.ElapsedTime() - start << endl;
00234 
00235   if (n<100) cout << "b1 = " << endl << b1 << endl;
00236   Epetra_Vector x2(View,transA->OperatorRangeMap(), &((*xexact)[0]));
00237   Epetra_Vector b2(transA->OperatorDomainMap());
00238   start = timer.ElapsedTime();
00239   transA->Multiply(false, x2, b2);
00240   if (verbose) cout << "\nTime to compute b2: matvec with transpose matrix                      = " << timer.ElapsedTime() - start << endl;
00241 
00242   if (n<100) cout << "b1 = " << endl << b1 << endl;
00243 
00244   double residual;
00245   Epetra_Vector resid(A->OperatorDomainMap());
00246 
00247   resid.Update(1.0, b1, -1.0, b2, 0.0);
00248   resid.Norm2(&residual);
00249   if (verbose) cout << "Norm of b1 - b2 = " << residual << endl;
00250 
00251   int ierr = 0;
00252 
00253   if (residual > 1.0e-10) ierr++;
00254 
00255   if (ierr!=0 && verbose) cerr << "Status: Test failed" << endl;
00256   else if (verbose) cerr << "Status: Test passed" << endl;
00257 
00258   return(ierr);
00259 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines