Actual source code: ex59.c

petsc-master 2019-06-18
Report Typos and Errors
  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:   PetscMalloc1(p+1,&glldata->zGL);
474:   PetscMemzero(glldata->zGL,(p+1)*sizeof(*glldata->zGL));

476:   glldata->zGL[0]=-1.0;
477:   glldata->zGL[p]= 1.0;
478:   if (p > 1) {
479:     if (p == 2) glldata->zGL[1]=0.0;
480:     else {
481:       PetscMalloc1(p-1,&M);
482:       for (i=0; i<p-1; i++) {
483:         si  = (PetscReal)(i+1.0);
484:         M[i]=0.5*PetscSqrtReal(si*(si+2.0)/((si+0.5)*(si+1.5)));
485:       }
486:       pm1  = p-1;
487:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
488:       PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("N",&pm1,&glldata->zGL[1],M,&x,&pm1,M,&lierr));
489:       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in STERF Lapack routine %d",(int)lierr);
490:       PetscFPTrapPop();
491:       PetscFree(M);
492:     }
493:   }

495:   /* Weights for 1D quadrature */
496:   PetscMalloc1(p+1,&glldata->rhoGL);

498:   glldata->rhoGL[0]=2.0/(PetscScalar)(p*(p+1.0));
499:   glldata->rhoGL[p]=glldata->rhoGL[0];
500:   z2 = -1;                      /* Dummy value to avoid -Wmaybe-initialized */
501:   for (i=1; i<p; i++) {
502:     x  = glldata->zGL[i];
503:     z0 = 1.0;
504:     z1 = x;
505:     for (n=1; n<p; n++) {
506:       z2 = x*z1*(2.0*n+1.0)/(n+1.0)-z0*(PetscScalar)(n/(n+1.0));
507:       z0 = z1;
508:       z1 = z2;
509:     }
510:     glldata->rhoGL[i]=2.0/(p*(p+1.0)*z2*z2);
511:   }

513:   /* Auxiliary mat for laplacian */
514:   PetscMalloc1(p+1,&glldata->A);
515:   PetscMalloc1((p+1)*(p+1),&glldata->A[0]);
516:   for (i=1; i<p+1; i++) glldata->A[i]=glldata->A[i-1]+p+1;

518:   for (j=1; j<p; j++) {
519:     x =glldata->zGL[j];
520:     z0=1.0;
521:     z1=x;
522:     for (n=1; n<p; n++) {
523:       z2=x*z1*(2.0*n+1.0)/(n+1.0)-z0*(PetscScalar)(n/(n+1.0));
524:       z0=z1;
525:       z1=z2;
526:     }
527:     Lpj=z2;
528:     for (r=1; r<p; r++) {
529:       if (r == j) {
530:         glldata->A[j][j]=2.0/(3.0*(1.0-glldata->zGL[j]*glldata->zGL[j])*Lpj*Lpj);
531:       } else {
532:         x  = glldata->zGL[r];
533:         z0 = 1.0;
534:         z1 = x;
535:         for (n=1; n<p; n++) {
536:           z2=x*z1*(2.0*n+1.0)/(n+1.0)-z0*(PetscScalar)(n/(n+1.0));
537:           z0=z1;
538:           z1=z2;
539:         }
540:         Lpr             = z2;
541:         glldata->A[r][j]=4.0/(p*(p+1.0)*Lpj*Lpr*(glldata->zGL[j]-glldata->zGL[r])*(glldata->zGL[j]-glldata->zGL[r]));
542:       }
543:     }
544:   }
545:   for (j=1; j<p+1; j++) {
546:     x  = glldata->zGL[j];
547:     z0 = 1.0;
548:     z1 = x;
549:     for (n=1; n<p; n++) {
550:       z2=x*z1*(2.0*n+1.0)/(n+1.0)-z0*(PetscScalar)(n/(n+1.0));
551:       z0=z1;
552:       z1=z2;
553:     }
554:     Lpj             = z2;
555:     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]));
556:     glldata->A[0][j]=glldata->A[j][0];
557:   }
558:   for (j=0; j<p; j++) {
559:     x  = glldata->zGL[j];
560:     z0 = 1.0;
561:     z1 = x;
562:     for (n=1; n<p; n++) {
563:       z2=x*z1*(2.0*n+1.0)/(n+1.0)-z0*(PetscScalar)(n/(n+1.0));
564:       z0=z1;
565:       z1=z2;
566:     }
567:     Lpj=z2;

569:     glldata->A[p][j]=4.0/(p*(p+1.0)*Lpj*(1.0-glldata->zGL[j])*(1.0-glldata->zGL[j]));
570:     glldata->A[j][p]=glldata->A[p][j];
571:   }
572:   glldata->A[0][0]=0.5+(p*(p+1.0)-2.0)/6.0;
573:   glldata->A[p][p]=glldata->A[0][0];

