Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
example/InverseIteration/cxx_main.cpp
Go to the documentation of this file.
1 //@HEADER
2 // ************************************************************************
3 //
4 // Epetra: Linear Algebra Services Package
5 // Copyright 2011 Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
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 #include <stdio.h>
43 #include <stdlib.h>
44 #include <ctype.h>
45 #include <assert.h>
46 #include <string.h>
47 #include <math.h>
48 #include "Epetra_Comm.h"
49 #include "Epetra_Map.h"
50 #include "Epetra_Time.h"
51 #include "Epetra_BlockMap.h"
52 #include "Epetra_MultiVector.h"
53 #include "Epetra_Vector.h"
54 #include "Epetra_CrsMatrix.h"
55 #ifdef EPETRA_MPI
56 #include "mpi.h"
57 #include "Epetra_MpiComm.h"
58 #else
59 #include "Epetra_SerialComm.h"
60 #endif
61 #include "Trilinos_Util.h"
62 #include "Ifpack_IlukGraph.h"
63 #include "Ifpack_CrsRiluk.h"
64 
65 // *************************************************************
66 // Function prototypes
67 // *************************************************************
68 
69 int invIteration(Epetra_CrsMatrix& A, double &lambda, bool verbose);
70 int applyInverseSetup(Epetra_CrsMatrix &A, Ifpack_CrsRiluk * & M);
71 int applyInverse(Epetra_CrsMatrix &A, Epetra_Vector & x, Epetra_Vector & b, Ifpack_CrsRiluk * M,
72  bool verbose);
73 int applyInverseDestroy(Ifpack_CrsRiluk * M);
75  Epetra_Vector &x,
76  Epetra_Vector &b,
77  Ifpack_CrsRiluk *M,
78  int Maxiter,
79  double Tolerance, bool verbose);
80 
81 // *************************************************************
82 // main program - This benchmark code reads a Harwell-Boeing data
83 // set and finds the minimal eigenvalue of the matrix
84 // using inverse iteration.
85 // *************************************************************
86 int main(int argc, char *argv[]) {
87 
88 #ifdef EPETRA_MPI
89  MPI_Init(&argc,&argv);
90  Epetra_MpiComm Comm (MPI_COMM_WORLD);
91 #else
92  Epetra_SerialComm Comm;
93 #endif
94 
95  cout << Comm << endl;
96 
97  int MyPID = Comm.MyPID();
98 
99  bool verbose = false;
100  if (MyPID==0) verbose = true; // Print out detailed results (turn off for best performance)
101 
102  if(argc != 2) {
103  if (verbose) cerr << "Usage: " << argv[0] << " HB_data_file" << endl;
104  exit(1); // Error
105  }
106 
107  // Define pointers that will be set by HB read function
108 
109  Epetra_Map * readMap;
110  Epetra_CrsMatrix * readA;
111  Epetra_Vector * readx;
112  Epetra_Vector * readb;
113  Epetra_Vector * readxexact;
114 
115  // Call function to read in HB problem
116  Trilinos_Util_ReadHb2Epetra(argv[1], Comm, readMap, readA, readx, readb, readxexact);
117 
118  // Not interested in x, b or xexact for an eigenvalue problem
119  delete readx;
120  delete readb;
121  delete readxexact;
122 
123 #ifdef EPETRA_MPI // If running in parallel, we need to distribute matrix across all PEs.
124 
125  // Create uniform distributed map
126  Epetra_Map map(readMap->NumGlobalElements(), 0, Comm);
127 
128  // Create Exporter to distribute read-in matrix and vectors
129 
130  Epetra_Export exporter(*readMap, map);
131  Epetra_CrsMatrix A(Copy, map, 0);
132 
133  A.Export(*readA, exporter, Add);
134  assert(A.FillComplete()==0);
135 
136  delete readA;
137  delete readMap;
138 
139 #else // If not running in parallel, we do not need to distribute the matrix
140  Epetra_CrsMatrix & A = *readA;
141 #endif
142 
143  // Create flop counter to collect all FLOPS
144  Epetra_Flops counter;
145  A.SetFlopCounter(counter);
146 
147  double lambda = 0; // Minimal eigenvalue returned here
148  // Call inverse iteration solver
149  Epetra_Time timer(Comm);
150  invIteration(A, lambda, verbose);
151  double elapsedTime = timer.ElapsedTime();
152  double totalFlops = counter.Flops();
153  double MFLOPS = totalFlops/elapsedTime/1000000.0;
154 
155 
156  cout << endl
157  << "*************************************************" << endl
158  << " Approximate smallest eigenvalue = " << lambda << endl
159  << " Total Time = " << elapsedTime << endl
160  << " Total FLOPS = " << totalFlops << endl
161  << " Total MFLOPS = " << MFLOPS << endl
162  << "*************************************************" << endl;
163 
164  // All done
165 #ifdef EPETRA_MPI
166  MPI_Finalize();
167 #else
168  delete readA;
169  delete readMap;
170 #endif
171 
172 return (0);
173 }
174 
175 // The following functions are used to solve the problem:
176 
177 
178 // *************************************************************
179 // Computes the smallest eigenvalue of the given matrix A.
180 // Returns result in lambda
181 // *************************************************************
182 
183 int invIteration(Epetra_CrsMatrix& A, double &lambda, bool verbose) {
184 
185  Ifpack_CrsRiluk * M;
186  applyInverseSetup(A, M);
187  Epetra_Vector q(A.RowMap());
188  Epetra_Vector z(A.RowMap());
189  Epetra_Vector resid(A.RowMap());
190 
191  Epetra_Flops * counter = A.GetFlopCounter();
192  if (counter!=0) {
193  q.SetFlopCounter(A);
194  z.SetFlopCounter(A);
195  resid.SetFlopCounter(A);
196  }
197 
198  // Fill z with random Numbers
199  z.Random();
200 
201  // variable needed for iteration
202  double normz, residual;
203 
204  int niters = 100;
205  double tolerance = 1.0E-6;
206  int ierr = 1;
207 
208  for (int iter = 0; iter < niters; iter++)
209  {
210  if (verbose)
211  cout << endl
212  << " ***** Performing step " << iter << " of inverse iteration ***** " << endl;
213 
214  z.Norm2(&normz); // Compute 2-norm of z
215  q.Scale(1.0/normz, z);
216  applyInverse(A, z, q, M, verbose); // Compute z such that Az = q
217  q.Dot(z, &lambda); // Approximate maximum eigenvalue
218  if (iter%10==0 || iter+1==niters)
219  {
220  resid.Update(1.0, z, -lambda, q, 0.0); // Compute A(inv)*q - lambda*q
221  resid.Norm2(&residual);
222  cout << endl
223  << "***** Inverse Iteration Step " << iter+1 << endl
224  << " Lambda = " << 1.0/lambda << endl
225  << " Residual of A(inv)*q - lambda*q = "
226  << residual << endl;
227  }
228  if (residual < tolerance) {
229  ierr = 0;
230  break;
231  }
232  }
233  // lambda is the largest eigenvalue of A(inv). 1/lambda is smallest eigenvalue of A.
234  lambda = 1.0/lambda;
235 
236  // Compute A*q - lambda*q explicitly
237  A.Multiply(false, q, z);
238  resid.Update(1.0, z, -lambda, q, 0.0); // Compute A*q - lambda*q
239  resid.Norm2(&residual);
240  cout << " Explicitly computed residual of A*q - lambda*q = " << residual << endl;
242  return(ierr);
243 }
244 
245 // Performs any setup that can be done once to reduce the cost of the applyInverse function
246 int applyInverseSetup(Epetra_CrsMatrix &A, Ifpack_CrsRiluk * & M) {
247  int LevelFill = 4;
248  int Overlap = 0;
249  Ifpack_IlukGraph * IlukGraph = new Ifpack_IlukGraph(A.Graph(), LevelFill, Overlap);
250  assert(IlukGraph->ConstructFilledGraph()==0);
251  M = new Ifpack_CrsRiluk(A, *IlukGraph);
252  M->SetFlopCounter(A);
253  assert(M->InitValues()==0);
254  assert(M->Factor()==0);
255  return(0);
256 }
257 
258 // *************************************************************
259 // Solves a problem Ax = b, for a given A and b.
260 // M is a preconditioner computed in applyInverseSetup.
261 // *************************************************************
262 
264  Epetra_Vector & x, Epetra_Vector & b, Ifpack_CrsRiluk * M, bool verbose) {
265 
266  int Maxiter = 1000;
267  double Tolerance = 1.0E-16;
268  BiCGSTAB(A, x, b, M, Maxiter, Tolerance, verbose);
269 
270  return(0);
271 }
272 
273 
274 // *************************************************************
275 // Releases any memory associated with the preconditioner M
276 // *************************************************************
277 
278 int applyInverseDestroy(Ifpack_CrsRiluk * M) {
279  Ifpack_IlukGraph & IlukGraph = (Ifpack_IlukGraph &) M->Graph();
280  delete M;
281  delete &IlukGraph;
282  return(0);
283 }
284 
285 // *************************************************************
286 // Uses the Bi-CGSTAB iterative method to solve Ax = b
287 // *************************************************************
288 
290  Epetra_Vector &x,
291  Epetra_Vector &b,
292  Ifpack_CrsRiluk *M,
293  int Maxiter,
294  double Tolerance, bool verbose) {
295 
296  // Allocate vectors needed for iterations
297  Epetra_Vector phat(x.Map()); phat.SetFlopCounter(x);
298  Epetra_Vector p(x.Map()); p.SetFlopCounter(x);
299  Epetra_Vector shat(x.Map()); shat.SetFlopCounter(x);
300  Epetra_Vector s(x.Map()); s.SetFlopCounter(x);
301  Epetra_Vector r(x.Map()); r.SetFlopCounter(x);
302  Epetra_Vector rtilde(x.Map()); rtilde.Random(); rtilde.SetFlopCounter(x);
303  Epetra_Vector v(x.Map()); v.SetFlopCounter(x);
304 
305 
306  A.Multiply(false, x, r); // r = A*x
307 
308  r.Update(1.0, b, -1.0); // r = b - A*x
309 
310  double r_norm, b_norm, scaled_r_norm, rhon, rhonm1 = 1.0;
311  double alpha = 1.0, omega = 1.0, sigma;
312  double omega_num, omega_den;
313  r.Norm2(&r_norm);
314  b.Norm2(&b_norm);
315  scaled_r_norm = r_norm/b_norm;
316  r.Dot(rtilde,&rhon);
317 
318  if (verbose) cout << "Initial residual = " << r_norm
319  << " Scaled residual = " << scaled_r_norm << endl;
320 
321 
322  for (int i=0; i<Maxiter; i++) { // Main iteration loop
323 
324  double beta = (rhon/rhonm1) * (alpha/omega);
325  rhonm1 = rhon;
326 
327  /* p = r + beta*(p - omega*v) */
328  /* phat = M^-1 p */
329  /* v = A phat */
330 
331  double dtemp = - beta*omega;
332 
333  p.Update(1.0, r, dtemp, v, beta);
334  if (M==0)
335  phat.Scale(1.0, p);
336  else
337  M->Solve(false, p, phat);
338  A.Multiply(false, phat, v);
339 
340 
341  rtilde.Dot(v,&sigma);
342  alpha = rhon/sigma;
343 
344  /* s = r - alpha*v */
345  /* shat = M^-1 s */
346  /* r = A shat (r is a tmp here for t ) */
347 
348  s.Update(-alpha, v, 1.0, r, 0.0);
349  if (M==0)
350  shat.Scale(1.0, s);
351  else
352  M->Solve(false, s, shat);
353  A.Multiply(false, shat, r);
354 
355  r.Dot(s, &omega_num);
356  r.Dot(r, &omega_den);
357  omega = omega_num / omega_den;
358 
359  /* x = x + alpha*phat + omega*shat */
360  /* r = s - omega*r */
361 
362  x.Update(alpha, phat, omega, shat, 1.0);
363  r.Update(1.0, s, -omega);
364 
365  r.Norm2(&r_norm);
366  scaled_r_norm = r_norm/b_norm;
367  r.Dot(rtilde,&rhon);
368 
369  if (verbose) cout << "Iter "<< i << " residual = " << r_norm
370  << " Scaled residual = " << scaled_r_norm << endl;
371 
372  if (scaled_r_norm < Tolerance) break;
373  }
374  return;
375 }
void SetFlopCounter(const Epetra_Flops &FlopCounter_in)
Set the internal Epetra_Flops() pointer.
int NumGlobalElements() const
Number of elements across all processors.
Epetra_Map: A class for partitioning vectors and matrices.
Definition: Epetra_Map.h:119
int Random()
Set multi-vector values to random numbers.
int Multiply(bool TransA, const Epetra_Vector &x, Epetra_Vector &y) const
Returns the result of a Epetra_CrsMatrix multiplied by a Epetra_Vector x in y.
double ElapsedTime(void) const
Epetra_Time elapsed time function.
void BiCGSTAB(Epetra_CrsMatrix &A, Epetra_Vector &x, Epetra_Vector &b, Ifpack_CrsRiluk *M, int Maxiter, double Tolerance, bool verbose)
Epetra_Export: This class builds an export object for efficient exporting of off-processor elements...
Definition: Epetra_Export.h:62
Epetra_Vector: A class for constructing and using dense vectors on a parallel computer.
int MyPID() const
Return my process ID.
Epetra_MpiComm: The Epetra MPI Communication Class.
int FillComplete(bool OptimizeDataStorage=true)
Signal that data entry is complete. Perform transformations to local index space. ...
Epetra_Time: The Epetra Timing Class.
Definition: Epetra_Time.h:75
const Epetra_Map & RowMap() const
Returns the Epetra_Map object associated with the rows of this matrix.
int applyInverse(Epetra_CrsMatrix &A, Epetra_Vector &x, Epetra_Vector &b, Ifpack_CrsRiluk *M, bool verbose)
int Export(const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
Exports an Epetra_DistObject using the Epetra_Import object.
int applyInverseDestroy(Ifpack_CrsRiluk *M)
Epetra_SerialComm: The Epetra Serial Communication Class.
Epetra_Flops: The Epetra Floating Point Operations Class.
Definition: Epetra_Flops.h:58
double Flops() const
Returns the number of floating point operations with this object and resets the count.
Definition: Epetra_Flops.h:74
int Update(double ScalarA, const Epetra_MultiVector &A, double ScalarThis)
Update multi-vector values with scaled values of A, this = ScalarThis*this + ScalarA*A.
Epetra_CrsMatrix: A class for constructing and using real-valued double-precision sparse compressed r...
int invIteration(Epetra_CrsMatrix &A, double &lambda, bool verbose)
int Norm2(double *Result) const
Compute 2-norm of each vector in multi-vector.
int main(int argc, char *argv[])
int applyInverseSetup(Epetra_CrsMatrix &A, Ifpack_CrsRiluk *&M)
const Epetra_CrsGraph & Graph() const
Returns a reference to the Epetra_CrsGraph object associated with this matrix.
Epetra_Flops * GetFlopCounter() const
Get the pointer to the Epetra_Flops() object associated with this object, returns 0 if none...
const Epetra_BlockMap & Map() const
Returns the address of the Epetra_BlockMap for this multi-vector.