Actual source code: ex55.c
petsc-3.8.4 2018-03-24
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,Pmat;
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 } };
36: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
37: comm = PETSC_COMM_WORLD;
38: MPI_Comm_rank(comm, &mype);
39: MPI_Comm_size(comm, &npe);
40: PetscOptionsGetInt(NULL,NULL,"-ne",&ne,NULL);
41: h = 1./ne;
42: /* ne*ne; number of global elements */
43: PetscOptionsGetReal(NULL,NULL,"-alpha",&soft_alpha,NULL);
44: PetscOptionsGetBool(NULL,NULL,"-use_coordinates",&use_coords,NULL);
45: M = 2*(ne+1)*(ne+1); /* global number of equations */
46: m = (ne+1)*(ne+1)/npe;
47: if (mype==npe-1) m = (ne+1)*(ne+1) - (npe-1)*m;
48: m *= 2;
49: /* create stiffness matrix */
50: MatCreate(comm,&Amat);
51: MatSetSizes(Amat,m,m,M,M);
52: MatSetBlockSize(Amat,2);
53: MatSetType(Amat,MATAIJ);
54: MatSeqAIJSetPreallocation(Amat,18,NULL);
55: MatMPIAIJSetPreallocation(Amat,18,NULL,18,NULL);
57: MatCreate(comm,&Pmat);
58: MatSetSizes(Pmat,m,m,M,M);
59: MatSetBlockSize(Pmat,2);
60: MatSetType(Pmat,MATAIJ);
61: MatSeqAIJSetPreallocation(Pmat,18,NULL);
62: MatMPIAIJSetPreallocation(Pmat,18,NULL,12,NULL);
64: MatGetOwnershipRange(Amat,&Istart,&Iend);
65: if (m != Iend - Istart) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"m %D does not equal Iend %D - Istart %D",m,Iend,Istart);
66: /* Generate vectors */
67: VecCreate(comm,&xx);
68: VecSetSizes(xx,m,M);
69: VecSetFromOptions(xx);
70: VecDuplicate(xx,&bb);
71: VecSet(bb,.0);
72: /* generate element matrices -- see ex56.c on how to use different data set */
73: {
74: DD[0][0] = 0.53333333333333321;
75: DD[0][1] = 0.20000000000000001;
76: DD[0][2] = -0.33333333333333331;
77: DD[0][3] = 0.0000000000000000;
78: DD[0][4] = -0.26666666666666666;
79: DD[0][5] = -0.20000000000000001;
80: DD[0][6] = 6.66666666666666796E-002;
81: DD[0][7] = 6.93889390390722838E-018;
82: DD[1][0] = 0.20000000000000001;
83: DD[1][1] = 0.53333333333333333;
84: DD[1][2] = 7.80625564189563192E-018;
85: DD[1][3] = 6.66666666666666935E-002;
86: DD[1][4] = -0.20000000000000001;
87: DD[1][5] = -0.26666666666666666;
88: DD[1][6] = -3.46944695195361419E-018;
89: DD[1][7] = -0.33333333333333331;
90: DD[2][0] = -0.33333333333333331;
91: DD[2][1] = 1.12757025938492461E-017;
92: DD[2][2] = 0.53333333333333333;
93: DD[2][3] = -0.20000000000000001;
94: DD[2][4] = 6.66666666666666935E-002;
95: DD[2][5] = -6.93889390390722838E-018;
96: DD[2][6] = -0.26666666666666666;
97: DD[2][7] = 0.19999999999999998;
98: DD[3][0] = 0.0000000000000000;
99: DD[3][1] = 6.66666666666666935E-002;
100: DD[3][2] = -0.20000000000000001;
101: DD[3][3] = 0.53333333333333333;
102: DD[3][4] = 4.33680868994201774E-018;
103: DD[3][5] = -0.33333333333333331;
104: DD[3][6] = 0.20000000000000001;
105: DD[3][7] = -0.26666666666666666;
106: DD[4][0] = -0.26666666666666666;
107: DD[4][1] = -0.20000000000000001;
108: DD[4][2] = 6.66666666666666935E-002;
109: DD[4][3] = 8.67361737988403547E-019;
110: DD[4][4] = 0.53333333333333333;
111: DD[4][5] = 0.19999999999999998;
112: DD[4][6] = -0.33333333333333331;
113: DD[4][7] = -3.46944695195361419E-018;
114: DD[5][0] = -0.20000000000000001;
115: DD[5][1] = -0.26666666666666666;
116: DD[5][2] = -1.04083408558608426E-017;
117: DD[5][3] = -0.33333333333333331;
118: DD[5][4] = 0.19999999999999998;
119: DD[5][5] = 0.53333333333333333;
120: DD[5][6] = 6.93889390390722838E-018;
121: DD[5][7] = 6.66666666666666519E-002;
122: DD[6][0] = 6.66666666666666796E-002;
123: DD[6][1] = -6.93889390390722838E-018;
124: DD[6][2] = -0.26666666666666666;
125: DD[6][3] = 0.19999999999999998;
126: DD[6][4] = -0.33333333333333331;
127: DD[6][5] = 6.93889390390722838E-018;
128: DD[6][6] = 0.53333333333333321;
129: DD[6][7] = -0.20000000000000001;
130: DD[7][0] = 6.93889390390722838E-018;
131: DD[7][1] = -0.33333333333333331;
132: DD[7][2] = 0.19999999999999998;
133: DD[7][3] = -0.26666666666666666;
134: DD[7][4] = 0.0000000000000000;
135: DD[7][5] = 6.66666666666666519E-002;
136: DD[7][6] = -0.20000000000000001;
137: DD[7][7] = 0.53333333333333321;
139: /* BC version of element */
140: for (i=0; i<8; i++) {
141: for (j=0; j<8; j++) {
142: if (i<4 || j < 4) {
143: if (i==j) DD2[i][j] = .1*DD1[i][j];
144: else DD2[i][j] = 0.0;
145: } else DD2[i][j] = DD1[i][j];
146: }
147: }
148: }
149: {
150: PetscReal *coords;
151: PetscMalloc1(m,&coords);
152: /* forms the element stiffness and coordinates */
153: for (Ii = Istart/2, ix = 0; Ii < Iend/2; Ii++, ix++) {
154: j = Ii/(ne+1); i = Ii%(ne+1);
155: /* coords */
156: x = h*(Ii % (ne+1)); y = h*(Ii/(ne+1));
157: coords[2*ix] = x; coords[2*ix+1] = y;
158: if (i<ne && j<ne) {
159: PetscInt jj,ii,idx[4];
160: /* radius */
161: PetscReal radius = PetscSqrtReal((x-.5+h/2)*(x-.5+h/2) + (y-.5+h/2)*(y-.5+h/2));
162: PetscReal alpha = 1.0;
163: if (radius < 0.25) alpha = soft_alpha;
165: idx[0] = Ii; idx[1] = Ii+1; idx[2] = Ii + (ne+1) + 1; idx[3] = Ii + (ne+1);
166: for (ii=0; ii<8; ii++) {
167: for (jj=0;jj<8;jj++) DD[ii][jj] = alpha*DD1[ii][jj];
168: }
169: MatSetValuesBlocked(Pmat,4,idx,4,idx,(const PetscScalar*)DD,ADD_VALUES);
170: if (j>0) {
171: MatSetValuesBlocked(Amat,4,idx,4,idx,(const PetscScalar*)DD,ADD_VALUES);
172: } else {
173: /* a BC */
174: for (ii=0; ii<8; ii++) {
175: for (jj=0;jj<8;jj++) DD[ii][jj] = alpha*DD2[ii][jj];
176: }
177: MatSetValuesBlocked(Amat,4,idx,4,idx,(const PetscScalar*)DD,ADD_VALUES);
178: }
179: }
180: if (j>0) {
181: PetscScalar v = h*h;
182: PetscInt jj = 2*Ii; /* load in x direction */
183: VecSetValues(bb,1,&jj,&v,INSERT_VALUES);
184: }
185: }
186: MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);
187: MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);
188: MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY);
189: MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY);
190: VecAssemblyBegin(bb);
191: VecAssemblyEnd(bb);
193: /* Setup solver */
194: KSPCreate(PETSC_COMM_WORLD,&ksp);
195: KSPSetFromOptions(ksp);
197: /* finish KSP/PC setup */
198: KSPSetOperators(ksp, Amat, Amat);
199: if (use_coords) {
200: PC pc;
201: KSPGetPC(ksp, &pc);
202: PCSetCoordinates(pc, 2, m/2, coords);
203: }
204: PetscFree(coords);
205: }
207: if (!PETSC_TRUE) {
208: PetscViewer viewer;
209: PetscViewerASCIIOpen(comm, "Amat.m", &viewer);
210: PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
211: MatView(Amat,viewer);
212: PetscViewerPopFormat(viewer);
213: PetscViewerDestroy(&viewer);
214: }
216: /* solve */
217: #if defined(PETSC_USE_LOG)
218: PetscLogStageRegister("Setup", &stage[0]);
219: PetscLogStageRegister("Solve", &stage[1]);
220: PetscLogStagePush(stage[0]);
221: #endif
222: KSPSetUp(ksp);
223: #if defined(PETSC_USE_LOG)
224: PetscLogStagePop();
225: #endif
227: VecSet(xx,.0);
229: #if defined(PETSC_USE_LOG)
230: PetscLogStagePush(stage[1]);
231: #endif
232: KSPSolve(ksp, bb, xx);
233: #if defined(PETSC_USE_LOG)
234: PetscLogStagePop();
235: #endif
237: KSPGetIterationNumber(ksp,&its);
239: if (0) {
240: PetscReal norm,norm2;
241: PetscViewer viewer;
242: Vec res;
244: PetscObjectGetComm((PetscObject)bb,&comm);
245: VecNorm(bb, NORM_2, &norm2);
247: VecDuplicate(xx, &res);
248: MatMult(Amat, xx, res);
249: VecAXPY(bb, -1.0, res);
250: VecDestroy(&res);
251: VecNorm(bb, NORM_2, &norm);
252: PetscPrintf(PETSC_COMM_WORLD,"[%d]%s |b-Ax|/|b|=%e, |b|=%e\n",0,PETSC_FUNCTION_NAME,norm/norm2,norm2);
253: PetscViewerASCIIOpen(comm, "residual.m", &viewer);
254: PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
255: VecView(bb,viewer);
256: PetscViewerPopFormat(viewer);
257: PetscViewerDestroy(&viewer);
258: }
260: /* Free work space */
261: KSPDestroy(&ksp);
262: VecDestroy(&xx);
263: VecDestroy(&bb);
264: MatDestroy(&Amat);
265: MatDestroy(&Pmat);
267: PetscFinalize();
268: return ierr;
269: }