575:   /* compute element matrix */
576:   xloc = p+1;
577:   yloc = p+1;
578:   zloc = p+1;
579:   if (dd.dim<2) yloc=1;
580:   if (dd.dim<3) zloc=1;
581:   xyloc  = xloc*yloc;
582:   xyzloc = xloc*yloc*zloc;

584:   MatCreate(PETSC_COMM_SELF,&glldata->elem_mat);
585:   MatSetSizes(glldata->elem_mat,xyzloc,xyzloc,xyzloc,xyzloc);
586:   MatSetType(glldata->elem_mat,MATSEQAIJ);
587:   MatSeqAIJSetPreallocation(glldata->elem_mat,xyzloc,NULL); /* overestimated */
588:   MatZeroEntries(glldata->elem_mat);
589:   MatSetOption(glldata->elem_mat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);

591:   for (k=0; k<zloc; k++) {
592:     if (dd.dim>2) rhoGLk=glldata->rhoGL[k];
593:     else rhoGLk=1.0;

595:     for (j=0; j<yloc; j++) {
596:       if (dd.dim>1) rhoGLj=glldata->rhoGL[j];
597:       else rhoGLj=1.0;

599:       for (i=0; i<xloc; i++) {
600:         ii = k*xyloc+j*xloc+i;
601:         s  = k;
602:         r  = j;
603:         for (q=0; q<xloc; q++) {
604:           jj   = s*xyloc+r*xloc+q;
605:           MatSetValue(glldata->elem_mat,jj,ii,glldata->A[i][q]*rhoGLj*rhoGLk,ADD_VALUES);
606:         }
607:         if (dd.dim>1) {
608:           s=k;
609:           q=i;
610:           for (r=0; r<yloc; r++) {
611:             jj   = s*xyloc+r*xloc+q;
612:             MatSetValue(glldata->elem_mat,jj,ii,glldata->A[j][r]*glldata->rhoGL[i]*rhoGLk,ADD_VALUES);
613:           }
614:         }
615:         if (dd.dim>2) {
616:           r=j;
617:           q=i;
618:           for (s=0; s<zloc; s++) {
619:             jj   = s*xyloc+r*xloc+q;
620:             MatSetValue(glldata->elem_mat,jj,ii,glldata->A[k][s]*rhoGLj*glldata->rhoGL[i],ADD_VALUES);
621:           }
622:         }
623:       }
624:     }
625:   }
626:   MatAssemblyBegin(glldata->elem_mat,MAT_FINAL_ASSEMBLY);
627:   MatAssemblyEnd  (glldata->elem_mat,MAT_FINAL_ASSEMBLY);
628: #if DEBUG
629:   {
630:     Vec       lvec,rvec;
631:     PetscReal norm;
632:     MatCreateVecs(glldata->elem_mat,&lvec,&rvec);
633:     VecSet(lvec,1.0);
634:     MatMult(glldata->elem_mat,lvec,rvec);
635:     VecNorm(rvec,NORM_INFINITY,&norm);
636:     VecDestroy(&lvec);
637:     VecDestroy(&rvec);
638:   }
639: #endif
640:   return(0);
641: }

