Actual source code: ex59.c
petsc-3.13.6 2020-09-29
1: static char help[] = "This example illustrates the use of PCBDDC/FETI-DP and its customization.\n\n\
2: Discrete system: 1D, 2D or 3D laplacian, discretized with spectral elements.\n\
3: Spectral degree can be specified by passing values to -p option.\n\
4: Global problem either with dirichlet boundary conditions on one side or in the pure neumann case (depending on runtime parameters).\n\
5: Domain is [-nex,nex]x[-ney,ney]x[-nez,nez]: ne_ number of elements in _ direction.\n\
6: Exaple usage: \n\
7: 1D: mpiexec -n 4 ex59 -nex 7\n\
8: 2D: mpiexec -n 4 ex59 -npx 2 -npy 2 -nex 2 -ney 2\n\
9: 3D: mpiexec -n 4 ex59 -npx 2 -npy 2 -npz 1 -nex 2 -ney 2 -nez 1\n\
10: Subdomain decomposition can be specified with -np_ parameters.\n\
11: Dirichlet boundaries on one side by default:\n\
12: it does not iterate on dirichlet nodes by default: if -usezerorows is passed in, it also iterates on Dirichlet nodes.\n\
13: Pure Neumann case can be requested by passing in -pureneumann.\n\
14: In the latter case, in order to avoid runtime errors during factorization, please specify also -coarse_redundant_pc_factor_zeropivot 0\n\n";
16: #include <petscksp.h>
17: #include <petscpc.h>
18: #include <petscdm.h>
19: #include <petscdmda.h>
20: #include <petscblaslapack.h>
21: #define DEBUG 0
23: /* structure holding domain data */
24: typedef struct {
25: /* communicator */
26: MPI_Comm gcomm;
27: /* space dimension */
28: PetscInt dim;
29: /* spectral degree */
30: PetscInt p;
31: /* subdomains per dimension */
32: PetscInt npx,npy,npz;
33: /* subdomain index in cartesian dimensions */
34: PetscInt ipx,ipy,ipz;
35: /* elements per dimension */
36: PetscInt nex,ney,nez;
37: /* local elements per dimension */
38: PetscInt nex_l,ney_l,nez_l;
39: /* global number of dofs per dimension */
40: PetscInt xm,ym,zm;
41: /* local number of dofs per dimension */
42: PetscInt xm_l,ym_l,zm_l;
43: /* starting global indexes for subdomain in lexicographic ordering */
44: PetscInt startx,starty,startz;
45: /* pure Neumann problem */
46: PetscBool pure_neumann;
47: /* Dirichlet BC implementation */
48: PetscBool DBC_zerorows;
49: /* Scaling factor for subdomain */
50: PetscScalar scalingfactor;
51: PetscBool testkspfetidp;
52: } DomainData;
54: /* structure holding GLL data */
55: typedef struct {
56: /* GLL nodes */
57: PetscReal *zGL;
58: /* GLL weights */
59: PetscScalar *rhoGL;
60: /* aux_mat */
61: PetscScalar **A;
62: /* Element matrix */
63: Mat elem_mat;
64: } GLLData;
67: static PetscErrorCode BuildCSRGraph(DomainData dd, PetscInt **xadj, PetscInt **adjncy)
68: {
70: PetscInt *xadj_temp,*adjncy_temp;
71: PetscInt i,j,k,ii,jj,kk,iindex,count_adj;
72: PetscInt istart_csr,iend_csr,jstart_csr,jend_csr,kstart_csr,kend_csr;
73: PetscBool internal_node;
75: /* first count dimension of adjncy */
76: count_adj=0;
77: for (k=0; k<dd.zm_l; k++) {
78: internal_node = PETSC_TRUE;
79: kstart_csr = k>0 ? k-1 : k;
80: kend_csr = k<dd.zm_l-1 ? k+2 : k+1;
81: if (k == 0 || k == dd.zm_l-1) {
82: internal_node = PETSC_FALSE;
83: kstart_csr = k;
84: kend_csr = k+1;
85: }
86: for (j=0; j<dd.ym_l; j++) {
87: jstart_csr = j>0 ? j-1 : j;
88: jend_csr = j<dd.ym_l-1 ? j+2 : j+1;
89: if (j == 0 || j == dd.ym_l-1) {
90: internal_node = PETSC_FALSE;
91: jstart_csr = j;
92: jend_csr = j+1;
93: }
94: for (i=0; i<dd.xm_l; i++) {
95: istart_csr = i>0 ? i-1 : i;
96: iend_csr = i<dd.xm_l-1 ? i+2 : i+1;
97: if (i == 0 || i == dd.xm_l-1) {
98: internal_node = PETSC_FALSE;
99: istart_csr = i;
100: iend_csr = i+1;
101: }
102: if (internal_node) {
103: istart_csr = i;
104: iend_csr = i+1;
105: jstart_csr = j;
106: jend_csr = j+1;
107: kstart_csr = k;
108: kend_csr = k+1;
109: }
110: for (kk=kstart_csr;kk<kend_csr;kk++) {
111: for (jj=jstart_csr;jj<jend_csr;jj++) {
112: for (ii=istart_csr;ii<iend_csr;ii++) {
113: count_adj=count_adj+1;
114: }
115: }
116: }
117: }
118: }
119: }
120: PetscMalloc1(dd.xm_l*dd.ym_l*dd.zm_l+1,&xadj_temp);
121: PetscMalloc1(count_adj,&adjncy_temp);
123: /* now fill CSR data structure */
124: count_adj=0;
125: for (k=0; k<dd.zm_l; k++) {
126: internal_node = PETSC_TRUE;
127: kstart_csr = k>0 ? k-1 : k;
128: kend_csr = k<dd.zm_l-1 ? k+2 : k+1;
129: if (k == 0 || k == dd.zm_l-1) {
130: internal_node = PETSC_FALSE;
131: kstart_csr = k;
132: kend_csr = k+1;
133: }
134: for (j=0; j<dd.ym_l; j++) {
135: jstart_csr = j>0 ? j-1 : j;
136: jend_csr = j<dd.ym_l-1 ? j+2 : j+1;
137: if (j == 0 || j == dd.ym_l-1) {
138: internal_node = PETSC_FALSE;
139: jstart_csr = j;
140: jend_csr = j+1;
141: }
142: for (i=0; i<dd.xm_l; i++) {
143: istart_csr = i>0 ? i-1 : i;
144: iend_csr = i<dd.xm_l-1 ? i+2 : i+1;
145: if (i == 0 || i == dd.xm_l-1) {
146: internal_node = PETSC_FALSE;
147: istart_csr = i;
148: iend_csr = i+1;
149: }
150: iindex = k*dd.xm_l*dd.ym_l+j*dd.xm_l+i;
151: xadj_temp[iindex] = count_adj;
152: if (internal_node) {
153: istart_csr = i;
154: iend_csr = i+1;
155: jstart_csr = j;
156: jend_csr = j+1;
157: kstart_csr = k;
158: kend_csr = k+1;
159: }
160: for (kk=kstart_csr; kk<kend_csr; kk++) {
161: for (jj=jstart_csr; jj<jend_csr; jj++) {
162: for (ii=istart_csr; ii<iend_csr; ii++) {
163: iindex = kk*dd.xm_l*dd.ym_l+jj*dd.xm_l+ii;
165: adjncy_temp[count_adj] = iindex;
166: count_adj = count_adj+1;
167: }
168: }
169: }
170: }
171: }
172: }
173: xadj_temp[dd.xm_l*dd.ym_l*dd.zm_l] = count_adj;
175: *xadj = xadj_temp;
176: *adjncy = adjncy_temp;
177: return(0);
178: }
180: static PetscErrorCode ComputeSpecialBoundaryIndices(DomainData dd,IS *dirichlet,IS *neumann)
181: {
183: IS temp_dirichlet=0,temp_neumann=0;
184: PetscInt localsize,i,j,k,*indices;
185: PetscBool *touched;
188: localsize = dd.xm_l*dd.ym_l*dd.zm_l;
190: PetscMalloc1(localsize,&indices);
191: PetscMalloc1(localsize,&touched);
192: for (i=0; i<localsize; i++) touched[i] = PETSC_FALSE;
194: if (dirichlet) {
195: i = 0;
196: /* west boundary */
197: if (dd.ipx == 0) {
198: for (k=0;k<dd.zm_l;k++) {
199: for (j=0;j<dd.ym_l;j++) {
200: indices[i]=k*dd.ym_l*dd.xm_l+j*dd.xm_l;
201: touched[indices[i]]=PETSC_TRUE;
202: i++;
203: }
204: }
205: }
206: ISCreateGeneral(dd.gcomm,i,indices,PETSC_COPY_VALUES,&temp_dirichlet);
207: }
208: if (neumann) {
209: i = 0;
210: /* west boundary */
211: if (dd.ipx == 0) {
212: for (k=0;k<dd.zm_l;k++) {
213: for (j=0;j<dd.ym_l;j++) {
214: indices[i]=k*dd.ym_l*dd.xm_l+j*dd.xm_l;
215: if (!touched[indices[i]]) {
216: touched[indices[i]]=PETSC_TRUE;
217: i++;
218: }
219: }
220: }
221: }
222: /* east boundary */
223: if (dd.ipx == dd.npx-1) {
224: for (k=0;k<dd.zm_l;k++) {
225: for (j=0;j<dd.ym_l;j++) {
226: indices[i]=k*dd.ym_l*dd.xm_l+j*dd.xm_l+dd.xm_l-1;
227: if (!touched[indices[i]]) {
228: touched[indices[i]]=PETSC_TRUE;
229: i++;
230: }
231: }
232: }
233: }
234: /* south boundary */
235: if (dd.ipy == 0 && dd.dim > 1) {
236: for (k=0;k<dd.zm_l;k++) {
237: for (j=0;j<dd.xm_l;j++) {
238: indices[i]=k*dd.ym_l*dd.xm_l+j;
239: if (!touched[indices[i]]) {
240: touched[indices[i]]=PETSC_TRUE;
241: i++;
242: }
243: }
244: }
245: }
246: /* north boundary */
247: if (dd.ipy == dd.npy-1 && dd.dim > 1) {
248: for (k=0;k<dd.zm_l;k++) {
249: for (j=0;j<dd.xm_l;j++) {
250: indices[i]=k*dd.ym_l*dd.xm_l+(dd.ym_l-1)*dd.xm_l+j;
251: if (!touched[indices[i]]) {
252: touched[indices[i]]=PETSC_TRUE;
253: i++;
254: }
255: }
256: }
257: }
258: /* bottom boundary */
259: if (dd.ipz == 0 && dd.dim > 2) {
260: for (k=0;k<dd.ym_l;k++) {
261: for (j=0;j<dd.xm_l;j++) {
262: indices[i]=k*dd.xm_l+j;
263: if (!touched[indices[i]]) {
264: touched[indices[i]]=PETSC_TRUE;
265: i++;
266: }
267: }
268: }
269: }
270: /* top boundary */
271: if (dd.ipz == dd.npz-1 && dd.dim > 2) {
272: for (k=0;k<dd.ym_l;k++) {
273: for (j=0;j<dd.xm_l;j++) {
274: indices[i]=(dd.zm_l-1)*dd.ym_l*dd.xm_l+k*dd.xm_l+j;
275: if (!touched[indices[i]]) {
276: touched[indices[i]]=PETSC_TRUE;
277: i++;
278: }
279: }
280: }
281: }
282: ISCreateGeneral(dd.gcomm,i,indices,PETSC_COPY_VALUES,&temp_neumann);
283: }
284: if (dirichlet) *dirichlet = temp_dirichlet;
285: if (neumann) *neumann = temp_neumann;
286: PetscFree(indices);
287: PetscFree(touched);
288: return(0);
289: }
291: static PetscErrorCode ComputeMapping(DomainData dd,ISLocalToGlobalMapping *isg2lmap)
292: {
293: PetscErrorCode ierr;
294: DM da;
295: AO ao;
296: DMBoundaryType bx = DM_BOUNDARY_NONE,by = DM_BOUNDARY_NONE, bz = DM_BOUNDARY_NONE;
297: DMDAStencilType stype = DMDA_STENCIL_BOX;
298: ISLocalToGlobalMapping temp_isg2lmap;
299: PetscInt i,j,k,ig,jg,kg,lindex,gindex,localsize;
300: PetscInt *global_indices;
303: /* Not an efficient mapping: this function computes a very simple lexicographic mapping
304: just to illustrate the creation of a MATIS object */
305: localsize = dd.xm_l*dd.ym_l*dd.zm_l;
306: PetscMalloc1(localsize,&global_indices);
307: for (k=0; k<dd.zm_l; k++) {
308: kg=dd.startz+k;
309: for (j=0; j<dd.ym_l; j++) {
310: jg=dd.starty+j;
311: for (i=0; i<dd.xm_l; i++) {
312: ig =dd.startx+i;
313: lindex =k*dd.xm_l*dd.ym_l+j*dd.xm_l+i;
314: gindex =kg*dd.xm*dd.ym+jg*dd.xm+ig;
315: global_indices[lindex]=gindex;
316: }
317: }
318: }
319: if (dd.dim==3) {
320: DMDACreate3d(dd.gcomm,bx,by,bz,stype,dd.xm,dd.ym,dd.zm,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,NULL,&da);
321: } else if (dd.dim==2) {
322: DMDACreate2d(dd.gcomm,bx,by,stype,dd.xm,dd.ym,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
323: } else {
324: DMDACreate1d(dd.gcomm,bx,dd.xm,1,1,NULL,&da);
325: }
326: DMSetFromOptions(da);
327: DMSetUp(da);
328: DMDASetAOType(da,AOMEMORYSCALABLE);
329: DMDAGetAO(da,&ao);
330: AOApplicationToPetsc(ao,dd.xm_l*dd.ym_l*dd.zm_l,global_indices);
331: ISLocalToGlobalMappingCreate(dd.gcomm,1,localsize,global_indices,PETSC_OWN_POINTER,&temp_isg2lmap);
332: DMDestroy(&da);
333: *isg2lmap = temp_isg2lmap;
334: return(0);
335: }
336: static PetscErrorCode ComputeSubdomainMatrix(DomainData dd, GLLData glldata, Mat *local_mat)
337: {
339: PetscInt localsize,zloc,yloc,xloc,auxnex,auxney,auxnez;
340: PetscInt ie,je,ke,i,j,k,ig,jg,kg,ii,ming;
341: PetscInt *indexg,*cols,*colsg;
342: PetscScalar *vals;
343: Mat temp_local_mat,elem_mat_DBC=0,*usedmat;
344: IS submatIS;
347: MatGetSize(glldata.elem_mat,&i,&j);
348: PetscMalloc1(i,&indexg);
349: PetscMalloc1(i,&colsg);
350: /* get submatrix of elem_mat without dirichlet nodes */
351: if (!dd.pure_neumann && !dd.DBC_zerorows && !dd.ipx) {
352: xloc = dd.p+1;
353: yloc = 1;
354: zloc = 1;
355: if (dd.dim>1) yloc = dd.p+1;
356: if (dd.dim>2) zloc = dd.p+1;
357: ii = 0;
358: for (k=0;k<zloc;k++) {
359: for (j=0;j<yloc;j++) {
360: for (i=1;i<xloc;i++) {
361: indexg[ii]=k*xloc*yloc+j*xloc+i;
362: ii++;
363: }
364: }
365: }
366: ISCreateGeneral(PETSC_COMM_SELF,ii,indexg,PETSC_COPY_VALUES,&submatIS);
367: MatCreateSubMatrix(glldata.elem_mat,submatIS,submatIS,MAT_INITIAL_MATRIX,&elem_mat_DBC);
368: ISDestroy(&submatIS);
369: }
371: /* Assemble subdomain matrix */
372: localsize = dd.xm_l*dd.ym_l*dd.zm_l;
373: MatCreate(PETSC_COMM_SELF,&temp_local_mat);
374: MatSetSizes(temp_local_mat,localsize,localsize,localsize,localsize);
375: MatSetOptionsPrefix(temp_local_mat,"subdomain_");
376: /* set local matrices type: here we use SEQSBAIJ primarily for testing purpose */
377: /* in order to avoid conversions inside the BDDC code, use SeqAIJ if possible */
378: if (dd.DBC_zerorows && !dd.ipx) { /* in this case, we need to zero out some of the rows, so use seqaij */
379: MatSetType(temp_local_mat,MATSEQAIJ);
380: } else {
381: MatSetType(temp_local_mat,MATSEQSBAIJ);
382: }
383: MatSetFromOptions(temp_local_mat);
385: i = PetscPowRealInt(3.0*(dd.p+1.0),dd.dim);
387: MatSeqAIJSetPreallocation(temp_local_mat,i,NULL); /* very overestimated */
388: MatSeqSBAIJSetPreallocation(temp_local_mat,1,i,NULL); /* very overestimated */
389: MatSeqBAIJSetPreallocation(temp_local_mat,1,i,NULL); /* very overestimated */
390: MatSetOption(temp_local_mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
392: yloc = dd.p+1;
393: zloc = dd.p+1;
394: if (dd.dim < 3) zloc = 1;
395: if (dd.dim < 2) yloc = 1;
397: auxnez = dd.nez_l;
398: auxney = dd.ney_l;
399: auxnex = dd.nex_l;
400: if (dd.dim < 3) auxnez = 1;
401: if (dd.dim < 2) auxney = 1;
403: for (ke=0; ke<auxnez; ke++) {
404: for (je=0; je<auxney; je++) {
405: for (ie=0; ie<auxnex; ie++) {
406: /* customize element accounting for BC */
407: xloc = dd.p+1;
408: ming = 0;
409: usedmat = &glldata.elem_mat;
410: if (!dd.pure_neumann && !dd.DBC_zerorows && !dd.ipx) {
411: if (ie == 0) {
412: xloc = dd.p;
413: usedmat = &elem_mat_DBC;
414: } else {
415: ming = -1;
416: usedmat = &glldata.elem_mat;
417: }
418: }
419: /* local to the element/global to the subdomain indexing */
420: for (k=0; k<zloc; k++) {
421: kg = ke*dd.p+k;
422: for (j=0; j<yloc; j++) {
423: jg = je*dd.p+j;
424: for (i=0; i<xloc; i++) {
425: ig = ie*dd.p+i+ming;
426: ii = k*xloc*yloc+j*xloc+i;
427: indexg[ii] = kg*dd.xm_l*dd.ym_l+jg*dd.xm_l+ig;
428: }
429: }
430: }
431: /* Set values */
432: for (i=0; i<xloc*yloc*zloc; i++) {
433: MatGetRow(*usedmat,i,&j,(const PetscInt**)&cols,(const PetscScalar**)&vals);
434: for (k=0; k<j; k++) colsg[k] = indexg[cols[k]];
435: MatSetValues(temp_local_mat,1,&indexg[i],j,colsg,vals,ADD_VALUES);
436: MatRestoreRow(*usedmat,i,&j,(const PetscInt**)&cols,(const PetscScalar**)&vals);
437: }
438: }
439: }
440: }
441: PetscFree(indexg);
442: PetscFree(colsg);
443: MatAssemblyBegin(temp_local_mat,MAT_FINAL_ASSEMBLY);
444: MatAssemblyEnd (temp_local_mat,MAT_FINAL_ASSEMBLY);
445: #if DEBUG
446: {
447: Vec lvec,rvec;
448: PetscReal norm;
449: MatCreateVecs(temp_local_mat,&lvec,&rvec);
450: VecSet(lvec,1.0);
451: MatMult(temp_local_mat,lvec,rvec);
452: VecNorm(rvec,NORM_INFINITY,&norm);
453: VecDestroy(&lvec);
454: VecDestroy(&rvec);
455: }
456: #endif
457: *local_mat = temp_local_mat;
458: MatDestroy(&elem_mat_DBC);
459: return(0);
460: }
462: static PetscErrorCode GLLStuffs(DomainData dd, GLLData *glldata)
463: {
465: PetscReal *M,si;
466: PetscScalar x,z0,z1,z2,Lpj,Lpr,rhoGLj,rhoGLk;
467: PetscBLASInt pm1,lierr;
468: PetscInt i,j,n,k,s,r,q,ii,jj,p=dd.p;
469: PetscInt xloc,yloc,zloc,xyloc,xyzloc;
472: /* Gauss-Lobatto-Legendre nodes zGL on [-1,1] */
473: PetscCalloc1(p+1,&glldata->zGL);
475: glldata->zGL[0]=-1.0;
476: glldata->zGL[p]= 1.0;
477: if (p > 1) {
478: if (p == 2) glldata->zGL[1]=0.0;
479: else {
480: PetscMalloc1(p-1,&M);
481: for (i=0; i<p-1; i++) {
482: si = (PetscReal)(i+1.0);
483: M[i]=0.5*PetscSqrtReal(si*(si+2.0)/((si+0.5)*(si+1.5)));
484: }
485: pm1 = p-1;
486: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
487: PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("N",&pm1,&glldata->zGL[1],M,&x,&pm1,M,&lierr));
488: if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in STERF Lapack routine %d",(int)lierr);
489: PetscFPTrapPop();
490: PetscFree(M);
491: }
492: }
494: /* Weights for 1D quadrature */
495: PetscMalloc1(p+1,&glldata->rhoGL);
497: glldata->rhoGL[0]=2.0/(PetscScalar)(p*(p+1.0));
498: glldata->rhoGL[p]=glldata->rhoGL[0];
499: z2 = -1; /* Dummy value to avoid -Wmaybe-initialized */
500: for (i=1; i<p; i++) {
501: x = glldata->zGL[i];
502: z0 = 1.0;
503: z1 = x;
504: for (n=1; n<p; n++) {
505: z2 = x*z1*(2.0*n+1.0)/(n+1.0)-z0*(PetscScalar)(n/(n+1.0));
506: z0 = z1;
507: z1 = z2;
508: }
509: glldata->rhoGL[i]=2.0/(p*(p+1.0)*z2*z2);
510: }
512: /* Auxiliary mat for laplacian */
513: PetscMalloc1(p+1,&glldata->A);
514: PetscMalloc1((p+1)*(p+1),&glldata->A[0]);
515: for (i=1; i<p+1; i++) glldata->A[i]=glldata->A[i-1]+p+1;
517: for (j=1; j<p; j++) {
518: x =glldata->zGL[j];
519: z0=1.0;
520: z1=x;
521: for (n=1; n<p; n++) {
522: z2=x*z1*(2.0*n+1.0)/(n+1.0)-z0*(PetscScalar)(n/(n+1.0));
523: z0=z1;
524: z1=z2;
525: }
526: Lpj=z2;
527: for (r=1; r<p; r++) {
528: if (r == j) {
529: glldata->A[j][j]=2.0/(3.0*(1.0-glldata->zGL[j]*glldata->zGL[j])*Lpj*Lpj);
530: } else {
531: x = glldata->zGL[r];
532: z0 = 1.0;
533: z1 = x;
534: for (n=1; n<p; n++) {
535: z2=x*z1*(2.0*n+1.0)/(n+1.0)-z0*(PetscScalar)(n/(n+1.0));
536: z0=z1;
537: z1=z2;
538: }
539: Lpr = z2;
540: glldata->A[r][j]=4.0/(p*(p+1.0)*Lpj*Lpr*(glldata->zGL[j]-glldata->zGL[r])*(glldata->zGL[j]-glldata->zGL[r]));
541: }
542: }
543: }
544: for (j=1; j<p+1; j++) {
545: x = glldata->zGL[j];
546: z0 = 1.0;
547: z1 = x;
548: for (n=1; n<p; n++) {
549: z2=x*z1*(2.0*n+1.0)/(n+1.0)-z0*(PetscScalar)(n/(n+1.0));
550: z0=z1;
551: z1=z2;
552: }
553: Lpj = z2;
554: glldata->A[j][0]=4.0*PetscPowRealInt(-1.0,p)/(p*(p+1.0)*Lpj*(1.0+glldata->zGL[j])*(1.0+glldata->zGL[j]));
555: glldata->A[0][j]=glldata->A[j][0];
556: }
557: for (j=0; j<p; j++) {
558: x = glldata->zGL[j];
559: z0 = 1.0;
560: z1 = x;
561: for (n=1; n<p; n++) {
562: z2=x*z1*(2.0*n+1.0)/(n+1.0)-z0*(PetscScalar)(n/(n+1.0));
563: z0=z1;
564: z1=z2;
565: }
566: Lpj=z2;
568: glldata->A[p][j]=4.0/(p*(p+1.0)*Lpj*(1.0-glldata->zGL[j])*(1.0-glldata->zGL[j]));
569: glldata->A[j][p]=glldata->A[p][j];
570: }
571: glldata->A[0][0]=0.5+(p*(p+1.0)-2.0)/6.0;
572: glldata->A[p][p]=glldata->A[0][0];
574: /* compute element matrix */
575: xloc = p+1;
576: yloc = p+1;
577: zloc = p+1;
578: if (dd.dim<2) yloc=1;
579: if (dd.dim<3) zloc=1;
580: xyloc = xloc*yloc;
581: xyzloc = xloc*yloc*zloc;
583: MatCreate(PETSC_COMM_SELF,&glldata->elem_mat);
584: MatSetSizes(glldata->elem_mat,xyzloc,xyzloc,xyzloc,xyzloc);
585: MatSetType(glldata->elem_mat,MATSEQAIJ);
586: MatSeqAIJSetPreallocation(glldata->elem_mat,xyzloc,NULL); /* overestimated */
587: MatZeroEntries(glldata->elem_mat);
588: MatSetOption(glldata->elem_mat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);
590: for (k=0; k<zloc; k++) {
591: if (dd.dim>2) rhoGLk=glldata->rhoGL[k];
592: else rhoGLk=1.0;
594: for (j=0; j<yloc; j++) {
595: if (dd.dim>1) rhoGLj=glldata->rhoGL[j];
596: else rhoGLj=1.0;
598: for (i=0; i<xloc; i++) {
599: ii = k*xyloc+j*xloc+i;
600: s = k;
601: r = j;
602: for (q=0; q<xloc; q++) {
603: jj = s*xyloc+r*xloc+q;
604: MatSetValue(glldata->elem_mat,jj,ii,glldata->A[i][q]*rhoGLj*rhoGLk,ADD_VALUES);
605: }
606: if (dd.dim>1) {
607: s=k;
608: q=i;
609: for (r=0; r<yloc; r++) {
610: jj = s*xyloc+r*xloc+q;
611: MatSetValue(glldata->elem_mat,jj,ii,glldata->A[j][r]*glldata->rhoGL[i]*rhoGLk,ADD_VALUES);
612: }
613: }
614: if (dd.dim>2) {
615: r=j;
616: q=i;
617: for (s=0; s<zloc; s++) {
618: jj = s*xyloc+r*xloc+q;
619: MatSetValue(glldata->elem_mat,jj,ii,glldata->A[k][s]*rhoGLj*glldata->rhoGL[i],ADD_VALUES);
620: }
621: }
622: }
623: }
624: }
625: MatAssemblyBegin(glldata->elem_mat,MAT_FINAL_ASSEMBLY);
626: MatAssemblyEnd (glldata->elem_mat,MAT_FINAL_ASSEMBLY);
627: #if DEBUG
628: {
629: Vec lvec,rvec;
630: PetscReal norm;
631: MatCreateVecs(glldata->elem_mat,&lvec,&rvec);
632: VecSet(lvec,1.0);
633: MatMult(glldata->elem_mat,lvec,rvec);
634: VecNorm(rvec,NORM_INFINITY,&norm);
635: VecDestroy(&lvec);
636: VecDestroy(&rvec);
637: }
638: #endif
639: return(0);
640: }
642: static PetscErrorCode DomainDecomposition(DomainData *dd)
643: {
644: PetscMPIInt rank;
645: PetscInt i,j,k;
648: /* Subdomain index in cartesian coordinates */
649: MPI_Comm_rank(dd->gcomm,&rank);
650: dd->ipx = rank%dd->npx;
651: if (dd->dim>1) dd->ipz = rank/(dd->npx*dd->npy);
652: else dd->ipz = 0;
654: dd->ipy = rank/dd->npx-dd->ipz*dd->npy;
655: /* number of local elements */
656: dd->nex_l = dd->nex/dd->npx;
657: if (dd->ipx < dd->nex%dd->npx) dd->nex_l++;
658: if (dd->dim>1) {
659: dd->ney_l = dd->ney/dd->npy;
660: if (dd->ipy < dd->ney%dd->npy) dd->ney_l++;
661: } else {
662: dd->ney_l = 0;
663: }
664: if (dd->dim>2) {
665: dd->nez_l = dd->nez/dd->npz;
666: if (dd->ipz < dd->nez%dd->npz) dd->nez_l++;
667: } else {
668: dd->nez_l = 0;
669: }
670: /* local and global number of dofs */
671: dd->xm_l = dd->nex_l*dd->p+1;
672: dd->xm = dd->nex*dd->p+1;
673: dd->ym_l = dd->ney_l*dd->p+1;
674: dd->ym = dd->ney*dd->p+1;
675: dd->zm_l = dd->nez_l*dd->p+1;
676: dd->zm = dd->nez*dd->p+1;
677: if (!dd->pure_neumann) {
678: if (!dd->DBC_zerorows && !dd->ipx) dd->xm_l--;
679: if (!dd->DBC_zerorows) dd->xm--;
680: }
682: /* starting global index for local dofs (simple lexicographic order) */
683: dd->startx = 0;
684: j = dd->nex/dd->npx;
685: for (i=0; i<dd->ipx; i++) {
686: k = j;
687: if (i<dd->nex%dd->npx) k++;
688: dd->startx = dd->startx+k*dd->p;
689: }
690: if (!dd->pure_neumann && !dd->DBC_zerorows && dd->ipx) dd->startx--;
692: dd->starty = 0;
693: if (dd->dim > 1) {
694: j = dd->ney/dd->npy;
695: for (i=0; i<dd->ipy; i++) {
696: k = j;
697: if (i<dd->ney%dd->npy) k++;
698: dd->starty = dd->starty+k*dd->p;
699: }
700: }
701: dd->startz = 0;
702: if (dd->dim > 2) {
703: j = dd->nez/dd->npz;
704: for (i=0; i<dd->ipz; i++) {
705: k = j;
706: if (i<dd->nez%dd->npz) k++;
707: dd->startz = dd->startz+k*dd->p;
708: }
709: }
710: return(0);
711: }
712: static PetscErrorCode ComputeMatrix(DomainData dd, Mat *A)
713: {
714: PetscErrorCode ierr;
715: GLLData gll;
716: Mat local_mat =0,temp_A=0;
717: ISLocalToGlobalMapping matis_map =0;
718: IS dirichletIS=0;
721: /* Compute some stuff of Gauss-Legendre-Lobatto quadrature rule */
722: GLLStuffs(dd,&gll);
723: /* Compute matrix of subdomain Neumann problem */
724: ComputeSubdomainMatrix(dd,gll,&local_mat);
725: /* Compute global mapping of local dofs */
726: ComputeMapping(dd,&matis_map);
727: /* Create MATIS object needed by BDDC */
728: MatCreateIS(dd.gcomm,1,PETSC_DECIDE,PETSC_DECIDE,dd.xm*dd.ym*dd.zm,dd.xm*dd.ym*dd.zm,matis_map,NULL,&temp_A);
729: /* Set local subdomain matrices into MATIS object */
730: MatScale(local_mat,dd.scalingfactor);
731: MatISSetLocalMat(temp_A,local_mat);
732: /* Call assembly functions */
733: MatAssemblyBegin(temp_A,MAT_FINAL_ASSEMBLY);
734: MatAssemblyEnd(temp_A,MAT_FINAL_ASSEMBLY);
736: if (dd.DBC_zerorows) {
737: PetscInt dirsize;
739: ComputeSpecialBoundaryIndices(dd,&dirichletIS,NULL);
740: MatSetOption(local_mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
741: MatZeroRowsColumnsLocalIS(temp_A,dirichletIS,1.0,NULL,NULL);
742: ISGetLocalSize(dirichletIS,&dirsize);
743: ISDestroy(&dirichletIS);
744: }
746: /* giving hints to local and global matrices could be useful for the BDDC */
747: MatSetOption(local_mat,MAT_SYMMETRIC,PETSC_TRUE);
748: MatSetOption(local_mat,MAT_SPD,PETSC_TRUE);
749: #if DEBUG
750: {
751: Vec lvec,rvec;
752: PetscReal norm;
753: MatCreateVecs(temp_A,&lvec,&rvec);
754: VecSet(lvec,1.0);
755: MatMult(temp_A,lvec,rvec);
756: VecNorm(rvec,NORM_INFINITY,&norm);
757: VecDestroy(&lvec);
758: VecDestroy(&rvec);
759: }
760: #endif
761: /* free allocated workspace */
762: PetscFree(gll.zGL);
763: PetscFree(gll.rhoGL);
764: PetscFree(gll.A[0]);
765: PetscFree(gll.A);
766: MatDestroy(&gll.elem_mat);
767: MatDestroy(&local_mat);
768: ISLocalToGlobalMappingDestroy(&matis_map);
769: /* give back the pointer to te MATIS object */
770: *A = temp_A;
771: return(0);
772: }
774: static PetscErrorCode ComputeKSPFETIDP(DomainData dd, KSP ksp_bddc, KSP *ksp_fetidp)
775: {
777: KSP temp_ksp;
778: PC pc;
781: KSPGetPC(ksp_bddc,&pc);
782: if (!dd.testkspfetidp) {
783: PC D;
784: Mat F;
786: PCBDDCCreateFETIDPOperators(pc,PETSC_TRUE,NULL,&F,&D);
787: KSPCreate(PetscObjectComm((PetscObject)F),&temp_ksp);
788: KSPSetOperators(temp_ksp,F,F);
789: KSPSetType(temp_ksp,KSPCG);
790: KSPSetPC(temp_ksp,D);
791: KSPSetComputeSingularValues(temp_ksp,PETSC_TRUE);
792: KSPSetOptionsPrefix(temp_ksp,"fluxes_");
793: KSPSetFromOptions(temp_ksp);
794: KSPSetUp(temp_ksp);
795: MatDestroy(&F);
796: PCDestroy(&D);
797: } else {
798: Mat A,Ap;
800: KSPCreate(PetscObjectComm((PetscObject)ksp_bddc),&temp_ksp);
801: KSPGetOperators(ksp_bddc,&A,&Ap);
802: KSPSetOperators(temp_ksp,A,Ap);
803: KSPSetOptionsPrefix(temp_ksp,"fluxes_");
804: KSPSetType(temp_ksp,KSPFETIDP);
805: KSPFETIDPSetInnerBDDC(temp_ksp,pc);
806: KSPSetComputeSingularValues(temp_ksp,PETSC_TRUE);
807: KSPSetFromOptions(temp_ksp);
808: KSPSetUp(temp_ksp);
809: }
810: *ksp_fetidp = temp_ksp;
811: return(0);
812: }
815: static PetscErrorCode ComputeKSPBDDC(DomainData dd,Mat A,KSP *ksp)
816: {
818: KSP temp_ksp;
819: PC pc;
820: IS primals,dirichletIS=0,neumannIS=0,*bddc_dofs_splitting;
821: PetscInt vidx[8],localsize,*xadj=NULL,*adjncy=NULL;
822: MatNullSpace near_null_space;
825: KSPCreate(dd.gcomm,&temp_ksp);
826: KSPSetOperators(temp_ksp,A,A);
827: KSPSetType(temp_ksp,KSPCG);
828: KSPGetPC(temp_ksp,&pc);
829: PCSetType(pc,PCBDDC);
831: localsize = dd.xm_l*dd.ym_l*dd.zm_l;
833: /* BDDC customization */
835: /* jumping coefficients case */
836: PCISSetSubdomainScalingFactor(pc,dd.scalingfactor);
838: /* Dofs splitting
839: Simple stride-1 IS
840: It is not needed since, by default, PCBDDC assumes a stride-1 split */
841: PetscMalloc1(1,&bddc_dofs_splitting);
842: #if 1
843: ISCreateStride(PETSC_COMM_WORLD,localsize,0,1,&bddc_dofs_splitting[0]);
844: PCBDDCSetDofsSplittingLocal(pc,1,bddc_dofs_splitting);
845: #else
846: /* examples for global ordering */
848: /* each process lists the nodes it owns */
849: PetscInt sr,er;
850: MatGetOwnershipRange(A,&sr,&er);
851: ISCreateStride(PETSC_COMM_WORLD,er-sr,sr,1,&bddc_dofs_splitting[0]);
852: PCBDDCSetDofsSplitting(pc,1,bddc_dofs_splitting);
853: /* Split can be passed in a more general way since any process can list any node */
854: #endif
855: ISDestroy(&bddc_dofs_splitting[0]);
856: PetscFree(bddc_dofs_splitting);
858: /* Primal constraints implemented by using a near null space attached to A -> now it passes in only the constants
859: (which in practice is not needed since by default PCBDDC builds the primal space using constants for quadrature formulas */
860: #if 0
861: Vec vecs[2];
862: PetscRandom rctx;
863: MatCreateVecs(A,&vecs[0],&vecs[1]);
864: PetscRandomCreate(dd.gcomm,&rctx);
865: VecSetRandom(vecs[0],rctx);
866: VecSetRandom(vecs[1],rctx);
867: MatNullSpaceCreate(dd.gcomm,PETSC_TRUE,2,vecs,&near_null_space);
868: VecDestroy(&vecs[0]);
869: VecDestroy(&vecs[1]);
870: PetscRandomDestroy(&rctx);
871: #else
872: MatNullSpaceCreate(dd.gcomm,PETSC_TRUE,0,NULL,&near_null_space);
873: #endif
874: MatSetNearNullSpace(A,near_null_space);
875: MatNullSpaceDestroy(&near_null_space);
877: /* CSR graph of subdomain dofs */
878: BuildCSRGraph(dd,&xadj,&adjncy);
879: PCBDDCSetLocalAdjacencyGraph(pc,localsize,xadj,adjncy,PETSC_OWN_POINTER);
881: /* Prescribe user-defined primal vertices: in this case we use the 8 corners in 3D (for 2D and 1D, the indices are duplicated) */
882: vidx[0] = 0*dd.xm_l+0;
883: vidx[1] = 0*dd.xm_l+dd.xm_l-1;
884: vidx[2] = (dd.ym_l-1)*dd.xm_l+0;
885: vidx[3] = (dd.ym_l-1)*dd.xm_l+dd.xm_l-1;
886: vidx[4] = (dd.zm_l-1)*dd.xm_l*dd.ym_l+0*dd.xm_l+0;
887: vidx[5] = (dd.zm_l-1)*dd.xm_l*dd.ym_l+0*dd.xm_l+dd.xm_l-1;
888: vidx[6] = (dd.zm_l-1)*dd.xm_l*dd.ym_l+(dd.ym_l-1)*dd.xm_l+0;
889: vidx[7] = (dd.zm_l-1)*dd.xm_l*dd.ym_l+(dd.ym_l-1)*dd.xm_l+dd.xm_l-1;
890: ISCreateGeneral(dd.gcomm,8,vidx,PETSC_COPY_VALUES,&primals);
891: PCBDDCSetPrimalVerticesLocalIS(pc,primals);
892: ISDestroy(&primals);
894: /* Neumann/Dirichlet indices on the global boundary */
895: if (dd.DBC_zerorows) {
896: /* Only in case you eliminate some rows matrix with zerorows function, you need to set dirichlet indices into PCBDDC data */
897: ComputeSpecialBoundaryIndices(dd,&dirichletIS,&neumannIS);
898: PCBDDCSetNeumannBoundariesLocal(pc,neumannIS);
899: PCBDDCSetDirichletBoundariesLocal(pc,dirichletIS);
900: } else {
901: if (dd.pure_neumann) {
902: /* In such a case, all interface nodes lying on the global boundary are neumann nodes */
903: ComputeSpecialBoundaryIndices(dd,NULL,&neumannIS);
904: PCBDDCSetNeumannBoundariesLocal(pc,neumannIS);
905: } else {
906: /* It is wrong setting dirichlet indices without having zeroed the corresponding rows in the global matrix */
907: /* But we can still compute them since nodes near the dirichlet boundaries does not need to be defined as neumann nodes */
908: ComputeSpecialBoundaryIndices(dd,&dirichletIS,&neumannIS);
909: PCBDDCSetNeumannBoundariesLocal(pc,neumannIS);
910: }
911: }
913: /* Pass local null space information to local matrices (needed when using approximate local solvers) */
914: if (dd.ipx || dd.pure_neumann) {
915: MatNullSpace nsp;
916: Mat local_mat;
918: MatISGetLocalMat(A,&local_mat);
919: MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nsp);
920: MatSetNullSpace(local_mat,nsp);
921: MatNullSpaceDestroy(&nsp);
922: }
923: KSPSetComputeSingularValues(temp_ksp,PETSC_TRUE);
924: KSPSetOptionsPrefix(temp_ksp,"physical_");
925: KSPSetFromOptions(temp_ksp);
926: KSPSetUp(temp_ksp);
927: *ksp = temp_ksp;
928: ISDestroy(&dirichletIS);
929: ISDestroy(&neumannIS);
930: return(0);
931: }
933: static PetscErrorCode InitializeDomainData(DomainData *dd)
934: {
936: PetscMPIInt sizes,rank;
937: PetscInt factor;
940: dd->gcomm = PETSC_COMM_WORLD;
941: MPI_Comm_size(dd->gcomm,&sizes);
942: MPI_Comm_rank(dd->gcomm,&rank);
943: /* Get informations from command line */
944: /* Processors/subdomains per dimension */
945: /* Default is 1d problem */
946: dd->npx = sizes;
947: dd->npy = 0;
948: dd->npz = 0;
949: dd->dim = 1;
950: PetscOptionsGetInt(NULL,NULL,"-npx",&dd->npx,NULL);
951: PetscOptionsGetInt(NULL,NULL,"-npy",&dd->npy,NULL);
952: PetscOptionsGetInt(NULL,NULL,"-npz",&dd->npz,NULL);
953: if (dd->npy) dd->dim++;
954: if (dd->npz) dd->dim++;
955: /* Number of elements per dimension */
956: /* Default is one element per subdomain */
957: dd->nex = dd->npx;
958: dd->ney = dd->npy;
959: dd->nez = dd->npz;
960: PetscOptionsGetInt(NULL,NULL,"-nex",&dd->nex,NULL);
961: PetscOptionsGetInt(NULL,NULL,"-ney",&dd->ney,NULL);
962: PetscOptionsGetInt(NULL,NULL,"-nez",&dd->nez,NULL);
963: if (!dd->npy) {
964: dd->ney=0;
965: dd->nez=0;
966: }
967: if (!dd->npz) dd->nez=0;
968: /* Spectral degree */
969: dd->p = 3;
970: PetscOptionsGetInt(NULL,NULL,"-p",&dd->p,NULL);
971: /* pure neumann problem? */
972: dd->pure_neumann = PETSC_FALSE;
973: PetscOptionsGetBool(NULL,NULL,"-pureneumann",&dd->pure_neumann,NULL);
975: /* How to enforce dirichlet boundary conditions (in case pureneumann has not been requested explicitly) */
976: dd->DBC_zerorows = PETSC_FALSE;
978: PetscOptionsGetBool(NULL,NULL,"-usezerorows",&dd->DBC_zerorows,NULL);
979: if (dd->pure_neumann) dd->DBC_zerorows = PETSC_FALSE;
980: dd->scalingfactor = 1.0;
982: factor = 0.0;
983: PetscOptionsGetInt(NULL,NULL,"-jump",&factor,NULL);
984: /* checkerboard pattern */
985: dd->scalingfactor = PetscPowScalar(10.0,(PetscScalar)factor*PetscPowScalar(-1.0,(PetscScalar)rank));
986: /* test data passed in */
987: if (dd->dim==1) {
988: if (sizes!=dd->npx) SETERRQ(dd->gcomm,PETSC_ERR_USER,"Number of mpi procs in 1D must be equal to npx");
989: if (dd->nex<dd->npx) SETERRQ(dd->gcomm,PETSC_ERR_USER,"Number of elements per dim must be greater/equal than number of procs per dim");
990: } else if (dd->dim==2) {
991: if (sizes!=dd->npx*dd->npy) SETERRQ(dd->gcomm,PETSC_ERR_USER,"Number of mpi procs in 2D must be equal to npx*npy");
992: if (dd->nex<dd->npx || dd->ney<dd->npy) SETERRQ(dd->gcomm,PETSC_ERR_USER,"Number of elements per dim must be greater/equal than number of procs per dim");
993: } else {
994: if (sizes!=dd->npx*dd->npy*dd->npz) SETERRQ(dd->gcomm,PETSC_ERR_USER,"Number of mpi procs in 3D must be equal to npx*npy*npz");
995: if (dd->nex<dd->npx || dd->ney<dd->npy || dd->nez<dd->npz) SETERRQ(dd->gcomm,PETSC_ERR_USER,"Number of elements per dim must be greater/equal than number of procs per dim");
996: }
997: return(0);
998: }
1000: int main(int argc,char **args)
1001: {
1002: PetscErrorCode ierr;
1003: DomainData dd;
1004: PetscReal norm,maxeig,mineig;
1005: PetscScalar scalar_value;
1006: PetscInt ndofs,its;
1007: Mat A = NULL,F = NULL;
1008: KSP KSPwithBDDC = NULL,KSPwithFETIDP = NULL;
1009: KSPConvergedReason reason;
1010: Vec exact_solution = NULL,bddc_solution = NULL,bddc_rhs = NULL;
1011: PetscBool testfetidp = PETSC_TRUE;
1013: /* Init PETSc */
1014: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
1015: /* Initialize DomainData */
1016: InitializeDomainData(&dd);
1017: /* Decompose domain */
1018: DomainDecomposition(&dd);
1019: #if DEBUG
1020: printf("Subdomain data\n");
1021: printf("IPS : %d %d %d\n",dd.ipx,dd.ipy,dd.ipz);
1022: printf("NEG : %d %d %d\n",dd.nex,dd.ney,dd.nez);
1023: printf("NEL : %d %d %d\n",dd.nex_l,dd.ney_l,dd.nez_l);
1024: printf("LDO : %d %d %d\n",dd.xm_l,dd.ym_l,dd.zm_l);
1025: printf("SIZES : %d %d %d\n",dd.xm,dd.ym,dd.zm);
1026: printf("STARTS: %d %d %d\n",dd.startx,dd.starty,dd.startz);
1027: #endif
1028: dd.testkspfetidp = PETSC_TRUE;
1029: PetscOptionsGetBool(NULL,NULL,"-testfetidp",&testfetidp,NULL);
1030: PetscOptionsGetBool(NULL,NULL,"-testkspfetidp",&dd.testkspfetidp,NULL);
1031: /* assemble global matrix */
1032: ComputeMatrix(dd,&A);
1033: /* get work vectors */
1034: MatCreateVecs(A,&bddc_solution,NULL);
1035: VecDuplicate(bddc_solution,&bddc_rhs);
1036: VecDuplicate(bddc_solution,&exact_solution);
1037: /* create and customize KSP/PC for BDDC */
1038: ComputeKSPBDDC(dd,A,&KSPwithBDDC);
1039: /* create KSP/PC for FETIDP */
1040: if (testfetidp) {
1041: ComputeKSPFETIDP(dd,KSPwithBDDC,&KSPwithFETIDP);
1042: }
1043: /* create random exact solution */
1044: #if defined(PETSC_USE_COMPLEX)
1045: VecSet(exact_solution,1.0 + PETSC_i);
1046: #else
1047: VecSetRandom(exact_solution,NULL);
1048: #endif
1049: VecShift(exact_solution,-0.5);
1050: VecScale(exact_solution,100.0);
1051: VecGetSize(exact_solution,&ndofs);
1052: if (dd.pure_neumann) {
1053: VecSum(exact_solution,&scalar_value);
1054: scalar_value = -scalar_value/(PetscScalar)ndofs;
1055: VecShift(exact_solution,scalar_value);
1056: }
1057: /* assemble BDDC rhs */
1058: MatMult(A,exact_solution,bddc_rhs);
1059: /* test ksp with BDDC */
1060: KSPSolve(KSPwithBDDC,bddc_rhs,bddc_solution);
1061: KSPGetIterationNumber(KSPwithBDDC,&its);
1062: KSPGetConvergedReason(KSPwithBDDC,&reason);
1063: KSPComputeExtremeSingularValues(KSPwithBDDC,&maxeig,&mineig);
1064: if (dd.pure_neumann) {
1065: VecSum(bddc_solution,&scalar_value);
1066: scalar_value = -scalar_value/(PetscScalar)ndofs;
1067: VecShift(bddc_solution,scalar_value);
1068: }
1069: /* check exact_solution and BDDC solultion */
1070: VecAXPY(bddc_solution,-1.0,exact_solution);
1071: VecNorm(bddc_solution,NORM_INFINITY,&norm);
1072: PetscPrintf(dd.gcomm,"---------------------BDDC stats-------------------------------\n");
1073: PetscPrintf(dd.gcomm,"Number of degrees of freedom : %8D\n",ndofs);
1074: if (reason < 0) {
1075: PetscPrintf(dd.gcomm,"Number of iterations : %8D\n",its);
1076: PetscPrintf(dd.gcomm,"Converged reason : %s\n",KSPConvergedReasons[reason]);
1077: }
1078: if (0.95 <= mineig && mineig <= 1.05) mineig = 1.0;
1079: PetscPrintf(dd.gcomm,"Eigenvalues preconditioned operator : %1.1e %1.1e\n",(double)PetscFloorReal(100.*mineig)/100.,(double)PetscCeilReal(100.*maxeig)/100.);
1080: if (norm > 1.e-1 || reason < 0) {
1081: PetscPrintf(dd.gcomm,"Error betweeen exact and computed solution : %1.2e\n",(double)norm);
1082: }
1083: PetscPrintf(dd.gcomm,"--------------------------------------------------------------\n");
1084: if (testfetidp) {
1085: Vec fetidp_solution_all = NULL,fetidp_solution = NULL,fetidp_rhs = NULL;
1087: VecDuplicate(bddc_solution,&fetidp_solution_all);
1088: if (!dd.testkspfetidp) {
1089: /* assemble fetidp rhs on the space of Lagrange multipliers */
1090: KSPGetOperators(KSPwithFETIDP,&F,NULL);
1091: MatCreateVecs(F,&fetidp_solution,&fetidp_rhs);
1092: PCBDDCMatFETIDPGetRHS(F,bddc_rhs,fetidp_rhs);
1093: VecSet(fetidp_solution,0.0);
1094: /* test ksp with FETIDP */
1095: KSPSolve(KSPwithFETIDP,fetidp_rhs,fetidp_solution);
1096: /* assemble fetidp solution on physical domain */
1097: PCBDDCMatFETIDPGetSolution(F,fetidp_solution,fetidp_solution_all);
1098: } else {
1099: KSP kspF;
1100: KSPSolve(KSPwithFETIDP,bddc_rhs,fetidp_solution_all);
1101: KSPFETIDPGetInnerKSP(KSPwithFETIDP,&kspF);
1102: KSPGetOperators(kspF,&F,NULL);
1103: }
1104: MatGetSize(F,&ndofs,NULL);
1105: KSPGetIterationNumber(KSPwithFETIDP,&its);
1106: KSPGetConvergedReason(KSPwithFETIDP,&reason);
1107: KSPComputeExtremeSingularValues(KSPwithFETIDP,&maxeig,&mineig);
1108: /* check FETIDP sol */
1109: if (dd.pure_neumann) {
1110: VecSum(fetidp_solution_all,&scalar_value);
1111: scalar_value = -scalar_value/(PetscScalar)ndofs;
1112: VecShift(fetidp_solution_all,scalar_value);
1113: }
1114: VecAXPY(fetidp_solution_all,-1.0,exact_solution);
1115: VecNorm(fetidp_solution_all,NORM_INFINITY,&norm);
1116: PetscPrintf(dd.gcomm,"------------------FETI-DP stats-------------------------------\n");
1117: PetscPrintf(dd.gcomm,"Number of degrees of freedom : %8D\n",ndofs);
1118: if (reason < 0) {
1119: PetscPrintf(dd.gcomm,"Number of iterations : %8D\n",its);
1120: PetscPrintf(dd.gcomm,"Converged reason : %D\n",reason);
1121: }
1122: if (0.95 <= mineig && mineig <= 1.05) mineig = 1.0;
1123: PetscPrintf(dd.gcomm,"Eigenvalues preconditioned operator : %1.1e %1.1e\n",(double)PetscFloorReal(100.*mineig)/100.,(double)PetscCeilReal(100.*maxeig)/100.);
1124: if (norm > 1.e-1 || reason < 0) {
1125: PetscPrintf(dd.gcomm,"Error betweeen exact and computed solution : %1.2e\n",(double)norm);
1126: }
1127: PetscPrintf(dd.gcomm,"--------------------------------------------------------------\n");
1128: VecDestroy(&fetidp_solution);
1129: VecDestroy(&fetidp_solution_all);
1130: VecDestroy(&fetidp_rhs);
1131: }
1132: KSPDestroy(&KSPwithFETIDP);
1133: VecDestroy(&exact_solution);
1134: VecDestroy(&bddc_solution);
1135: VecDestroy(&bddc_rhs);
1136: MatDestroy(&A);
1137: KSPDestroy(&KSPwithBDDC);
1138: /* Quit PETSc */
1139: PetscFinalize();
1140: return ierr;
1141: }
1143: /*TEST
1145: testset:
1146: nsize: 4
1147: args: -nex 7 -physical_pc_bddc_coarse_eqs_per_proc 3 -physical_pc_bddc_switch_static
1148: output_file: output/ex59_bddc_fetidp_1.out
1149: test:
1150: suffix: bddc_fetidp_1
1151: test:
1152: requires: viennacl
1153: suffix: bddc_fetidp_1_viennacl
1154: args: -subdomain_mat_type aijviennacl
1155: test:
1156: requires: cuda
1157: suffix: bddc_fetidp_1_cuda
1158: args: -subdomain_mat_type aijcusparse -physical_pc_bddc_dirichlet_pc_factor_mat_solver_type cusparse
1160: testset:
1161: nsize: 4
1162: args: -npx 2 -npy 2 -nex 6 -ney 6 -fluxes_ksp_max_it 10 -physical_ksp_max_it 10
1163: requires: !single
1164: test:
1165: suffix: bddc_fetidp_2
1166: test:
1167: suffix: bddc_fetidp_3
1168: args: -npz 1 -nez 1
1169: test:
1170: suffix: bddc_fetidp_4
1171: args: -npz 1 -nez 1 -physical_pc_bddc_use_change_of_basis -physical_sub_schurs_mat_solver_type petsc -physical_pc_bddc_use_deluxe_scaling -physical_pc_bddc_deluxe_singlemat -fluxes_fetidp_ksp_type cg
1173: testset:
1174: nsize: 8
1175: suffix: bddc_fetidp_approximate
1176: args: -npx 2 -npy 2 -npz 2 -p 2 -nex 8 -ney 7 -nez 9 -physical_ksp_max_it 20 -subdomain_mat_type aij -physical_pc_bddc_switch_static -physical_pc_bddc_dirichlet_approximate -physical_pc_bddc_neumann_approximate -physical_pc_bddc_dirichlet_pc_type gamg -physical_pc_bddc_dirichlet_mg_levels_ksp_max_it 3 -physical_pc_bddc_neumann_pc_type sor -physical_pc_bddc_neumann_approximate_scale -testfetidp 0
1178: testset:
1179: nsize: 4
1180: args: -npx 2 -npy 2 -nex 6 -ney 6 -fluxes_ksp_max_it 10 -physical_ksp_max_it 10 -physical_ksp_view -physical_pc_bddc_levels 1
1181: filter: grep -v "variant HERMITIAN"
1182: requires: !single
1183: test:
1184: suffix: bddc_fetidp_ml_1
1185: args: -physical_pc_bddc_coarsening_ratio 1
1186: test:
1187: suffix: bddc_fetidp_ml_2
1188: args: -physical_pc_bddc_coarsening_ratio 2 -mat_partitioning_type average
1189: test:
1190: suffix: bddc_fetidp_ml_3
1191: args: -physical_pc_bddc_coarsening_ratio 4
1193: testset:
1194: nsize: 9
1195: args: -npx 3 -npy 3 -p 2 -nex 6 -ney 6 -physical_pc_bddc_deluxe_singlemat -physical_sub_schurs_mat_solver_type petsc -physical_pc_bddc_use_deluxe_scaling -physical_pc_bddc_graph_maxcount 1 -physical_pc_bddc_levels 3 -physical_pc_bddc_coarsening_ratio 2 -physical_pc_bddc_coarse_ksp_type gmres
1196: output_file: output/ex59_bddc_fetidp_ml_eqlimit.out
1197: test:
1198: suffix: bddc_fetidp_ml_eqlimit_1
1199: args: -physical_pc_bddc_coarse_eqs_limit 31 -mat_partitioning_type average -physical_pc_bddc_coarse_pc_bddc_graph_maxcount 1
1200: test:
1201: suffix: bddc_fetidp_ml_eqlimit_2
1202: args: -physical_pc_bddc_coarse_eqs_limit 46
1205: TEST*/