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