643: static PetscErrorCode DomainDecomposition(DomainData *dd)
644: {
645:   PetscMPIInt rank;
646:   PetscInt    i,j,k;

649:   /* Subdomain index in cartesian coordinates */
650:   MPI_Comm_rank(dd->gcomm,&rank);
651:   dd->ipx = rank%dd->npx;
652:   if (dd->dim>1) dd->ipz = rank/(dd->npx*dd->npy);
653:   else dd->ipz = 0;

655:   dd->ipy = rank/dd->npx-dd->ipz*dd->npy;
656:   /* number of local elements */
657:   dd->nex_l = dd->nex/dd->npx;
658:   if (dd->ipx < dd->nex%dd->npx) dd->nex_l++;
659:   if (dd->dim>1) {
660:     dd->ney_l = dd->ney/dd->npy;
661:     if (dd->ipy < dd->ney%dd->npy) dd->ney_l++;
662:   } else {
663:     dd->ney_l = 0;
664:   }
665:   if (dd->dim>2) {
666:     dd->nez_l = dd->nez/dd->npz;
667:     if (dd->ipz < dd->nez%dd->npz) dd->nez_l++;
668:   } else {
669:     dd->nez_l = 0;
670:   }
671:   /* local and global number of dofs */
672:   dd->xm_l = dd->nex_l*dd->p+1;
673:   dd->xm   = dd->nex*dd->p+1;
674:   dd->ym_l = dd->ney_l*dd->p+1;
675:   dd->ym   = dd->ney*dd->p+1;
676:   dd->zm_l = dd->nez_l*dd->p+1;
677:   dd->zm   = dd->nez*dd->p+1;
678:   if (!dd->pure_neumann) {
679:     if (!dd->DBC_zerorows && !dd->ipx) dd->xm_l--;
680:     if (!dd->DBC_zerorows) dd->xm--;
681:   }

683:   /* starting global index for local dofs (simple lexicographic order) */
684:   dd->startx = 0;
685:   j          = dd->nex/dd->npx;
686:   for (i=0; i<dd->ipx; i++) {
687:     k = j;
688:     if (i<dd->nex%dd->npx) k++;
689:     dd->startx = dd->startx+k*dd->p;
690:   }
691:   if (!dd->pure_neumann && !dd->DBC_zerorows && dd->ipx) dd->startx--;

693:   dd->starty = 0;
694:   if (dd->dim > 1) {
695:     j = dd->ney/dd->npy;
696:     for (i=0; i<dd->ipy; i++) {
697:       k = j;
698:       if (i<dd->ney%dd->npy) k++;
699:       dd->starty = dd->starty+k*dd->p;
700:     }
701:   }
702:   dd->startz = 0;
703:   if (dd->dim > 2) {
704:     j = dd->nez/dd->npz;
705:     for (i=0; i<dd->ipz; i++) {
706:       k = j;
707:       if (i<dd->nez%dd->npz) k++;
708:       dd->startz = dd->startz+k*dd->p;
709:     }
710:   }
711:   return(0);
712: }
713: static PetscErrorCode ComputeMatrix(DomainData dd, Mat *A)
714: {
715:   PetscErrorCode         ierr;
716:   GLLData                gll;
717:   Mat                    local_mat  =0,temp_A=0;
718:   ISLocalToGlobalMapping matis_map  =0;
719:   IS                     dirichletIS=0;

722:   /* Compute some stuff of Gauss-Legendre-Lobatto quadrature rule */
723:   GLLStuffs(dd,&gll);
724:   /* Compute matrix of subdomain Neumann problem */
725:   ComputeSubdomainMatrix(dd,gll,&local_mat);
726:   /* Compute global mapping of local dofs */
727:   ComputeMapping(dd,&matis_map);
728:   /* Create MATIS object needed by BDDC */
729:   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);
730:   /* Set local subdomain matrices into MATIS object */
731:   MatScale(local_mat,dd.scalingfactor);
732:   MatISSetLocalMat(temp_A,local_mat);
733:   /* Call assembly functions */
734:   MatAssemblyBegin(temp_A,MAT_FINAL_ASSEMBLY);
735:   MatAssemblyEnd(temp_A,MAT_FINAL_ASSEMBLY);

737:   if (dd.DBC_zerorows) {
738:     PetscInt dirsize;

740:     ComputeSpecialBoundaryIndices(dd,&dirichletIS,NULL);
741:     MatSetOption(local_mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
742:     MatZeroRowsColumnsLocalIS(temp_A,dirichletIS,1.0,NULL,NULL);
743:     ISGetLocalSize(dirichletIS,&dirsize);
744:     ISDestroy(&dirichletIS);
745:   }

747:   /* giving hints to local and global matrices could be useful for the BDDC */
748:   MatSetOption(local_mat,MAT_SYMMETRIC,PETSC_TRUE);
749:   MatSetOption(local_mat,MAT_SPD,PETSC_TRUE);
750: #if DEBUG
751:   {
752:     Vec       lvec,rvec;
753:     PetscReal norm;
754:     MatCreateVecs(temp_A,&lvec,&rvec);
755:     VecSet(lvec,1.0);
756:     MatMult(temp_A,lvec,rvec);
757:     VecNorm(rvec,NORM_INFINITY,&norm);
758:     VecDestroy(&lvec);
759:     VecDestroy(&rvec);
760:   }
761: #endif
762:   /* free allocated workspace */
763:   PetscFree(gll.zGL);
764:   PetscFree(gll.rhoGL);
765:   PetscFree(gll.A[0]);
766:   PetscFree(gll.A);
767:   MatDestroy(&gll.elem_mat);
768:   MatDestroy(&local_mat);
769:   ISLocalToGlobalMappingDestroy(&matis_map);
770:   /* give back the pointer to te MATIS object */
771:   *A = temp_A;
772:   return(0);
773: }

