Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Test_SuperLU_DIST/cxx_main.cpp
Go to the documentation of this file.
1 #include "Amesos_ConfigDefs.h"
2 
3 #ifdef HAVE_MPI
4 #include "mpi.h"
5 #include "Epetra_MpiComm.h"
6 #else
7 #include "Epetra_SerialComm.h"
8 #endif
9 #include "Epetra_Map.h"
10 #include "Epetra_Vector.h"
11 #include "Epetra_Time.h"
12 #include "Epetra_Util.h"
13 #include "Epetra_CrsMatrix.h"
14 #include "Amesos_Superludist.h"
15 #include "Amesos_TestRowMatrix.h"
17 
18 //=============================================================================
19 bool CheckError(bool verbose,
20  const Epetra_RowMatrix& A,
21  const Epetra_MultiVector& x,
22  const Epetra_MultiVector& b,
23  const Epetra_MultiVector& x_exact)
24 {
25  std::vector<double> Norm;
26  int NumVectors = x.NumVectors();
27  Norm.resize(NumVectors);
28  Epetra_MultiVector Ax(x);
29  A.Multiply(false,x,Ax);
30  Ax.Update(1.0,b,-1.0);
31  Ax.Norm2(&Norm[0]);
32  bool TestPassed = false;
33  double TotalNorm = 0.0;
34  for (int i = 0 ; i < NumVectors ; ++i) {
35  TotalNorm += Norm[i];
36  }
37  if (verbose && A.Comm().MyPID() == 0)
38  std::cout << "||Ax - b|| = " << TotalNorm << std::endl;
39  if (TotalNorm < 1e-5 )
40  TestPassed = true;
41  else
42  TestPassed = false;
43 
44  Ax.Update (1.0,x,-1.0,x_exact,0.0);
45  Ax.Norm2(&Norm[0]);
46  for (int i = 0 ; i < NumVectors ; ++i) {
47  TotalNorm += Norm[i];
48  }
49  if (verbose && A.Comm().MyPID() == 0)
50  std::cout << "||Ax - b|| = " << TotalNorm << std::endl;
51  if (TotalNorm < 1e-5 )
52  TestPassed = true;
53  else
54  TestPassed = false;
55 
56  return(TestPassed);
57 }
58 
59 
60 //=============================================================================
61 int main(int argc, char *argv[]) {
62 
63 #ifdef HAVE_MPI
64  MPI_Init(&argc, &argv);
65  Epetra_MpiComm Comm(MPI_COMM_WORLD);
66 #else
67  Epetra_SerialComm Comm;
68 #endif
69 
70  int NumGlobalElements = 1000; // kludge see bug #1142 - when set to 7, and Min_jGlobal is set to zero,
71  // Test_SuperLU_DIST.exe dies during Amesos_Superludist destruction
72  int NumVectors = 7;
73 
74  bool verbose = true ;
75  if ( argc > 1 && argv[1][0] == '-' && argv[1][1] == 'q' ) verbose = false ;
76 
77  // =================== //
78  // create a random map //
79  // =================== //
80 
81  int* part = new int[NumGlobalElements];
82 
83  if (Comm.MyPID() == 0) {
84  Epetra_Util Util;
85 
86  for( int i=0 ; i<NumGlobalElements ; ++i ) {
87  unsigned int r = Util.RandomInt();
88  part[i] = r%(Comm.NumProc());
89  }
90  }
91 
92  Comm.Broadcast(part,NumGlobalElements,0);
93 
94  // count the elements assigned to this proc
95  int NumMyElements = 0;
96  for (int i = 0 ; i < NumGlobalElements ; ++i) {
97  if (part[i] == Comm.MyPID())
98  NumMyElements++;
99  }
100 
101  // get the loc2global list
102  int* MyGlobalElements = new int[NumMyElements];
103  int count = 0;
104  for (int i = 0 ; i < NumGlobalElements ; ++i) {
105  if (part[i] == Comm.MyPID() )
106  MyGlobalElements[count++] = i;
107  }
108 
109  Epetra_Map Map(NumGlobalElements,NumMyElements,MyGlobalElements,
110  0,Comm);
111 
112  delete [] part;
113 
114  // ===================== //
115  // Create a dense matrix //
116  // ===================== //
117 
118  Epetra_CrsMatrix Matrix(Copy,Map,NumGlobalElements);
119 
120  int* Indices = new int[NumGlobalElements];
121  double* Values = new double[NumGlobalElements];
122 
123  for (int i = 0 ; i < NumGlobalElements ; ++i)
124  Indices[i] = i;
125 
126  for (int i = 0 ; i < NumMyElements ; ++i)
127  {
128  int iGlobal = MyGlobalElements[i];
129  const int MakeNotDense = 1; // kludge see bug #1142 - set to zero to demonstrate bug #1142 on atlantis
130  int Min_jGlobal = std::min(i,MakeNotDense );
131  for (int jGlobal = Min_jGlobal ; jGlobal < NumGlobalElements ; ++jGlobal) {
132  if (iGlobal == jGlobal)
133  Values[jGlobal-Min_jGlobal] = 1.0 * (NumGlobalElements + 1 ) *
134  (NumGlobalElements + 1);
135  else if (iGlobal > jGlobal)
136  Values[jGlobal-Min_jGlobal] = -1.0*(jGlobal+1);
137  else
138  Values[jGlobal-Min_jGlobal] = 1.0*(iGlobal+1);
139  }
140  Matrix.InsertGlobalValues(MyGlobalElements[i],
141  NumGlobalElements-MakeNotDense, Values, &Indices[Min_jGlobal]);
142 
143  }
144 
145  Matrix.FillComplete();
146 
147  delete[] MyGlobalElements;
148  delete[] Indices;
149  delete[] Values;
150 
151  // ======================== //
152  // other data for this test //
153  // ======================== //
154 
155  Amesos_TestRowMatrix A(&Matrix);
156  Epetra_MultiVector x(Map,NumVectors);
157  Epetra_MultiVector x_exact(Map,NumVectors);
158  Epetra_MultiVector b(Map,NumVectors);
159  x_exact.Random();
160  A.Multiply(false,x_exact,b);
161 
162  // =========== //
163  // AMESOS PART //
164  // =========== //
165 
166  Epetra_LinearProblem Problem(&Matrix, &x, &b);
168 
170  List.set("MaxProcs", Comm.NumProc());
171 
172  Solver->SetParameters(List);
173 
176  AMESOS_CHK_ERR(Solver->Solve());
177 
178  if (!CheckError(verbose,A,x,b,x_exact))
179  AMESOS_CHK_ERR(-1);
180 
181  delete Solver;
182 
183 #ifdef HAVE_MPI
184  MPI_Finalize();
185 #endif
186  return(EXIT_SUCCESS);
187 }
int Solve()
Solves A X = B (or AT x = B)
Amesos_Superludist: An object-oriented wrapper for Superludist.
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
virtual int InsertGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)
static bool verbose
Definition: Amesos.cpp:67
Amesos_TestRowMatrix: a class to test Epetra_RowMatrix based codes.
int SymbolicFactorization()
Performs SymbolicFactorization on the matrix A.
int MyPID() const
virtual int MyPID() const =0
int FillComplete(bool OptimizeDataStorage=true)
int SetParameters(Teuchos::ParameterList &ParameterList)
Updates internal variables.
unsigned int RandomInt()
virtual int Multiply(bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Epetra_RowMatrix multiplied by a Epetra_MultiVector X in Y.
virtual const Epetra_Comm & Comm() const =0
int NumericFactorization()
Performs NumericFactorization on the matrix A.
#define AMESOS_CHK_ERR(a)
int main(int argc, char *argv[])
int NumProc() const
virtual int Multiply(bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0
bool CheckError(const std::string SolverType, const std::string Descriptor, const Epetra_RowMatrix &A, const Epetra_MultiVector &x, const Epetra_MultiVector &b, const Epetra_MultiVector &x_exact)
int Broadcast(double *MyVals, int Count, int Root) const