Amesos Package Browser (Single Doxygen Collection) Development
Test_SuperLU/cxx_main.cpp
Go to the documentation of this file.
00001 #include "Amesos_ConfigDefs.h"
00002 
00003 #ifdef HAVE_MPI
00004 #include "mpi.h"
00005 #include "Epetra_MpiComm.h"
00006 #else
00007 #include "Epetra_SerialComm.h"
00008 #endif
00009 #include "Epetra_Map.h"
00010 #include "Epetra_Vector.h"
00011 #include "Epetra_Util.h"
00012 #include "Epetra_CrsMatrix.h"
00013 #include "Amesos_Superlu.h"
00014 #include "Amesos_TestRowMatrix.h"
00015 #include "Teuchos_ParameterList.hpp"
00016 
00017 //=============================================================================
00018 bool CheckError(const Epetra_RowMatrix& A,
00019     const Epetra_MultiVector& x,
00020     const Epetra_MultiVector& b,
00021     const Epetra_MultiVector& x_exact)
00022 {
00023   std::vector<double> Norm;
00024   int NumVectors = x.NumVectors();
00025   Norm.resize(NumVectors);
00026   Epetra_MultiVector Ax(x);
00027   A.Multiply(false,x,Ax);
00028   Ax.Update(1.0,b,-1.0);
00029   Ax.Norm2(&Norm[0]);
00030   bool TestPassed = false;
00031   double TotalNorm = 0.0;
00032   for (int i = 0 ; i < NumVectors ; ++i) {
00033     TotalNorm += Norm[i];
00034   }
00035   if (A.Comm().MyPID() == 0)
00036     std::cout << "||Ax - b||  = " << TotalNorm << std::endl;
00037   if (TotalNorm < 1e-5 )
00038     TestPassed = true;
00039   else
00040     TestPassed = false;
00041 
00042   Ax.Update (1.0,x,-1.0,x_exact,0.0);
00043   Ax.Norm2(&Norm[0]);
00044   for (int i = 0 ; i < NumVectors ; ++i) {
00045     TotalNorm += Norm[i];
00046   }
00047   if (A.Comm().MyPID() == 0)
00048     std::cout << "||Ax - b||  = " << TotalNorm << std::endl;
00049   if (TotalNorm < 1e-5 )
00050     TestPassed = true;
00051   else
00052     TestPassed = false;
00053 
00054   return(TestPassed);
00055 }
00056 
00057 //=============================================================================
00058 int main(int argc, char *argv[]) 
00059 {
00060 #ifdef HAVE_MPI
00061   MPI_Init(&argc, &argv);
00062   Epetra_MpiComm Comm(MPI_COMM_WORLD);
00063 #else
00064   Epetra_SerialComm Comm;
00065 #endif
00066 
00067   int NumGlobalElements = 1000;
00068   int NumVectors = 7;
00069   
00070   // =================== //
00071   // create a random map //
00072   // =================== //
00073  
00074   int* part = new int[NumGlobalElements];
00075 
00076   if (Comm.MyPID() == 0) {
00077     Epetra_Util Util;
00078 
00079     for( int i=0 ; i<NumGlobalElements ; ++i ) {
00080       unsigned int r = Util.RandomInt();  
00081       part[i] = r%(Comm.NumProc());
00082     }
00083   }
00084 
00085   Comm.Broadcast(part,NumGlobalElements,0);
00086 
00087   // count the elements assigned to this proc
00088   int NumMyElements = 0;
00089   for (int i = 0 ; i < NumGlobalElements ; ++i) {
00090     if (part[i] == Comm.MyPID()) 
00091       NumMyElements++;
00092   }
00093 
00094   // get the loc2global list
00095   int* MyGlobalElements = new int[NumMyElements];
00096   int count = 0;
00097   for (int i = 0 ; i < NumGlobalElements ; ++i) {
00098     if (part[i] == Comm.MyPID() ) 
00099       MyGlobalElements[count++] = i;
00100   }
00101 
00102   Epetra_Map Map(NumGlobalElements,NumMyElements,MyGlobalElements,
00103      0,Comm);
00104 
00105   delete [] part;
00106 
00107   // ===================== //
00108   // Create a dense matrix //
00109   // ===================== //
00110  
00111   Epetra_CrsMatrix Matrix(Copy,Map,NumGlobalElements);
00112 
00113   int* Indices = new int[NumGlobalElements];
00114   double* Values = new double[NumGlobalElements];
00115 
00116   for (int i = 0 ; i < NumGlobalElements ; ++i) 
00117     Indices[i] = i;
00118 
00119   for (int i = 0 ; i < NumMyElements ; ++i) {
00120     int iGlobal = MyGlobalElements[i];
00121     for (int jGlobal = 0 ; jGlobal < NumGlobalElements ; ++jGlobal) {
00122       if (iGlobal == jGlobal) 
00123   Values[jGlobal] = 1.0 * (NumGlobalElements + 1 ) *
00124     (NumGlobalElements + 1);
00125       else if (iGlobal > jGlobal)
00126   Values[jGlobal] = -1.0*(jGlobal+1);
00127       else
00128   Values[jGlobal] = 1.0*(iGlobal+1);
00129     }
00130     Matrix.InsertGlobalValues(MyGlobalElements[i],
00131                               NumGlobalElements, Values, Indices);
00132 
00133   }
00134 
00135   Matrix.FillComplete();
00136 
00137   delete [] MyGlobalElements;
00138   delete [] Indices;
00139   delete [] Values;
00140  
00141   // ======================== //
00142   // other data for this test //
00143   // ======================== //
00144 
00145   Amesos_TestRowMatrix A(&Matrix);
00146   Epetra_MultiVector x(Map,NumVectors);
00147   Epetra_MultiVector x_exact(Map,NumVectors);
00148   Epetra_MultiVector b(Map,NumVectors);
00149   x_exact.Random();
00150   A.Multiply(false,x_exact,b);
00151 
00152   // =========== //
00153   // AMESOS PART //
00154   // =========== //
00155 
00156   Epetra_LinearProblem Problem;
00157   Amesos_Superlu Solver(Problem);
00158 
00159   Problem.SetOperator(&A);
00160   Problem.SetLHS(&x);
00161   Problem.SetRHS(&b);
00162 
00163   AMESOS_CHK_ERR(Solver.SymbolicFactorization());
00164   AMESOS_CHK_ERR(Solver.NumericFactorization());
00165   AMESOS_CHK_ERR(Solver.Solve());
00166 
00167   bool TestPassed = true;
00168 
00169   TestPassed = TestPassed &&
00170     CheckError(A,x,b,x_exact);
00171 
00172   AMESOS_CHK_ERR( ! TestPassed ) ; 
00173 
00174 #ifdef HAVE_MPI
00175   MPI_Finalize();
00176 #endif
00177 
00178 
00179   return(EXIT_SUCCESS);
00180 
00181 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines