download the original source code.
1 /*
2 Example 11
3
4 Interface: Linear-Algebraic (IJ)
5
6 Compile with: make ex11
7
8 Sample run: mpirun -np 4 ex11
9
10 Description: This example solves the 2-D Laplacian eigenvalue
11 problem with zero boundary conditions on an nxn grid.
12 The number of unknowns is N=n^2. The standard 5-point
13 stencil is used, and we solve for the interior nodes
14 only.
15
16 We use the same matrix as in Examples 3 and 5.
17 The eigensolver is LOBPCG with AMG preconditioner.
18 */
19
20 #include <math.h>
21 #include "_hypre_utilities.h"
22 #include "krylov.h"
23 #include "HYPRE.h"
24 #include "HYPRE_parcsr_ls.h"
25
26 /* lobpcg stuff */
27 #include "HYPRE_lobpcg.h"
28 #include "interpreter.h"
29 #include "HYPRE_MatvecFunctions.h"
30 #include "temp_multivector.h"
31 #include "_hypre_parcsr_mv.h"
32
33 #include "vis.c"
34
35 int main (int argc, char *argv[])
36 {
37 int i;
38 int myid, num_procs;
39 int N, n;
40 int blockSize;
41
42 int ilower, iupper;
43 int local_size, extra;
44
45 int vis;
46
47 double h, h2;
48
49 HYPRE_IJMatrix A;
50 HYPRE_ParCSRMatrix parcsr_A;
51 HYPRE_IJVector b;
52 HYPRE_ParVector par_b;
53 HYPRE_IJVector x;
54 HYPRE_ParVector par_x;
55 HYPRE_ParVector* pvx;
56
57 HYPRE_Solver precond, lobpcg_solver;
58 mv_InterfaceInterpreter* interpreter;
59 HYPRE_MatvecFunctions matvec_fn;
60
61 /* Initialize MPI */
62 MPI_Init(&argc, &argv);
63 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
64 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
65
66 /* Default problem parameters */
67 n = 33;
68 blockSize = 10;
69 vis = 0;
70
71 /* Parse command line */
72 {
73 int arg_index = 0;
74 int print_usage = 0;
75
76 while (arg_index < argc)
77 {
78 if ( strcmp(argv[arg_index], "-n") == 0 )
79 {
80 arg_index++;
81 n = atoi(argv[arg_index++]);
82 }
83 else if ( strcmp(argv[arg_index], "-blockSize") == 0 )
84 {
85 arg_index++;
86 blockSize = atoi(argv[arg_index++]);
87 }
88 else if ( strcmp(argv[arg_index], "-vis") == 0 )
89 {
90 arg_index++;
91 vis = 1;
92 }
93 else if ( strcmp(argv[arg_index], "-help") == 0 )
94 {
95 print_usage = 1;
96 break;
97 }
98 else
99 {
100 arg_index++;
101 }
102 }
103
104 if ((print_usage) && (myid == 0))
105 {
106 printf("\n");
107 printf("Usage: %s [<options>]\n", argv[0]);
108 printf("\n");
109 printf(" -n <n> : problem size in each direction (default: 33)\n");
110 printf(" -blockSize <n> : eigenproblem block size (default: 10)\n");
111 printf(" -vis : save the solution for GLVis visualization\n");
112 printf("\n");
113 }
114
115 if (print_usage)
116 {
117 MPI_Finalize();
118 return (0);
119 }
120 }
121
122 /* Preliminaries: want at least one processor per row */
123 if (n*n < num_procs) n = sqrt(num_procs) + 1;
124 N = n*n; /* global number of rows */
125 h = 1.0/(n+1); /* mesh size*/
126 h2 = h*h;
127
128 /* Each processor knows only of its own rows - the range is denoted by ilower
129 and iupper. Here we partition the rows. We account for the fact that
130 N may not divide evenly by the number of processors. */
131 local_size = N/num_procs;
132 extra = N - local_size*num_procs;
133
134 ilower = local_size*myid;
135 ilower += hypre_min(myid, extra);
136
137 iupper = local_size*(myid+1);
138 iupper += hypre_min(myid+1, extra);
139 iupper = iupper - 1;
140
141 /* How many rows do I have? */
142 local_size = iupper - ilower + 1;
143
144 /* Create the matrix.
145 Note that this is a square matrix, so we indicate the row partition
146 size twice (since number of rows = number of cols) */
147 HYPRE_IJMatrixCreate(MPI_COMM_WORLD, ilower, iupper, ilower, iupper, &A);
148
149 /* Choose a parallel csr format storage (see the User's Manual) */
150 HYPRE_IJMatrixSetObjectType(A, HYPRE_PARCSR);
151
152 /* Initialize before setting coefficients */
153 HYPRE_IJMatrixInitialize(A);
154
155 /* Now go through my local rows and set the matrix entries.
156 Each row has at most 5 entries. For example, if n=3:
157
158 A = [M -I 0; -I M -I; 0 -I M]
159 M = [4 -1 0; -1 4 -1; 0 -1 4]
160
161 Note that here we are setting one row at a time, though
162 one could set all the rows together (see the User's Manual).
163 */
164 {
165 int nnz;
166 double values[5];
167 int cols[5];
168
169 for (i = ilower; i <= iupper; i++)
170 {
171 nnz = 0;
172
173 /* The left identity block:position i-n */
174 if ((i-n)>=0)
175 {
176 cols[nnz] = i-n;
177 values[nnz] = -1.0;
178 nnz++;
179 }
180
181 /* The left -1: position i-1 */
182 if (i%n)
183 {
184 cols[nnz] = i-1;
185 values[nnz] = -1.0;
186 nnz++;
187 }
188
189 /* Set the diagonal: position i */
190 cols[nnz] = i;
191 values[nnz] = 4.0;
192 nnz++;
193
194 /* The right -1: position i+1 */
195 if ((i+1)%n)
196 {
197 cols[nnz] = i+1;
198 values[nnz] = -1.0;
199 nnz++;
200 }
201
202 /* The right identity block:position i+n */
203 if ((i+n)< N)
204 {
205 cols[nnz] = i+n;
206 values[nnz] = -1.0;
207 nnz++;
208 }
209
210 /* Set the values for row i */
211 HYPRE_IJMatrixSetValues(A, 1, &nnz, &i, cols, values);
212 }
213 }
214
215 /* Assemble after setting the coefficients */
216 HYPRE_IJMatrixAssemble(A);
217 /* Get the parcsr matrix object to use */
218 HYPRE_IJMatrixGetObject(A, (void**) &parcsr_A);
219
220 /* Create sample rhs and solution vectors */
221 HYPRE_IJVectorCreate(MPI_COMM_WORLD, ilower, iupper,&b);
222 HYPRE_IJVectorSetObjectType(b, HYPRE_PARCSR);
223 HYPRE_IJVectorInitialize(b);
224 HYPRE_IJVectorAssemble(b);
225 HYPRE_IJVectorGetObject(b, (void **) &par_b);
226
227 HYPRE_IJVectorCreate(MPI_COMM_WORLD, ilower, iupper,&x);
228 HYPRE_IJVectorSetObjectType(x, HYPRE_PARCSR);
229 HYPRE_IJVectorInitialize(x);
230 HYPRE_IJVectorAssemble(x);
231 HYPRE_IJVectorGetObject(x, (void **) &par_x);
232
233 /* Create a preconditioner and solve the eigenproblem */
234
235 /* AMG preconditioner */
236 {
237 HYPRE_BoomerAMGCreate(&precond);
238 HYPRE_BoomerAMGSetPrintLevel(precond, 1); /* print amg solution info */
239 HYPRE_BoomerAMGSetCoarsenType(precond, 6);
240 HYPRE_BoomerAMGSetRelaxType(precond, 6); /* Sym G.S./Jacobi hybrid */
241 HYPRE_BoomerAMGSetNumSweeps(precond, 1);
242 HYPRE_BoomerAMGSetTol(precond, 0.0); /* conv. tolerance zero */
243 HYPRE_BoomerAMGSetMaxIter(precond, 1); /* do only one iteration! */
244 }
245
246 /* LOBPCG eigensolver */
247 {
248 int time_index;
249
250 int maxIterations = 100; /* maximum number of iterations */
251 int pcgMode = 1; /* use rhs as initial guess for inner pcg iterations */
252 int verbosity = 1; /* print iterations info */
253 double tol = 1.e-8; /* absolute tolerance (all eigenvalues) */
254 int lobpcgSeed = 775; /* random seed */
255
256 mv_MultiVectorPtr eigenvectors = NULL;
257 mv_MultiVectorPtr constraints = NULL;
258 double *eigenvalues = NULL;
259
260 if (myid != 0)
261 verbosity = 0;
262
263 /* define an interpreter for the ParCSR interface */
264 interpreter = hypre_CTAlloc(mv_InterfaceInterpreter,1);
265 HYPRE_ParCSRSetupInterpreter(interpreter);
266 HYPRE_ParCSRSetupMatvec(&matvec_fn);
267
268 /* eigenvectors - create a multivector */
269 eigenvectors =
270 mv_MultiVectorCreateFromSampleVector(interpreter, blockSize, par_x);
271 mv_MultiVectorSetRandom (eigenvectors, lobpcgSeed);
272
273 /* eigenvectors - get a pointer */
274 {
275 mv_TempMultiVector* tmp = mv_MultiVectorGetData(eigenvectors);
276 pvx = (HYPRE_ParVector*)(tmp -> vector);
277 }
278
279 /* eigenvalues - allocate space */
280 eigenvalues = (double*) calloc( blockSize, sizeof(double) );
281
282 HYPRE_LOBPCGCreate(interpreter, &matvec_fn, &lobpcg_solver);
283 HYPRE_LOBPCGSetMaxIter(lobpcg_solver, maxIterations);
284 HYPRE_LOBPCGSetPrecondUsageMode(lobpcg_solver, pcgMode);
285 HYPRE_LOBPCGSetTol(lobpcg_solver, tol);
286 HYPRE_LOBPCGSetPrintLevel(lobpcg_solver, verbosity);
287
288 /* use a preconditioner */
289 HYPRE_LOBPCGSetPrecond(lobpcg_solver,
290 (HYPRE_PtrToSolverFcn) HYPRE_BoomerAMGSolve,
291 (HYPRE_PtrToSolverFcn) HYPRE_BoomerAMGSetup,
292 precond);
293
294 HYPRE_LOBPCGSetup(lobpcg_solver, (HYPRE_Matrix)parcsr_A,
295 (HYPRE_Vector)par_b, (HYPRE_Vector)par_x);
296
297 time_index = hypre_InitializeTiming("LOBPCG Solve");
298 hypre_BeginTiming(time_index);
299
300 HYPRE_LOBPCGSolve(lobpcg_solver, constraints, eigenvectors, eigenvalues );
301
302 hypre_EndTiming(time_index);
303 hypre_PrintTiming("Solve phase times", MPI_COMM_WORLD);
304 hypre_FinalizeTiming(time_index);
305 hypre_ClearTiming();
306
307 /* clean-up */
308 HYPRE_BoomerAMGDestroy(precond);
309 HYPRE_LOBPCGDestroy(lobpcg_solver);
310 hypre_TFree(eigenvalues);
311 hypre_TFree(interpreter);
312 }
313
314 /* Save the solution for GLVis visualization, see vis/glvis-ex11.sh */
315 if (vis)
316 {
317 FILE *file;
318 char filename[255];
319
320 int nvalues = local_size;
321 double *values;
322
323 /* get the local solution */
324 values = hypre_VectorData(hypre_ParVectorLocalVector(
325 (hypre_ParVector*)pvx[blockSize-1]));
326
327 sprintf(filename, "%s.%06d", "vis/ex11.sol", myid);
328 if ((file = fopen(filename, "w")) == NULL)
329 {
330 printf("Error: can't open output file %s\n", filename);
331 MPI_Finalize();
332 exit(1);
333 }
334
335 /* save solution */
336 for (i = 0; i < nvalues; i++)
337 fprintf(file, "%.14e\n", values[i]);
338
339 fflush(file);
340 fclose(file);
341
342 /* save global finite element mesh */
343 if (myid == 0)
344 GLVis_PrintGlobalSquareMesh("vis/ex11.mesh", n-1);
345 }
346
347 /* Clean up */
348 HYPRE_IJMatrixDestroy(A);
349 HYPRE_IJVectorDestroy(b);
350 HYPRE_IJVectorDestroy(x);
351
352 /* Finalize MPI*/
353 MPI_Finalize();
354
355 return(0);
356 }
syntax highlighted by Code2HTML, v. 0.9.1