Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/CrsRiluk/cxx_main.cpp
Go to the documentation of this file.
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack: Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2002) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 */
42 
43 #include "Ifpack_ConfigDefs.h"
44 #include <stdio.h>
45 #include <stdlib.h>
46 #include <ctype.h>
47 #include <assert.h>
48 #include <string.h>
49 #include <math.h>
50 #include "Epetra_Comm.h"
51 #include "Epetra_Map.h"
52 #include "Epetra_Time.h"
53 #include "Epetra_CrsMatrix.h"
54 #include "Epetra_Vector.h"
55 #include "Epetra_Object.h"
56 #include "Ifpack_Version.h"
57 #ifdef EPETRA_MPI
58 #include "mpi.h"
59 #include "Epetra_MpiComm.h"
60 #else
61 #include "Epetra_SerialComm.h"
62 #endif
63 
64 
65 // prototype
67  Epetra_Vector& q,
68  Epetra_Vector& z,
69  Epetra_Vector& resid,
70  double * lambda, int niters, double tolerance,
71  bool verbose);
72 
73 
74 int main(int argc, char *argv[])
75 {
76  using std::cout;
77  using std::endl;
78 
79  int ierr = 0, i;
80 
81 #ifdef EPETRA_MPI
82 
83  // Initialize MPI
84 
85  MPI_Init(&argc,&argv);
86  int rank; // My process ID
87 
88  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
89  Epetra_MpiComm Comm( MPI_COMM_WORLD );
90 
91 #else
92 
93  int rank = 0;
94  Epetra_SerialComm Comm;
95 
96 #endif
97 
98  bool verbose = false;
99 
100  // Check if we should print results to standard out
101  if (argc>1) if (argv[1][0]=='-' && argv[1][1]=='v') verbose = true;
102 
103  // char tmp;
104  // if (rank==0) cout << "Press any key to continue..."<< endl;
105  // if (rank==0) cin >> tmp;
106  // Comm.Barrier();
107 
108  int MyPID = Comm.MyPID();
109  int NumProc = Comm.NumProc();
110 
111  if (verbose && MyPID==0)
112  cout << Ifpack_Version() << endl << endl;
113 
114  if (verbose) cout << "Processor "<<MyPID<<" of "<< NumProc
115  << " is alive."<<endl;
116 
117  bool verbose1 = verbose;
118 
119  // Redefine verbose to only print on PE 0
120  if (verbose && rank!=0) verbose = false;
121 
122  int NumMyEquations = 10000;
123  int NumGlobalEquations = NumMyEquations*NumProc+EPETRA_MIN(NumProc,3);
124  if (MyPID < 3) NumMyEquations++;
125 
126  // Construct a Map that puts approximately the same Number of equations on each processor
127 
128  Epetra_Map Map(NumGlobalEquations, NumMyEquations, 0, Comm);
129 
130  // Get update list and number of local equations from newly created Map
131  int * MyGlobalElements = new int[Map.NumMyElements()];
132  Map.MyGlobalElements(MyGlobalElements);
133 
134  // Create an integer vector NumNz that is used to build the Petra Matrix.
135  // NumNz[i] is the Number of OFF-DIAGONAL term for the ith global equation on this processor
136 
137  int * NumNz = new int[NumMyEquations];
138 
139  // We are building a tridiagonal matrix where each row has (-1 2 -1)
140  // So we need 2 off-diagonal terms (except for the first and last equation)
141 
142  for (i=0; i<NumMyEquations; i++)
143  if (MyGlobalElements[i]==0 || MyGlobalElements[i] == NumGlobalEquations-1)
144  NumNz[i] = 1;
145  else
146  NumNz[i] = 2;
147 
148  // Create a Epetra_Matrix
149 
150  Epetra_CrsMatrix A(Copy, Map, NumNz);
151 
152  // Add rows one-at-a-time
153  // Need some vectors to help
154  // Off diagonal Values will always be -1
155 
156 
157  double *Values = new double[2];
158  Values[0] = -1.0; Values[1] = -1.0;
159  int *Indices = new int[2];
160  double two = 2.0;
161  int NumEntries;
162 
163  for (i=0; i<NumMyEquations; i++)
164  {
165  if (MyGlobalElements[i]==0)
166  {
167  Indices[0] = 1;
168  NumEntries = 1;
169  }
170  else if (MyGlobalElements[i] == NumGlobalEquations-1)
171  {
172  Indices[0] = NumGlobalEquations-2;
173  NumEntries = 1;
174  }
175  else
176  {
177  Indices[0] = MyGlobalElements[i]-1;
178  Indices[1] = MyGlobalElements[i]+1;
179  NumEntries = 2;
180  }
181  int ierr2;
182  ierr2 = A.InsertGlobalValues(MyGlobalElements[i], NumEntries, Values, Indices);
183  IFPACK_CHK_ERR(ierr2);
184  ierr2 = A.InsertGlobalValues(MyGlobalElements[i], 1, &two, MyGlobalElements+i); // Put in the diagonal entry
185  IFPACK_CHK_ERR(ierr2);
186  }
187 
188  // Finish up
189  A.FillComplete();
190 
191  // Create vectors for Power method
192 
193  Epetra_Vector q(Map);
194  Epetra_Vector z(Map);
195  Epetra_Vector resid(Map);
196 
197  // variable needed for iteration
198  double lambda = 0.0;
199  int niters = 10000;
200  double tolerance = 1.0e-3;
201 
202  // Iterate
203  Epetra_Time timer(Comm);
204  ierr += power_method(A, q, z, resid, &lambda, niters, tolerance, verbose);
205  double elapsed_time = timer.ElapsedTime();
206  double total_flops = A.Flops() + q.Flops() + z.Flops() + resid.Flops();
207  double MFLOPs = total_flops/elapsed_time/1000000.0;
208 
209  if (verbose) cout << "\n\nTotal MFLOPs for first solve = " << MFLOPs << endl<< endl;
210 
211  // Increase diagonal dominance
212  if (verbose) cout << "\n\nIncreasing the magnitude of first diagonal term and solving again\n\n"
213  << endl;
214 
215 
216  if (A.MyGlobalRow(0)) {
217  int numvals = A.NumGlobalEntries(0);
218  double * Rowvals = new double [numvals];
219  int * Rowinds = new int [numvals];
220  A.ExtractGlobalRowCopy(0, numvals, numvals, Rowvals, Rowinds); // Get A[0,0]
221 
222  for (i=0; i<numvals; i++) if (Rowinds[i] == 0) Rowvals[i] *= 10.0;
223 
224  A.ReplaceGlobalValues(0, numvals, Rowvals, Rowinds);
225 
226  delete [] Rowvals;
227  delete [] Rowinds;
228  }
229  // Iterate (again)
230  lambda = 0.0;
231  timer.ResetStartTime();
232  A.ResetFlops(); q.ResetFlops(); z.ResetFlops(); resid.ResetFlops();
233  ierr += power_method(A, q, z, resid, &lambda, niters, tolerance, verbose);
234  elapsed_time = timer.ElapsedTime();
235  total_flops = A.Flops() + q.Flops() + z.Flops() + resid.Flops();
236  MFLOPs = total_flops/elapsed_time/1000000.0;
237 
238  if (verbose) cout << "\n\nTotal MFLOPs for second solve = " << MFLOPs << endl<< endl;
239 
240 
241  // Release all objects
242  delete [] NumNz;
243  delete [] Values;
244  delete [] Indices;
245  delete [] MyGlobalElements;
246 
247 
248 
249  if (verbose1) {
250  // Test ostream << operator (if verbose1)
251  // Construct a Map that puts 2 equations on each PE
252 
253  int NumMyElements1 = 2;
254  int NumMyEquations1 = NumMyElements1;
255  int NumGlobalEquations1 = NumMyEquations1*NumProc;
256 
257  Epetra_Map Map1(-1, NumMyElements1, 0, Comm);
258 
259  // Get update list and number of local equations from newly created Map
260  int * MyGlobalElements1 = new int[Map1.NumMyElements()];
261  Map1.MyGlobalElements(MyGlobalElements1);
262 
263  // Create an integer vector NumNz that is used to build the Petra Matrix.
264  // NumNz[i] is the Number of OFF-DIAGONAL term for the ith global equation on this processor
265 
266  int * NumNz1 = new int[NumMyEquations1];
267 
268  // We are building a tridiagonal matrix where each row has (-1 2 -1)
269  // So we need 2 off-diagonal terms (except for the first and last equation)
270 
271  for (i=0; i<NumMyEquations1; i++)
272  if (MyGlobalElements1[i]==0 || MyGlobalElements1[i] == NumGlobalEquations1-1)
273  NumNz1[i] = 1;
274  else
275  NumNz1[i] = 2;
276 
277  // Create a Epetra_Matrix
278 
279  Epetra_CrsMatrix A1(Copy, Map1, NumNz1);
280 
281  // Add rows one-at-a-time
282  // Need some vectors to help
283  // Off diagonal Values will always be -1
284 
285 
286  double *Values1 = new double[2];
287  Values1[0] = -1.0; Values1[1] = -1.0;
288  int *Indices1 = new int[2];
289  double two1 = 2.0;
290  int NumEntries1;
291 
292  for (i=0; i<NumMyEquations1; i++)
293  {
294  if (MyGlobalElements1[i]==0)
295  {
296  Indices1[0] = 1;
297  NumEntries1 = 1;
298  }
299  else if (MyGlobalElements1[i] == NumGlobalEquations1-1)
300  {
301  Indices1[0] = NumGlobalEquations1-2;
302  NumEntries1 = 1;
303  }
304  else
305  {
306  Indices1[0] = MyGlobalElements1[i]-1;
307  Indices1[1] = MyGlobalElements1[i]+1;
308  NumEntries1 = 2;
309  }
310  int ierr2;
311  ierr2 = A1.InsertGlobalValues(MyGlobalElements1[i], NumEntries1, Values1, Indices1);
312  IFPACK_CHK_ERR(ierr2);
313  ierr2 = A1.InsertGlobalValues(MyGlobalElements1[i], 1, &two1, MyGlobalElements1+i); // Put in the diagonal entry
314  IFPACK_CHK_ERR(ierr2);
315  }
316 
317  // Finish up
318  A1.FillComplete();
319 
320  if (verbose) cout << "\n\nPrint out tridiagonal matrix, each part on each processor.\n\n" << endl;
321  cout << A1 << endl;
322 
323  // Release all objects
324  delete [] NumNz1;
325  delete [] Values1;
326  delete [] Indices1;
327  delete [] MyGlobalElements1;
328 
329  }
330 
331 #ifdef EPETRA_MPI
332  MPI_Finalize() ;
333 #endif
334 
335 /* end main
336 */
337 return ierr ;
338 }
339 
341  Epetra_Vector& q,
342  Epetra_Vector& z,
343  Epetra_Vector& resid,
344  double * lambda, int niters, double tolerance,
345  bool verbose) {
346  using std::cout;
347  using std::endl;
348 
349  // Fill z with random Numbers
350  z.Random();
351 
352  // variable needed for iteration
353  double normz, residual;
354 
355  int ierr = 1;
356 
357  for (int iter = 0; iter < niters; iter++)
358  {
359  z.Norm2(&normz); // Compute 2-norm of z
360  q.Scale(1.0/normz, z);
361  A.Multiply(false, q, z); // Compute z = A*q
362  q.Dot(z, lambda); // Approximate maximum eigenvaluE
363  if (iter%100==0 || iter+1==niters)
364  {
365  resid.Update(1.0, z, -(*lambda), q, 0.0); // Compute A*q - lambda*q
366  resid.Norm2(&residual);
367  if (verbose) cout << "Iter = " << iter << " Lambda = " << *lambda
368  << " Residual of A*q - lambda*q = " << residual << endl;
369  }
370  if (residual < tolerance) {
371  ierr = 0;
372  break;
373  }
374  }
375  return(ierr);
376 }
bool MyGlobalRow(int GID) const
int NumGlobalEntries(long long Row) const
int power_method(Epetra_CrsMatrix &A, Epetra_Vector &q, Epetra_Vector &z, Epetra_Vector &resid, double *lambda, int niters, double tolerance, bool verbose)
int ExtractGlobalRowCopy(int_type Row, int Length, int &NumEntries, double *values, int_type *Indices) const
int Multiply(bool TransA, const Epetra_Vector &x, Epetra_Vector &y) const
int MyGlobalElements(int *MyGlobalElementList) const
double ElapsedTime(void) const
virtual int InsertGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)
#define EPETRA_MIN(x, y)
int MyPID() const
int FillComplete(bool OptimizeDataStorage=true)
int NumMyElements() const
std::string Ifpack_Version()
int main(int argc, char *argv[])
int NumProc() const
TransListIter iter
#define IFPACK_CHK_ERR(ifpack_err)
void ResetStartTime(void)
virtual int ReplaceGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)