Actual source code: ex55.c
1: static char help[] = "2D, bi-linear quadrilateral (Q1), displacement finite element formulation\n\
2: of plain strain linear elasticity. E=1.0, nu=0.25.\n\
3: Unit square domain with Dirichelet boundary condition on the y=0 side only.\n\
4: Load of 1.0 in x direction on all nodes (not a true uniform load).\n\
5: -ne <size> : number of (square) quadrilateral elements in each dimension\n\
6: -alpha <v> : scaling of material coeficient in embedded circle\n\n";
8: #include <petscksp.h>
10: int main(int argc,char **args)
11: {
12: Mat Amat;
14: PetscInt i,m,M,its,Istart,Iend,j,Ii,ix,ne=4;
15: PetscReal x,y,h;
16: Vec xx,bb;
17: KSP ksp;
18: PetscReal soft_alpha = 1.e-3;
19: MPI_Comm comm;
20: PetscBool use_coords = PETSC_FALSE;
21: PetscMPIInt npe,mype;
22: PetscScalar DD[8][8],DD2[8][8];
23: #if defined(PETSC_USE_LOG)
24: PetscLogStage stage[2];
25: #endif
26: PetscScalar DD1[8][8] = { {5.333333333333333E-01, 2.0000E-01, -3.333333333333333E-01, 0.0000E+00, -2.666666666666667E-01, -2.0000E-01, 6.666666666666667E-02, 0.0000E-00 },
27: {2.0000E-01, 5.333333333333333E-01, 0.0000E-00, 6.666666666666667E-02, -2.0000E-01, -2.666666666666667E-01, 0.0000E-00, -3.333333333333333E-01 },
28: {-3.333333333333333E-01, 0.0000E-00, 5.333333333333333E-01, -2.0000E-01, 6.666666666666667E-02, 0.0000E-00, -2.666666666666667E-01, 2.0000E-01 },
29: {0.0000E+00, 6.666666666666667E-02, -2.0000E-01, 5.333333333333333E-01, 0.0000E-00, -3.333333333333333E-01, 2.0000E-01, -2.666666666666667E-01 },
30: {-2.666666666666667E-01, -2.0000E-01, 6.666666666666667E-02, 0.0000E-00, 5.333333333333333E-01, 2.0000E-01, -3.333333333333333E-01, 0.0000E+00 },
31: {-2.0000E-01, -2.666666666666667E-01, 0.0000E-00, -3.333333333333333E-01, 2.0000E-01, 5.333333333333333E-01, 0.0000E-00, 6.666666666666667E-02 },
32: {6.666666666666667E-02, 0.0000E-00, -2.666666666666667E-01, 2.0000E-01, -3.333333333333333E-01, 0.0000E-00, 5.333333333333333E-01, -2.0000E-01 },
33: {0.0000E-00, -3.333333333333333E-01, 2.0000E-01, -2.666666666666667E-01, 0.0000E-00, 6.666666666666667E-02, -2.0000E-01, 5.333333333333333E-01 } };
35: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
36: comm = PETSC_COMM_WORLD;
37: MPI_Comm_rank(comm, &mype);
38: MPI_Comm_size(comm, &npe);
39: PetscOptionsGetInt(NULL,NULL,"-ne",&ne,NULL);
40: h = 1./ne;
41: /* ne*ne; number of global elements */
42: PetscOptionsGetReal(NULL,NULL,"-alpha",&soft_alpha,NULL);
43: PetscOptionsGetBool(NULL,NULL,"-use_coordinates",&use_coords,NULL);
44: M = 2*(ne+1)*(ne+1); /* global number of equations */
45: m = (ne+1)*(ne+1)/npe;
46: if (mype==npe-1) m = (ne+1)*(ne+1) - (npe-1)*m;
47: m *= 2;
48: /* create stiffness matrix */
49: MatCreate(comm,&Amat);
50: MatSetSizes(Amat,m,m,M,M);
51: MatSetType(Amat,MATAIJ);
52: MatSetOption(Amat,MAT_SPD,PETSC_TRUE);
53: MatSetFromOptions(Amat);
54: MatSetBlockSize(Amat,2);
55: MatSeqAIJSetPreallocation(Amat,18,NULL);
56: MatMPIAIJSetPreallocation(Amat,18,NULL,18,NULL);
57: #if defined(PETSC_HAVE_HYPRE)
58: MatHYPRESetPreallocation(Amat,18,NULL,18,NULL);
59: #endif
61: MatGetOwnershipRange(Amat,&Istart,&Iend);
62: if (m != Iend - Istart) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"m %D does not equal Iend %D - Istart %D",m,Iend,Istart);
63: /* Generate vectors */
64: MatCreateVecs(Amat,&xx,&bb);
65: VecSet(bb,.0);
66: /* generate element matrices -- see ex56.c on how to use different data set */
67: {
68: DD[0][0] = 0.53333333333333321;
69: DD[0][1] = 0.20000000000000001;
70: DD[0][2] = -0.33333333333333331;
71: DD[0][3] = 0.0000000000000000;
72: DD[0][4] = -0.26666666666666666;
73: DD[0][5] = -0.20000000000000001;
74: DD[0][6] = 6.66666666666666796E-002;
75: DD[0][7] = 6.93889390390722838E-018;
76: DD[1][0] = 0.20000000000000001;
77: DD[1][1] = 0.53333333333333333;
78: DD[1][2] = 7.80625564189563192E-018;
79: DD[1][3] = 6.66666666666666935E-002;
80: DD[1][4] = -0.20000000000000001;
81: DD[1][5] = -0.26666666666666666;
82: DD[1][6] = -3.46944695195361419E-018;
83: DD[1][7] = -0.33333333333333331;
84: DD[2][0] = -0.33333333333333331;
85: DD[2][1] = 1.12757025938492461E-017;
86: DD[2][2] = 0.53333333333333333;
87: DD[2][3] = -0.20000000000000001;
88: DD[2][4] = 6.66666666666666935E-002;
89: DD[2][5] = -6.93889390390722838E-018;
90: DD[2][6] = -0.26666666666666666;
91: DD[2][7] = 0.19999999999999998;
92: DD[3][0] = 0.0000000000000000;
93: DD[3][1] = 6.66666666666666935E-002;
94: DD[3][2] = -0.20000000000000001;
95: DD[3][3] = 0.53333333333333333;
96: DD[3][4] = 4.33680868994201774E-018;
97: DD[3][5] = -0.33333333333333331;
98: DD[3][6] = 0.20000000000000001;
99: DD[3][7] = -0.26666666666666666;
100: DD[4][0] = -0.26666666666666666;
101: DD[4][1] = -0.20000000000000001;
102: DD[4][2] = 6.66666666666666935E-002;
103: DD[4][3] = 8.67361737988403547E-019;
104: DD[4][4] = 0.53333333333333333;
105: DD[4][5] = 0.19999999999999998;
106: DD[4][6] = -0.33333333333333331;
107: DD[4][7] = -3.46944695195361419E-018;
108: DD[5][0] = -0.20000000000000001;
109: DD[5][1] = -0.26666666666666666;
110: DD[5][2] = -1.04083408558608426E-017;
111: DD[5][3] = -0.33333333333333331;
112: DD[5][4] = 0.19999999999999998;
113: DD[5][5] = 0.53333333333333333;
114: DD[5][6] = 6.93889390390722838E-018;
115: DD[5][7] = 6.66666666666666519E-002;
116: DD[6][0] = 6.66666666666666796E-002;
117: DD[6][1] = -6.93889390390722838E-018;
118: DD[6][2] = -0.26666666666666666;
119: DD[6][3] = 0.19999999999999998;
120: DD[6][4] = -0.33333333333333331;
121: DD[6][5] = 6.93889390390722838E-018;
122: DD[6][6] = 0.53333333333333321;
123: DD[6][7] = -0.20000000000000001;
124: DD[7][0] = 6.93889390390722838E-018;
125: DD[7][1] = -0.33333333333333331;
126: DD[7][2] = 0.19999999999999998;
127: DD[7][3] = -0.26666666666666666;
128: DD[7][4] = 0.0000000000000000;
129: DD[7][5] = 6.66666666666666519E-002;
130: DD[7][6] = -0.20000000000000001;
131: DD[7][7] = 0.53333333333333321;
133: /* BC version of element */
134: for (i=0; i<8; i++) {
135: for (j=0; j<8; j++) {
136: if (i<4 || j < 4) {
137: if (i==j) DD2[i][j] = .1*DD1[i][j];
138: else DD2[i][j] = 0.0;
139: } else DD2[i][j] = DD1[i][j];
140: }
141: }
142: }
143: {
144: PetscReal *coords;
145: PetscMalloc1(m,&coords);
146: /* forms the element stiffness and coordinates */
147: for (Ii = Istart/2, ix = 0; Ii < Iend/2; Ii++, ix++) {
148: j = Ii/(ne+1); i = Ii%(ne+1);
149: /* coords */
150: x = h*(Ii % (ne+1)); y = h*(Ii/(ne+1));
151: coords[2*ix] = x; coords[2*ix+1] = y;
152: if (i<ne && j<ne) {
153: PetscInt jj,ii,idx[4];
154: /* radius */
155: PetscReal radius = PetscSqrtReal((x-.5+h/2)*(x-.5+h/2) + (y-.5+h/2)*(y-.5+h/2));
156: PetscReal alpha = 1.0;
157: if (radius < 0.25) alpha = soft_alpha;
159: idx[0] = Ii; idx[1] = Ii+1; idx[2] = Ii + (ne+1) + 1; idx[3] = Ii + (ne+1);
160: for (ii=0; ii<8; ii++) {
161: for (jj=0;jj<8;jj++) DD[ii][jj] = alpha*DD1[ii][jj];
162: }
163: if (j>0) {
164: MatSetValuesBlocked(Amat,4,idx,4,idx,(const PetscScalar*)DD,ADD_VALUES);
165: } else {
166: /* a BC */
167: for (ii=0; ii<8; ii++) {
168: for (jj=0;jj<8;jj++) DD[ii][jj] = alpha*DD2[ii][jj];
169: }
170: MatSetValuesBlocked(Amat,4,idx,4,idx,(const PetscScalar*)DD,ADD_VALUES);
171: }
172: }
173: if (j>0) {
174: PetscScalar v = h*h;
175: PetscInt jj = 2*Ii; /* load in x direction */
176: VecSetValues(bb,1,&jj,&v,INSERT_VALUES);
177: }
178: }
179: MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);
180: MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);
181: VecAssemblyBegin(bb);
182: VecAssemblyEnd(bb);
184: /* Setup solver */
185: KSPCreate(PETSC_COMM_WORLD,&ksp);
186: KSPSetFromOptions(ksp);
188: /* finish KSP/PC setup */
189: KSPSetOperators(ksp, Amat, Amat);
190: if (use_coords) {
191: PC pc;
193: KSPGetPC(ksp, &pc);
194: PCSetCoordinates(pc, 2, m/2, coords);
195: }
196: PetscFree(coords);
197: }
199: if (!PETSC_TRUE) {
200: PetscViewer viewer;
201: PetscViewerASCIIOpen(comm, "Amat.m", &viewer);
202: PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
203: MatView(Amat,viewer);
204: PetscViewerPopFormat(viewer);
205: PetscViewerDestroy(&viewer);
206: }
208: /* solve */
209: #if defined(PETSC_USE_LOG)
210: PetscLogStageRegister("Setup", &stage[0]);
211: PetscLogStageRegister("Solve", &stage[1]);
212: PetscLogStagePush(stage[0]);
213: #endif
214: KSPSetUp(ksp);
215: #if defined(PETSC_USE_LOG)
216: PetscLogStagePop();
217: #endif
219: VecSet(xx,.0);
221: #if defined(PETSC_USE_LOG)
222: PetscLogStagePush(stage[1]);
223: #endif
224: KSPSolve(ksp, bb, xx);
225: #if defined(PETSC_USE_LOG)
226: PetscLogStagePop();
227: #endif
229: KSPGetIterationNumber(ksp,&its);
231: if (0) {
232: PetscReal norm,norm2;
233: PetscViewer viewer;
234: Vec res;
236: PetscObjectGetComm((PetscObject)bb,&comm);
237: VecNorm(bb, NORM_2, &norm2);
239: VecDuplicate(xx, &res);
240: MatMult(Amat, xx, res);
241: VecAXPY(bb, -1.0, res);
242: VecDestroy(&res);
243: VecNorm(bb, NORM_2, &norm);
244: PetscPrintf(PETSC_COMM_WORLD,"[%d]%s |b-Ax|/|b|=%e, |b|=%e\n",0,PETSC_FUNCTION_NAME,norm/norm2,norm2);
245: PetscViewerASCIIOpen(comm, "residual.m", &viewer);
246: PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
247: VecView(bb,viewer);
248: PetscViewerPopFormat(viewer);
249: PetscViewerDestroy(&viewer);
250: }
252: /* Free work space */
253: KSPDestroy(&ksp);
254: VecDestroy(&xx);
255: VecDestroy(&bb);
256: MatDestroy(&Amat);
258: PetscFinalize();
259: return ierr;
260: }
262: /*TEST
264: test:
265: nsize: 4
266: args: -ne 29 -alpha 1.e-3 -ksp_type cg -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -use_coordinates -ksp_converged_reason -pc_gamg_esteig_ksp_max_it 5 -ksp_rtol 1.e-3 -ksp_monitor_short -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.2
267: output_file: output/ex55_sa.out
269: test:
270: suffix: Classical
271: nsize: 4
272: args: -ne 29 -alpha 1.e-3 -ksp_type cg -pc_type gamg -pc_gamg_type classical -mg_levels_ksp_max_it 5 -ksp_converged_reason
273: output_file: output/ex55_classical.out
275: test:
276: suffix: NC
277: nsize: 4
278: args: -ne 29 -alpha 1.e-3 -ksp_type cg -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -ksp_converged_reason -pc_gamg_esteig_ksp_max_it 10 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.2
280: test:
281: suffix: geo
282: nsize: 4
283: args: -ne 29 -alpha 1.e-3 -ksp_type cg -pc_type gamg -pc_gamg_type geo -use_coordinates -ksp_monitor_short -ksp_type cg -ksp_norm_type unpreconditioned -mg_levels_ksp_max_it 3
284: output_file: output/ex55_0.out
285: requires: triangle
287: test:
288: suffix: hypre
289: nsize: 4
290: requires: hypre !complex !defined(PETSC_HAVE_HYPRE_DEVICE)
291: args: -ne 29 -alpha 1.e-3 -ksp_type cg -pc_type hypre -pc_hypre_type boomeramg -ksp_monitor_short
293: # command line options match GPU defaults
294: test:
295: suffix: hypre_device
296: nsize: 4
297: requires: hypre !complex
298: args: -mat_type hypre -ksp_view -ne 29 -alpha 1.e-3 -ksp_type cg -pc_type hypre -pc_hypre_type boomeramg -ksp_monitor_short -pc_hypre_boomeramg_relax_type_all l1scaled-Jacobi -pc_hypre_boomeramg_interp_type ext+i -pc_hypre_boomeramg_coarsen_type PMIS -pc_hypre_boomeramg_no_CF
300: TEST*/