775: static PetscErrorCode ComputeKSPFETIDP(DomainData dd, KSP ksp_bddc, KSP *ksp_fetidp)
776: {
778:   KSP            temp_ksp;
779:   PC             pc;

782:   KSPGetPC(ksp_bddc,&pc);
783:   if (!dd.testkspfetidp) {
784:     PC  D;
785:     Mat F;

787:     PCBDDCCreateFETIDPOperators(pc,PETSC_TRUE,NULL,&F,&D);
788:     KSPCreate(PetscObjectComm((PetscObject)F),&temp_ksp);
789:     KSPSetOperators(temp_ksp,F,F);
790:     KSPSetType(temp_ksp,KSPCG);
791:     KSPSetPC(temp_ksp,D);
792:     KSPSetComputeSingularValues(temp_ksp,PETSC_TRUE);
793:     KSPSetOptionsPrefix(temp_ksp,"fluxes_");
794:     KSPSetFromOptions(temp_ksp);
795:     KSPSetUp(temp_ksp);
796:     MatDestroy(&F);
797:     PCDestroy(&D);
798:   } else {
799:     Mat A,Ap;

801:     KSPCreate(PetscObjectComm((PetscObject)ksp_bddc),&temp_ksp);
802:     KSPGetOperators(ksp_bddc,&A,&Ap);
803:     KSPSetOperators(temp_ksp,A,Ap);
804:     KSPSetOptionsPrefix(temp_ksp,"fluxes_");
805:     KSPSetType(temp_ksp,KSPFETIDP);
806:     KSPFETIDPSetInnerBDDC(temp_ksp,pc);
807:     KSPSetComputeSingularValues(temp_ksp,PETSC_TRUE);
808:     KSPSetFromOptions(temp_ksp);
809:     KSPSetUp(temp_ksp);
810:   }
811:   *ksp_fetidp = temp_ksp;
812:   return(0);
813: }


