download the original source code.
1 /*
2 Example 10
3
4 Interface: Finite Element Interface (FEI)
5
6 Compile with: make ex10
7
8 Sample run: mpirun -np 4 ex10 -n 120 -solver 2
9
10 To see options: ex10 -help
11
12 Description: This code solves a system corresponding to a discretization
13 of the Laplace equation -Delta u = 1 with zero boundary
14 conditions on the unit square. The domain is split into
15 a n x n grid of quadrilateral elements and each processors
16 owns a horizontal strip of size m x n, where m = n/nprocs. We
17 use bilinear finite element discretization, so there are
18 nodes (vertices) that are shared between neighboring
19 processors. The Finite Element Interface is used to assemble
20 the matrix and solve the problem. Nine different solvers are
21 available.
22 */
23
24 #include <math.h>
25 #include <iostream.h>
26 #include <fstream.h>
27 #include "_hypre_utilities.h"
28 #include "LLNL_FEI_Impl.h"
29
30 int main(int argc, char *argv[])
31 {
32 int i, j, k;
33
34 int nprocs, mypid;
35
36 int n, m, offset;
37 double h;
38
39 int solverID;
40 int print_solution;
41
42 // Initialize MPI
43 MPI_Init(&argc, &argv);
44 MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
45 MPI_Comm_rank(MPI_COMM_WORLD, &mypid);
46
47 // Set default parameters
48 n = 4*nprocs;
49 solverID = 2;
50 print_solution = 0;
51
52 // Parse command line
53 {
54 int arg_index = 0;
55 int print_usage = 0;
56
57 while (arg_index < argc)
58 {
59 if ( strcmp(argv[arg_index], "-n") == 0 )
60 {
61 arg_index++;
62 n = atoi(argv[arg_index++]);
63 }
64 else if ( strcmp(argv[arg_index], "-solver") == 0 )
65 {
66 arg_index++;
67 solverID = atoi(argv[arg_index++]);
68 }
69 else if ( strcmp(argv[arg_index], "-print_solution") == 0 )
70 {
71 arg_index++;
72 print_solution = 1;
73 }
74 else if ( strcmp(argv[arg_index], "-help") == 0 )
75 {
76 print_usage = 1;
77 break;
78 }
79 else
80 {
81 arg_index++;
82 }
83 }
84
85 if ((print_usage) && (mypid == 0))
86 {
87 printf("\n");
88 printf("Usage: %s [<options>]\n", argv[0]);
89 printf("\n");
90 printf(" -n <n> : problem size per processor (default: %d)\n", 4*nprocs);
91 printf(" -solver <ID> : solver ID\n");
92 printf(" 0 - DS-PCG\n");
93 printf(" 1 - ParaSails-PCG\n");
94 printf(" 2 - AMG-PCG (default)\n");
95 printf(" 3 - AMGSA-PCG\n");
96 printf(" 4 - Euclid-PCG\n");
97 printf(" 5 - DS-GMRES\n");
98 printf(" 6 - AMG-GMRES\n");
99 printf(" 7 - AMGSA-GMRES\n");
100 printf(" 8 - Euclid-GMRES\n");
101 printf(" -print_solution : print the solution vector\n");
102 printf("\n");
103 }
104
105 if (print_usage)
106 {
107 MPI_Finalize();
108 return (0);
109 }
110 }
111
112 // Each processor owns a m x n grid of quadrilateral finite elements.
113 // The unknowns are located in the nodes (vertices of the mesh) and
114 // are numbered globally starting from the lower left corner and moving
115 // row-wise to the upper right corner.
116 m = n / nprocs;
117 offset = mypid*(m*(n+1));
118
119 h = 1.0 / n; // mesh size
120
121 // 1. FEI initialization phase
122
123 // Instantiate the FEI object
124 LLNL_FEI_Impl *feiPtr = new LLNL_FEI_Impl(MPI_COMM_WORLD);
125
126 // Set the matrix storage type to HYPRE
127 {
128 char **paramStrings = new char*[1];
129 paramStrings[0] = new char[100];
130 strcpy(paramStrings[0], "externalSolver HYPRE");
131 feiPtr->parameters(1, paramStrings);
132 delete paramStrings[0];
133 delete [] paramStrings;
134 }
135
136 // The unknowns in FEI are called fields. Each field has an
137 // identifier (fieldID) and rank (fieldSize).
138 int nFields = 1;
139 int *fieldSizes = new int[nFields]; fieldSizes[0] = 1;
140 int *fieldIDs = new int[nFields]; fieldIDs[0] = 0;
141
142 // Pass the field information to the FEI
143 feiPtr->initFields(nFields, fieldSizes, fieldIDs);
144
145 // Elements are grouped into blocks (in this case one block), and we
146 // have to describe the number of elements in the block (nElems) as
147 // well as the fields (unknowns) per element.
148 int elemBlkID = 0;
149 int nElems = m*n;
150 int elemNNodes = 4; // number of (shared) nodes per element
151 int *nodeNFields = new int[elemNNodes]; // fields per node
152 int **nodeFieldIDs = new int*[elemNNodes]; // node-fields IDs
153 int elemNFields = 0; // number of (non-shared) fields per element
154 int *elemFieldIDs = NULL; // element-fields IDs
155 for (i = 0; i < elemNNodes; i++)
156 {
157 nodeNFields[i] = 1;
158 nodeFieldIDs[i] = new int[nodeNFields[i]];
159 nodeFieldIDs[i][0] = fieldIDs[0];
160 }
161
162 // Pass the block information to the FEI. The interleave parameter
163 // controls how different fields are ordered in the element matrices.
164 int interleave = 0;
165 feiPtr->initElemBlock(elemBlkID, nElems, elemNNodes, nodeNFields,
166 nodeFieldIDs, elemNFields, elemFieldIDs, interleave);
167
168 // List the global indexes (IDs) of the nodes in each element
169 int **elemConn = new int*[nElems];
170 for (i = 0; i < m; i++)
171 for (j = 0; j < n; j++)
172 {
173 elemConn[i*n+j] = new int[elemNNodes]; // element with coordinates (i,j)
174 elemConn[i*n+j][0] = offset + i*(n+1)+j; // node in the lower left
175 elemConn[i*n+j][1] = elemConn[i*n+j][0]+1; // node in the lower right
176 elemConn[i*n+j][2] = elemConn[i*n+j][1]+n+1; // node in the upper right
177 elemConn[i*n+j][3] = elemConn[i*n+j][2]-1; // node in the upper left
178 }
179
180 // Pass the element topology information to the FEI
181 for (i = 0; i < nElems; i++)
182 feiPtr->initElem(elemBlkID, i, elemConn[i]);
183
184 // List the global indexes of nodes that are shared between processors
185 int nShared, *SharedIDs, *SharedLengs, **SharedProcs;
186 if (mypid == 0)
187 {
188 // Nodes in the top row are shared
189 nShared = n+1;
190 SharedIDs = new int[nShared];
191 for (i = 0; i < nShared; i++)
192 SharedIDs[i] = offset + m*(n+1) + i;
193 SharedLengs = new int[nShared];
194 for (i = 0; i < nShared; i++)
195 SharedLengs[i] = 2;
196 SharedProcs = new int*[nShared];
197 for (i = 0; i < nShared; i++)
198 {
199 SharedProcs[i] = new int[SharedLengs[i]];
200 SharedProcs[i][0] = mypid;
201 SharedProcs[i][1] = mypid+1;
202 }
203 }
204 else if (mypid == nprocs-1)
205 {
206 // Nodes in the bottom row are shared
207 nShared = n+1;
208 SharedIDs = new int[nShared];
209 for (i = 0; i < nShared; i++)
210 SharedIDs[i] = offset + i;
211 SharedLengs = new int[nShared];
212 for (i = 0; i < nShared; i++)
213 SharedLengs[i] = 2;
214 SharedProcs = new int*[nShared];
215 for (i = 0; i < nShared; i++)
216 {
217 SharedProcs[i] = new int[SharedLengs[i]];
218 SharedProcs[i][0] = mypid-1;
219 SharedProcs[i][1] = mypid;
220 }
221 }
222 else
223 {
224 // Nodes in the top and bottom rows are shared
225 nShared = 2*(n+1);
226 SharedIDs = new int[nShared];
227 for (i = 0; i < n+1; i++)
228 {
229 SharedIDs[i] = offset + i;
230 SharedIDs[n+1+i] = offset + m*(n+1) + i;
231 }
232 SharedLengs = new int[nShared];
233 for (i = 0; i < nShared; i++)
234 SharedLengs[i] = 2;
235 SharedProcs = new int*[nShared];
236 for (i = 0; i < n+1; i++)
237 {
238 SharedProcs[i] = new int[SharedLengs[i]];
239 SharedProcs[i][0] = mypid-1;
240 SharedProcs[i][1] = mypid;
241
242 SharedProcs[n+1+i] = new int[SharedLengs[n+1+i]];
243 SharedProcs[n+1+i][0] = mypid;
244 SharedProcs[n+1+i][1] = mypid+1;
245 }
246 }
247
248 // Pass the shared nodes information to the FEI
249 if (nprocs != 1 && nShared > 0)
250 feiPtr->initSharedNodes(nShared, SharedIDs, SharedLengs, SharedProcs);
251
252 // Finish the FEI initialization phase
253 feiPtr->initComplete();
254
255 // 2. FEI load phase
256
257 // Specify the boundary conditions
258 int nBCs, *BCEqn;
259 double **alpha, **beta, **gamma;
260 if (mypid == 0)
261 {
262 // Nodes in the bottom row and left and right columns
263 nBCs = n+1 + 2*m;
264 BCEqn = new int[nBCs];
265 for (i = 0; i < n+1; i++)
266 BCEqn[i] = offset + i;
267 for (i = 0; i < m; i++)
268 {
269 BCEqn[n+1+2*i] = offset + (i+1)*(n+1);
270 BCEqn[n+2+2*i] = offset + (i+1)*(n+1)+n;
271 }
272 }
273 else if (mypid == nprocs-1)
274 {
275 // Nodes in the top row and left and right columns
276 nBCs = n+1 + 2*m;
277 BCEqn = new int[nBCs];
278 for (i = 0; i < n+1; i++)
279 BCEqn[i] = offset + m*(n+1) + i;
280 for (i = 0; i < m; i++)
281 {
282 BCEqn[n+1+2*i] = offset + i*(n+1);
283 BCEqn[n+2+2*i] = offset + i*(n+1)+n;
284 }
285 }
286 else
287 {
288 // Nodes in the left and right columns
289 nBCs = 2*(m+1);
290 BCEqn = new int[nBCs];
291 for (i = 0; i < m+1; i++)
292 {
293 BCEqn[2*i] = offset + i*(n+1);
294 BCEqn[2*i+1] = offset + i*(n+1)+n;
295 }
296 }
297
298 // The arrays alpha, beta and gamma specify the type of boundary
299 // condition (essential, natural, mixed). The most general form
300 // for Laplace problems is alpha U + beta dU/dn = gamma. In this
301 // example we impose zero Dirichlet boundary conditions.
302 alpha = new double*[nBCs];
303 beta = new double*[nBCs];
304 gamma = new double*[nBCs];
305 for (i = 0; i < nBCs; i++)
306 {
307 alpha[i] = new double[1]; alpha[i][0] = 1.0;
308 beta[i] = new double[1]; beta[i][0] = 0.0;
309 gamma[i] = new double[1]; gamma[i][0] = 0.0;
310 }
311
312 // Pass the boundary condition information to the FEI
313 feiPtr->loadNodeBCs(nBCs, BCEqn, fieldIDs[0], alpha, beta, gamma);
314
315 // Specify element stiffness matrices
316 double ***elemStiff = new double**[nElems];
317 for (i = 0; i < m; i++)
318 for (j = 0; j < n; j++)
319 {
320 // Element with coordinates (i,j)
321 elemStiff[i*n+j] = new double*[elemNNodes];
322 for (k = 0; k < elemNNodes; k++)
323 elemStiff[i*n+j][k] = new double[elemNNodes];
324
325 // Stiffness matrix for the reference square
326 // 3 +---+ 2
327 // | |
328 // 0 +---+ 1
329
330 double **A = elemStiff[i*n+j];
331
332 for (k = 0; k < 4; k++)
333 A[k][k] = 2/3.;
334
335 A[0][1] = A[1][0] = -1/6.;
336 A[0][2] = A[2][0] = -1/3.;
337 A[0][3] = A[3][0] = -1/6.;
338 A[1][2] = A[2][1] = -1/6.;
339 A[1][3] = A[3][1] = -1/3.;
340 A[2][3] = A[3][2] = -1/6.;
341 }
342
343 // Specify element load vectors
344 double *elemLoad = new double[nElems*elemNNodes];
345 for (i = 0; i < nElems*elemNNodes; i++)
346 elemLoad[i] = h*h/4;
347
348 // Assemble the matrix. The elemFormat parameter describes
349 // the storage (symmetric/non-symmetric, row/column-wise)
350 // of the element stiffness matrices.
351 int elemFormat = 0;
352 for (i = 0; i < nElems; i++)
353 feiPtr->sumInElem(elemBlkID, i, elemConn[i], elemStiff[i],
354 &(elemLoad[i*elemNNodes]), elemFormat);
355
356 // Finish the FEI load phase
357 feiPtr->loadComplete();
358
359 // Clean up
360 for (i = 0; i < nElems; i++) delete [] elemConn[i];
361 delete [] elemConn;
362 for (i = 0; i < nElems; i++)
363 {
364 for (j = 0; j < elemNNodes; j++) delete [] elemStiff[i][j];
365 delete [] elemStiff[i];
366 }
367 delete [] elemStiff;
368 delete [] elemLoad;
369
370 delete [] BCEqn;
371 for (i = 0; i < nBCs; i++)
372 {
373 delete [] alpha[i];
374 delete [] beta[i];
375 delete [] gamma[i];
376 }
377 delete [] alpha;
378 delete [] beta;
379 delete [] gamma;
380
381 if (nShared > 0)
382 {
383 delete [] SharedIDs;
384 delete [] SharedLengs;
385 for (i = 0; i < nShared; i++) delete [] SharedProcs[i];
386 delete [] SharedProcs;
387 }
388
389 delete [] nodeNFields;
390 for (i = 0; i < elemNNodes; i++) delete [] nodeFieldIDs[i];
391 delete [] nodeFieldIDs;
392
393 delete [] fieldSizes;
394 delete [] fieldIDs;
395
396 // 3. Set up problem parameters and pass them to the FEI
397 {
398 int nParams = 19;
399 char **paramStrings = new char*[nParams];
400 for (i = 0; i < nParams; i++)
401 paramStrings[i] = new char[100];
402
403 strcpy(paramStrings[0], "outputLevel 2");
404 switch(solverID)
405 {
406 case 0:
407 strcpy(paramStrings[1], "solver cg");
408 strcpy(paramStrings[2], "preconditioner diagonal");
409 break;
410 case 1:
411 strcpy(paramStrings[1], "solver cg");
412 strcpy(paramStrings[2], "preconditioner parasails");
413 break;
414 default:
415 case 2:
416 strcpy(paramStrings[1], "solver cg");
417 strcpy(paramStrings[2], "preconditioner boomeramg");
418 break;
419 case 3:
420 strcpy(paramStrings[1], "solver cg");
421 strcpy(paramStrings[2], "preconditioner mli");
422 break;
423 case 4:
424 strcpy(paramStrings[1], "solver cg");
425 strcpy(paramStrings[2], "preconditioner euclid");
426 break;
427 case 5:
428 strcpy(paramStrings[1], "solver gmres");
429 strcpy(paramStrings[2], "preconditioner diagonal");
430 break;
431 case 6:
432 strcpy(paramStrings[1], "solver gmres");
433 strcpy(paramStrings[2], "preconditioner boomeramg");
434 break;
435 case 7:
436 strcpy(paramStrings[1], "solver gmres");
437 strcpy(paramStrings[2], "preconditioner mli");
438 break;
439 case 8:
440 strcpy(paramStrings[1], "solver gmres");
441 strcpy(paramStrings[2], "preconditioner euclid");
442 break;
443 }
444 strcpy(paramStrings[3], "maxIterations 100");
445 strcpy(paramStrings[4], "tolerance 1e-6");
446 strcpy(paramStrings[5], "gmresDim 30");
447 strcpy(paramStrings[6], "amgNumSweeps 1");
448 strcpy(paramStrings[7], "amgCoarsenType falgout");
449 strcpy(paramStrings[8], "amgRelaxType hybridsym");
450 strcpy(paramStrings[9], "amgSystemSize 1");
451 strcpy(paramStrings[10], "amgStrongThreshold 0.25");
452 strcpy(paramStrings[11], "MLI smoother HSGS");
453 strcpy(paramStrings[12], "MLI numSweeps 1");
454 strcpy(paramStrings[13], "MLI smootherWeight 1.0");
455 strcpy(paramStrings[14], "MLI nodeDOF 1");
456 strcpy(paramStrings[15], "MLI nullSpaceDim 1");
457 strcpy(paramStrings[16], "MLI minCoarseSize 50");
458 strcpy(paramStrings[17], "MLI outputLevel 0");
459 strcpy(paramStrings[18], "parasailsSymmetric outputLevel 0");
460
461 feiPtr->parameters(nParams, paramStrings);
462
463 for (i = 0; i < nParams; i++)
464 delete [] paramStrings[i];
465 delete [] paramStrings;
466 }
467
468 // 4. Solve the system
469 int status;
470 feiPtr->solve(&status);
471
472 // 5. Save the solution and destroy the FEI object
473 if (print_solution)
474 {
475 int numNodes, *nodeIDList, *solnOffsets;
476 double *solnValues;
477
478 // Get the number of nodes in the element block
479 feiPtr->getNumBlockActNodes(elemBlkID, &numNodes);
480
481 // Get their global IDs
482 nodeIDList = new int[numNodes];
483 feiPtr->getBlockNodeIDList(elemBlkID, numNodes, nodeIDList);
484
485 // Get the values corresponding to nodeIDList
486 solnOffsets = new int[numNodes];
487 solnValues = new double[numNodes];
488 feiPtr->getBlockNodeSolution(elemBlkID, numNodes, nodeIDList,
489 solnOffsets, solnValues);
490
491 // Find the location of the ith local node
492 for (i = 0; i < numNodes; i++)
493 solnOffsets[nodeIDList[i]-offset] = i;
494
495 // Save the ordered nodal values to a file
496 char sol_out[20];
497 sprintf(sol_out, "fei-test.sol_%d", mypid);
498 ofstream sol(sol_out);
499 sol << sol_out << endl;
500 for (i = 0; i < numNodes; i++)
501 sol << solnValues[solnOffsets[i]] << endl;
502
503 // Clean up
504 delete [] solnValues;
505 delete [] solnOffsets;
506 delete [] nodeIDList;
507 }
508 delete feiPtr;
509
510 // Finalize MPI
511 MPI_Finalize();
512
513 return (0);
514 }
syntax highlighted by Code2HTML, v. 0.9.1