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