816: static PetscErrorCode ComputeKSPBDDC(DomainData dd,Mat A,KSP *ksp)
817: {
819:   KSP            temp_ksp;
820:   PC             pc;
821:   IS             primals,dirichletIS=0,neumannIS=0,*bddc_dofs_splitting;
822:   PetscInt       vidx[8],localsize,*xadj=NULL,*adjncy=NULL;
823:   MatNullSpace   near_null_space;

826:   KSPCreate(dd.gcomm,&temp_ksp);
827:   KSPSetOperators(temp_ksp,A,A);
828:   KSPSetType(temp_ksp,KSPCG);
829:   KSPGetPC(temp_ksp,&pc);
830:   PCSetType(pc,PCBDDC);

832:   localsize = dd.xm_l*dd.ym_l*dd.zm_l;

834:   /* BDDC customization */

836:   /* jumping coefficients case */
837:   PCISSetSubdomainScalingFactor(pc,dd.scalingfactor);

839:   /* Dofs splitting
840:      Simple stride-1 IS
841:      It is not needed since, by default, PCBDDC assumes a stride-1 split */
842:   PetscMalloc1(1,&bddc_dofs_splitting);
843: #if 1
844:   ISCreateStride(PETSC_COMM_WORLD,localsize,0,1,&bddc_dofs_splitting[0]);
845:   PCBDDCSetDofsSplittingLocal(pc,1,bddc_dofs_splitting);
846: #else
847:   /* examples for global ordering */

849:   /* each process lists the nodes it owns */
850:   PetscInt sr,er;
851:   MatGetOwnershipRange(A,&sr,&er);
852:   ISCreateStride(PETSC_COMM_WORLD,er-sr,sr,1,&bddc_dofs_splitting[0]);
853:   PCBDDCSetDofsSplitting(pc,1,bddc_dofs_splitting);
854:   /* Split can be passed in a more general way since any process can list any node */
855: #endif
856:   ISDestroy(&bddc_dofs_splitting[0]);
857:   PetscFree(bddc_dofs_splitting);

859:   /* Primal constraints implemented by using a near null space attached to A -> now it passes in only the constants
860:     (which in practice is not needed since by default PCBDDC builds the primal space using constants for quadrature formulas */
861: #if 0
862:   Vec vecs[2];
863:   PetscRandom rctx;
864:   MatCreateVecs(A,&vecs[0],&vecs[1]);
865:   PetscRandomCreate(dd.gcomm,&rctx);
866:   VecSetRandom(vecs[0],rctx);
867:   VecSetRandom(vecs[1],rctx);
868:   MatNullSpaceCreate(dd.gcomm,PETSC_TRUE,2,vecs,&near_null_space);
869:   VecDestroy(&vecs[0]);
870:   VecDestroy(&vecs[1]);
871:   PetscRandomDestroy(&rctx);
872: #else
873:   MatNullSpaceCreate(dd.gcomm,PETSC_TRUE,0,NULL,&near_null_space);
874: #endif
875:   MatSetNearNullSpace(A,near_null_space);
876:   MatNullSpaceDestroy(&near_null_space);

878:   /* CSR graph of subdomain dofs */
879:   BuildCSRGraph(dd,&xadj,&adjncy);
880:   PCBDDCSetLocalAdjacencyGraph(pc,localsize,xadj,adjncy,PETSC_OWN_POINTER);

882:   /* Prescribe user-defined primal vertices: in this case we use the 8 corners in 3D (for 2D and 1D, the indices are duplicated) */
883:   vidx[0] = 0*dd.xm_l+0;
884:   vidx[1] = 0*dd.xm_l+dd.xm_l-1;
885:   vidx[2] = (dd.ym_l-1)*dd.xm_l+0;
886:   vidx[3] = (dd.ym_l-1)*dd.xm_l+dd.xm_l-1;
887:   vidx[4] = (dd.zm_l-1)*dd.xm_l*dd.ym_l+0*dd.xm_l+0;
888:   vidx[5] = (dd.zm_l-1)*dd.xm_l*dd.ym_l+0*dd.xm_l+dd.xm_l-1;
889:   vidx[6] = (dd.zm_l-1)*dd.xm_l*dd.ym_l+(dd.ym_l-1)*dd.xm_l+0;
890:   vidx[7] = (dd.zm_l-1)*dd.xm_l*dd.ym_l+(dd.ym_l-1)*dd.xm_l+dd.xm_l-1;
891:   ISCreateGeneral(dd.gcomm,8,vidx,PETSC_COPY_VALUES,&primals);
892:   PCBDDCSetPrimalVerticesLocalIS(pc,primals);
893:   ISDestroy(&primals);

895:   /* Neumann/Dirichlet indices on the global boundary */
896:   if (dd.DBC_zerorows) {
897:     /* Only in case you eliminate some rows matrix with zerorows function, you need to set dirichlet indices into PCBDDC data */
898:     ComputeSpecialBoundaryIndices(dd,&dirichletIS,&neumannIS);
899:     PCBDDCSetNeumannBoundariesLocal(pc,neumannIS);
900:     PCBDDCSetDirichletBoundariesLocal(pc,dirichletIS);
901:   } else {
902:     if (dd.pure_neumann) {
903:       /* In such a case, all interface nodes lying on the global boundary are neumann nodes */
904:       ComputeSpecialBoundaryIndices(dd,NULL,&neumannIS);
905:       PCBDDCSetNeumannBoundariesLocal(pc,neumannIS);
906:     } else {
907:       /* It is wrong setting dirichlet indices without having zeroed the corresponding rows in the global matrix */
908:       /* But we can still compute them since nodes near the dirichlet boundaries does not need to be defined as neumann nodes */
909:       ComputeSpecialBoundaryIndices(dd,&dirichletIS,&neumannIS);
910:       PCBDDCSetNeumannBoundariesLocal(pc,neumannIS);
911:     }
912:   }

914:   /* Pass local null space information to local matrices (needed when using approximate local solvers) */
915:   if (dd.ipx || dd.pure_neumann) {
916:     MatNullSpace nsp;
917:     Mat          local_mat;

919:     MatISGetLocalMat(A,&local_mat);
920:     MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nsp);
921:     MatSetNullSpace(local_mat,nsp);
922:     MatNullSpaceDestroy(&nsp);
923:   }
924:   KSPSetComputeSingularValues(temp_ksp,PETSC_TRUE);
925:   KSPSetOptionsPrefix(temp_ksp,"physical_");
926:   KSPSetFromOptions(temp_ksp);
927:   KSPSetUp(temp_ksp);
928:   *ksp = temp_ksp;
929:   ISDestroy(&dirichletIS);
930:   ISDestroy(&neumannIS);
931:   return(0);
932: }

