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