934: static PetscErrorCode InitializeDomainData(DomainData *dd)
935: {
937:   PetscMPIInt    sizes,rank;
938:   PetscInt       factor;

941:   dd->gcomm = PETSC_COMM_WORLD;
942:   MPI_Comm_size(dd->gcomm,&sizes);
943:   MPI_Comm_rank(dd->gcomm,&rank);
944:   /* Get informations from command line */
945:   /* Processors/subdomains per dimension */
946:   /* Default is 1d problem */
947:   dd->npx = sizes;
948:   dd->npy = 0;
949:   dd->npz = 0;
950:   dd->dim = 1;
951:   PetscOptionsGetInt(NULL,NULL,"-npx",&dd->npx,NULL);
952:   PetscOptionsGetInt(NULL,NULL,"-npy",&dd->npy,NULL);
953:   PetscOptionsGetInt(NULL,NULL,"-npz",&dd->npz,NULL);
954:   if (dd->npy) dd->dim++;
955:   if (dd->npz) dd->dim++;
956:   /* Number of elements per dimension */
957:   /* Default is one element per subdomain */
958:   dd->nex = dd->npx;
959:   dd->ney = dd->npy;
960:   dd->nez = dd->npz;
961:   PetscOptionsGetInt(NULL,NULL,"-nex",&dd->nex,NULL);
962:   PetscOptionsGetInt(NULL,NULL,"-ney",&dd->ney,NULL);
963:   PetscOptionsGetInt(NULL,NULL,"-nez",&dd->nez,NULL);
964:   if (!dd->npy) {
965:     dd->ney=0;
966:     dd->nez=0;
967:   }
968:   if (!dd->npz) dd->nez=0;
969:   /* Spectral degree */
970:   dd->p = 3;
971:   PetscOptionsGetInt(NULL,NULL,"-p",&dd->p,NULL);
972:   /* pure neumann problem? */
973:   dd->pure_neumann = PETSC_FALSE;
974:   PetscOptionsGetBool(NULL,NULL,"-pureneumann",&dd->pure_neumann,NULL);

976:   /* How to enforce dirichlet boundary conditions (in case pureneumann has not been requested explicitly) */
977:   dd->DBC_zerorows = PETSC_FALSE;

979:   PetscOptionsGetBool(NULL,NULL,"-usezerorows",&dd->DBC_zerorows,NULL);
980:   if (dd->pure_neumann) dd->DBC_zerorows = PETSC_FALSE;
981:   dd->scalingfactor = 1.0;

983:   factor = 0.0;
984:   PetscOptionsGetInt(NULL,NULL,"-jump",&factor,NULL);
985:   /* checkerboard pattern */
986:   dd->scalingfactor = PetscPowScalar(10.0,(PetscScalar)factor*PetscPowScalar(-1.0,(PetscScalar)rank));
987:   /* test data passed in */
988:   if (dd->dim==1) {
989:     if (sizes!=dd->npx) SETERRQ(dd->gcomm,PETSC_ERR_USER,"Number of mpi procs in 1D must be equal to npx");
990:     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");
991:   } else if (dd->dim==2) {
992:     if (sizes!=dd->npx*dd->npy) SETERRQ(dd->gcomm,PETSC_ERR_USER,"Number of mpi procs in 2D must be equal to npx*npy");
993:     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");
994:   } else {
995:     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");
996:     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");
997:   }
998:   return(0);
999: }

1001: int main(int argc,char **args)
1002: {
1003:   PetscErrorCode     ierr;
1004:   DomainData         dd;
1005:   PetscReal          norm,maxeig,mineig;
1006:   PetscScalar        scalar_value;
1007:   PetscInt           ndofs,its;
1008:   Mat                A = NULL,F = NULL;
1009:   KSP                KSPwithBDDC = NULL,KSPwithFETIDP = NULL;
1010:   KSPConvergedReason reason;
1011:   Vec                exact_solution = NULL,bddc_solution = NULL,bddc_rhs = NULL;
1012:   PetscBool          testfetidp = PETSC_TRUE;

1014:   /* Init PETSc */
1015:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
1016:   /* Initialize DomainData */
1017:   InitializeDomainData(&dd);
1018:   /* Decompose domain */
1019:   DomainDecomposition(&dd);
1020: #if DEBUG
1021:   printf("Subdomain data\n");
1022:   printf("IPS   : %d %d %d\n",dd.ipx,dd.ipy,dd.ipz);
1023:   printf("NEG   : %d %d %d\n",dd.nex,dd.ney,dd.nez);
1024:   printf("NEL   : %d %d %d\n",dd.nex_l,dd.ney_l,dd.nez_l);
1025:   printf("LDO   : %d %d %d\n",dd.xm_l,dd.ym_l,dd.zm_l);
1026:   printf("SIZES : %d %d %d\n",dd.xm,dd.ym,dd.zm);
1027:   printf("STARTS: %d %d %d\n",dd.startx,dd.starty,dd.startz);
1028: #endif
1029:   dd.testkspfetidp = PETSC_TRUE;
1030:   PetscOptionsGetBool(NULL,NULL,"-testfetidp",&testfetidp,NULL);
1031:   PetscOptionsGetBool(NULL,NULL,"-testkspfetidp",&dd.testkspfetidp,NULL);
1032:   /* assemble global matrix */
1033:   ComputeMatrix(dd,&A);
1034:   /* get work vectors */
1035:   MatCreateVecs(A,&bddc_solution,NULL);
1036:   VecDuplicate(bddc_solution,&bddc_rhs);
1037:   VecDuplicate(bddc_solution,&exact_solution);
1038:   /* create and customize KSP/PC for BDDC */
1039:   ComputeKSPBDDC(dd,A,&KSPwithBDDC);
1040:   /* create KSP/PC for FETIDP */
1041:   if (testfetidp) {
1042:     ComputeKSPFETIDP(dd,KSPwithBDDC,&KSPwithFETIDP);
1043:   }
1044:   /* create random exact solution */
1045: #if defined(PETSC_USE_COMPLEX)
1046:   VecSet(exact_solution,1.0 + PETSC_i);
1047: #else
1048:   VecSetRandom(exact_solution,NULL);
1049: #endif
1050:   VecShift(exact_solution,-0.5);
1051:   VecScale(exact_solution,100.0);
1052:   VecGetSize(exact_solution,&ndofs);
1053:   if (dd.pure_neumann) {
1054:     VecSum(exact_solution,&scalar_value);
1055:     scalar_value = -scalar_value/(PetscScalar)ndofs;
1056:     VecShift(exact_solution,scalar_value);
1057:   }
1058:   /* assemble BDDC rhs */
1059:   MatMult(A,exact_solution,bddc_rhs);
1060:   /* test ksp with BDDC */
1061:   KSPSolve(KSPwithBDDC,bddc_rhs,bddc_solution);
1062:   KSPGetIterationNumber(KSPwithBDDC,&its);
1063:   KSPGetConvergedReason(KSPwithBDDC,&reason);
1064:   KSPComputeExtremeSingularValues(KSPwithBDDC,&maxeig,&mineig);
1065:   if (dd.pure_neumann) {
1066:     VecSum(bddc_solution,&scalar_value);
1067:     scalar_value = -scalar_value/(PetscScalar)ndofs;
1068:     VecShift(bddc_solution,scalar_value);
1069:   }
1070:   /* check exact_solution and BDDC solultion */
1071:   VecAXPY(bddc_solution,-1.0,exact_solution);
1072:   VecNorm(bddc_solution,NORM_INFINITY,&norm);
1073:   PetscPrintf(dd.gcomm,"---------------------BDDC stats-------------------------------\n");
1074:   PetscPrintf(dd.gcomm,"Number of degrees of freedom               : %8D\n",ndofs);
1075:   if (reason < 0) {
1076:     PetscPrintf(dd.gcomm,"Number of iterations                       : %8D\n",its);
1077:     PetscPrintf(dd.gcomm,"Converged reason                           : %s\n",KSPConvergedReasons[reason]);
1078:   }
1079:   if (0.95 <= mineig && mineig <= 1.05) mineig = 1.0;
1080:   PetscPrintf(dd.gcomm,"Eigenvalues preconditioned operator        : %1.1e %1.1e\n",(double)PetscFloorReal(100.*mineig)/100.,(double)PetscCeilReal(100.*maxeig)/100.);
1081:   if (norm > 1.e-1 || reason < 0) {
1082:     PetscPrintf(dd.gcomm,"Error betweeen exact and computed solution : %1.2e\n",(double)norm);
1083:   }
1084:   PetscPrintf(dd.gcomm,"--------------------------------------------------------------\n");
1085:   if (testfetidp) {
1086:     Vec fetidp_solution_all = NULL,fetidp_solution = NULL,fetidp_rhs = NULL;

1088:     VecDuplicate(bddc_solution,&fetidp_solution_all);
1089:     if (!dd.testkspfetidp) {
1090:       /* assemble fetidp rhs on the space of Lagrange multipliers */
1091:       KSPGetOperators(KSPwithFETIDP,&F,NULL);
1092:       MatCreateVecs(F,&fetidp_solution,&fetidp_rhs);
1093:       PCBDDCMatFETIDPGetRHS(F,bddc_rhs,fetidp_rhs);
1094:       VecSet(fetidp_solution,0.0);
1095:       /* test ksp with FETIDP */
1096:       KSPSolve(KSPwithFETIDP,fetidp_rhs,fetidp_solution);
1097:       /* assemble fetidp solution on physical domain */
1098:       PCBDDCMatFETIDPGetSolution(F,fetidp_solution,fetidp_solution_all);
1099:     } else {
1100:       KSP kspF;
1101:       KSPSolve(KSPwithFETIDP,bddc_rhs,fetidp_solution_all);
1102:       KSPFETIDPGetInnerKSP(KSPwithFETIDP,&kspF);
1103:       KSPGetOperators(kspF,&F,NULL);
1104:     }
1105:     MatGetSize(F,&ndofs,NULL);
1106:     KSPGetIterationNumber(KSPwithFETIDP,&its);
1107:     KSPGetConvergedReason(KSPwithFETIDP,&reason);
1108:     KSPComputeExtremeSingularValues(KSPwithFETIDP,&maxeig,&mineig);
1109:     /* check FETIDP sol */
1110:     if (dd.pure_neumann) {
1111:       VecSum(fetidp_solution_all,&scalar_value);
1112:       scalar_value = -scalar_value/(PetscScalar)ndofs;
1113:       VecShift(fetidp_solution_all,scalar_value);
1114:     }
1115:     VecAXPY(fetidp_solution_all,-1.0,exact_solution);
1116:     VecNorm(fetidp_solution_all,NORM_INFINITY,&norm);
1117:     PetscPrintf(dd.gcomm,"------------------FETI-DP stats-------------------------------\n");
1118:     PetscPrintf(dd.gcomm,"Number of degrees of freedom               : %8D\n",ndofs);
1119:     if (reason < 0) {
1120:       PetscPrintf(dd.gcomm,"Number of iterations                       : %8D\n",its);
1121:       PetscPrintf(dd.gcomm,"Converged reason                           : %D\n",reason);
1122:     }
1123:     if (0.95 <= mineig && mineig <= 1.05) mineig = 1.0;
1124:     PetscPrintf(dd.gcomm,"Eigenvalues preconditioned operator        : %1.1e %1.1e\n",(double)PetscFloorReal(100.*mineig)/100.,(double)PetscCeilReal(100.*maxeig)/100.);
1125:     if (norm > 1.e-1 || reason < 0) {
1126:       PetscPrintf(dd.gcomm,"Error betweeen exact and computed solution : %1.2e\n",(double)norm);
1127:     }
1128:     PetscPrintf(dd.gcomm,"--------------------------------------------------------------\n");
1129:     VecDestroy(&fetidp_solution);
1130:     VecDestroy(&fetidp_solution_all);
1131:     VecDestroy(&fetidp_rhs);
1132:   }
1133:   KSPDestroy(&KSPwithFETIDP);
1134:   VecDestroy(&exact_solution);
1135:   VecDestroy(&bddc_solution);
1136:   VecDestroy(&bddc_rhs);
1137:   MatDestroy(&A);
1138:   KSPDestroy(&KSPwithBDDC);
1139:   /* Quit PETSc */
1140:   PetscFinalize();
1141:   return ierr;
1142: }

1144: /*TEST

1146:  testset:
1147:    nsize: 4
1148:    args: -nex 7 -physical_pc_bddc_coarse_eqs_per_proc 3 -physical_pc_bddc_switch_static
1149:    output_file: output/ex59_bddc_fetidp_1.out
1150:    test:
1151:      suffix: bddc_fetidp_1
1152:    test:
1153:      requires: viennacl
1154:      suffix: bddc_fetidp_1_viennacl
1155:      args: -subdomain_mat_type aijviennacl
1156:    test:
1157:      requires: cuda
1158:      suffix: bddc_fetidp_1_cuda
1159:      args: -subdomain_mat_type aijcusparse -physical_pc_bddc_dirichlet_pc_factor_mat_solver_type cusparse

1161:  testset:
1162:    nsize: 4
1163:    args: -npx 2 -npy 2 -nex 6 -ney 6 -fluxes_ksp_max_it 10 -physical_ksp_max_it 10
1164:    requires: !single
1165:    test:
1166:      suffix: bddc_fetidp_2
1167:    test:
1168:      suffix: bddc_fetidp_3
1169:      args: -npz 1 -nez 1
1170:    test:
1171:      suffix: bddc_fetidp_4
1172:      args: -npz 1 -nez 1 -physical_pc_bddc_use_change_of_basis -physical_pc_bddc_use_deluxe_scaling -physical_pc_bddc_deluxe_singlemat -fluxes_fetidp_ksp_type cg

1174:  testset:
1175:    nsize: 8
1176:    suffix: bddc_fetidp_approximate
1177:    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

1179:  testset:
1180:    nsize: 4
1181:    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
1182:    filter: grep -v "variant HERMITIAN"
1183:    requires: !single
1184:    test:
1185:      suffix: bddc_fetidp_ml_1
1186:      args: -physical_pc_bddc_coarsening_ratio 1
1187:    test:
1188:      suffix: bddc_fetidp_ml_2
1189:      args: -physical_pc_bddc_coarsening_ratio 2 -mat_partitioning_type average
1190:    test:
1191:      suffix: bddc_fetidp_ml_3
1192:      args: -physical_pc_bddc_coarsening_ratio 4

1194:  testset:
1195:    nsize: 9
1196:    args: -npx 3 -npy 3 -p 2 -nex 6 -ney 6 -physical_pc_bddc_deluxe_singlemat -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
1197:    output_file: output/ex59_bddc_fetidp_ml_eqlimit.out
1198:    test:
1199:      suffix: bddc_fetidp_ml_eqlimit_1
1200:      args: -physical_pc_bddc_coarse_eqs_limit 31 -mat_partitioning_type average -physical_pc_bddc_coarse_pc_bddc_graph_maxcount 1
1201:    test:
1202:      suffix: bddc_fetidp_ml_eqlimit_2
1203:      args: -physical_pc_bddc_coarse_eqs_limit 46


1206: TEST*/