Actual source code: test17.c

slepc-3.8.3 2018-04-03
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Test interface functions of spectrum-slicing Krylov-Schur.\n\n"
 12:   "This is based on ex12.c. The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 14:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";

 16: #include <slepceps.h>

 18: int main(int argc,char **argv)
 19: {
 20:   Mat            A,B;         /* matrices */
 21:   Mat            As,Bs;       /* matrices distributed in subcommunicators */
 22:   Mat            Au;          /* matrix used to modify A on subcommunicators */
 23:   EPS            eps;         /* eigenproblem solver context */
 24:   ST             st;          /* spectral transformation context */
 25:   KSP            ksp;
 26:   PC             pc;
 27:   Vec            v;
 28:   PetscMPIInt    size,rank;
 29:   PetscInt       N,n=35,m,Istart,Iend,II,nev,ncv,mpd,i,j,k,*inertias,npart,nval,start,nloc,nlocs,mlocs;
 30:   PetscBool      flag,showinertia=PETSC_TRUE,lock,detect;
 31:   PetscReal      int0,int1,*shifts,keep,*subint;
 32:   PetscScalar    eval;
 33:   size_t         count;
 34:   char           vlist[4000];

 37:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 38:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 39:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 41:   PetscOptionsGetBool(NULL,NULL,"-showinertia",&showinertia,NULL);
 42:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 43:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 44:   if (!flag) m=n;
 45:   N = n*m;
 46:   PetscPrintf(PETSC_COMM_WORLD,"\nSpectrum-slicing test, N=%D (%Dx%D grid)\n\n",N,n,m);

 48:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 49:      Compute the matrices that define the eigensystem, Ax=kBx
 50:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 52:   MatCreate(PETSC_COMM_WORLD,&A);
 53:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 54:   MatSetFromOptions(A);
 55:   MatSetUp(A);

 57:   MatCreate(PETSC_COMM_WORLD,&B);
 58:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
 59:   MatSetFromOptions(B);
 60:   MatSetUp(B);

 62:   MatGetOwnershipRange(A,&Istart,&Iend);
 63:   for (II=Istart;II<Iend;II++) {
 64:     i = II/n; j = II-i*n;
 65:     if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
 66:     if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
 67:     if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
 68:     if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
 69:     MatSetValue(A,II,II,4.0,INSERT_VALUES);
 70:     MatSetValue(B,II,II,2.0,INSERT_VALUES);
 71:   }
 72:   if (Istart==0) {
 73:     MatSetValue(B,0,0,6.0,INSERT_VALUES);
 74:     MatSetValue(B,0,1,-1.0,INSERT_VALUES);
 75:     MatSetValue(B,1,0,-1.0,INSERT_VALUES);
 76:     MatSetValue(B,1,1,1.0,INSERT_VALUES);
 77:   }

 79:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 80:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 81:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 82:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 84:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 85:                 Create the eigensolver and set various options
 86:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 88:   EPSCreate(PETSC_COMM_WORLD,&eps);
 89:   EPSSetOperators(eps,A,B);
 90:   EPSSetProblemType(eps,EPS_GHEP);
 91:   EPSSetType(eps,EPSKRYLOVSCHUR);

 93:   /*
 94:      Set interval and other settings for spectrum slicing
 95:   */
 96:   EPSSetWhichEigenpairs(eps,EPS_ALL);
 97:   int0 = 1.1; int1 = 1.3;
 98:   EPSSetInterval(eps,int0,int1);
 99:   EPSGetST(eps,&st);
100:   STSetType(st,STSINVERT);
101:   STGetKSP(st,&ksp);
102:   KSPGetPC(ksp,&pc);
103:   KSPSetType(ksp,KSPPREONLY);
104:   PCSetType(pc,PCCHOLESKY);

106:   /*
107:      Test interface functions of Krylov-Schur solver
108:   */
109:   EPSKrylovSchurGetRestart(eps,&keep);
110:   PetscPrintf(PETSC_COMM_WORLD," Restart parameter before changing = %g",(double)keep);
111:   EPSKrylovSchurSetRestart(eps,0.4);
112:   EPSKrylovSchurGetRestart(eps,&keep);
113:   PetscPrintf(PETSC_COMM_WORLD," ... changed to %g\n",(double)keep);

115:   EPSKrylovSchurGetDetectZeros(eps,&detect);
116:   PetscPrintf(PETSC_COMM_WORLD," Detect zeros before changing = %d",(int)detect);
117:   EPSKrylovSchurSetDetectZeros(eps,PETSC_TRUE);
118:   EPSKrylovSchurGetDetectZeros(eps,&detect);
119:   PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)detect);

121:   EPSKrylovSchurGetLocking(eps,&lock);
122:   PetscPrintf(PETSC_COMM_WORLD," Locking flag before changing = %d",(int)lock);
123:   EPSKrylovSchurSetLocking(eps,PETSC_FALSE);
124:   EPSKrylovSchurGetLocking(eps,&lock);
125:   PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)lock);

127:   EPSKrylovSchurGetDimensions(eps,&nev,&ncv,&mpd);
128:   PetscPrintf(PETSC_COMM_WORLD," Sub-solve dimensions before changing = [%D,%D,%D]",nev,ncv,mpd);
129:   EPSKrylovSchurSetDimensions(eps,30,60,60);
130:   EPSKrylovSchurGetDimensions(eps,&nev,&ncv,&mpd);
131:   PetscPrintf(PETSC_COMM_WORLD," ... changed to [%D,%D,%D]\n",nev,ncv,mpd);

133:   if (size>1) {
134:     EPSKrylovSchurSetPartitions(eps,size);
135:     EPSKrylovSchurGetPartitions(eps,&npart);
136:     PetscPrintf(PETSC_COMM_WORLD," Using %D partitions\n",npart);

138:     PetscMalloc1(npart+1,&subint);
139:     subint[0] = int0;
140:     subint[npart] = int1;
141:     for (i=1;i<npart;i++) subint[i] = int0+i*(int1-int0)/npart;
142:     EPSKrylovSchurSetSubintervals(eps,subint);
143:     PetscFree(subint);
144:     EPSKrylovSchurGetSubintervals(eps,&subint);
145:     PetscPrintf(PETSC_COMM_WORLD," Using sub-interval separations = ");
146:     for (i=1;i<npart;i++) { PetscPrintf(PETSC_COMM_WORLD," %g",(double)subint[i]); }
147:     PetscFree(subint);
148:     PetscPrintf(PETSC_COMM_WORLD,"\n");
149:   }

151:   EPSSetFromOptions(eps);

153:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154:            Compute all eigenvalues in interval and display info
155:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

157:   EPSSetUp(eps);
158:   EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias);
159:   PetscPrintf(PETSC_COMM_WORLD," Inertias after EPSSetUp:\n");
160:   for (i=0;i<k;i++) {
161:     PetscPrintf(PETSC_COMM_WORLD," .. %g (%D)\n",(double)shifts[i],inertias[i]);
162:   }
163:   PetscFree(shifts);
164:   PetscFree(inertias);

166:   EPSSolve(eps);
167:   EPSGetDimensions(eps,&nev,NULL,NULL);
168:   EPSGetInterval(eps,&int0,&int1);
169:   PetscPrintf(PETSC_COMM_WORLD," Found %D eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1);

171:   if (showinertia) {
172:     EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias);
173:     PetscPrintf(PETSC_COMM_WORLD," Used %D shifts (inertia):\n",k);
174:     for (i=0;i<k;i++) {
175:       PetscPrintf(PETSC_COMM_WORLD," .. %g (%D)\n",(double)shifts[i],inertias[i]);
176:     }
177:     PetscFree(shifts);
178:     PetscFree(inertias);
179:   }

181:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);

183:   if (size>1) {
184:     EPSKrylovSchurGetSubcommInfo(eps,&k,&nval,&v);
185:     start = 0;
186:     for (i=0;i<nval;i++) {
187:       EPSKrylovSchurGetSubcommPairs(eps,i,&eval,v);
188:       PetscSNPrintfCount(vlist+start,4000-start,"%g ",&count,(double)PetscRealPart(eval));
189:       start += count;
190:     }
191:     PetscSynchronizedPrintf(PETSC_COMM_WORLD," Process %d has worked in sub-interval %D, containing %D eigenvalues: %s\n",(int)rank,k,nval,vlist);
192:     PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
193:     VecDestroy(&v);

195:     EPSKrylovSchurGetSubcommMats(eps,&As,&Bs);
196:     MatGetLocalSize(A,&nloc,NULL);
197:     MatGetLocalSize(As,&nlocs,&mlocs);
198:     PetscSynchronizedPrintf(PETSC_COMM_WORLD," Process %d owns %D rows of the global matrices, and %D rows in the subcommunicator\n",(int)rank,nloc,nlocs);
199:     PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);

201:     /* modify A on subcommunicators */
202:     MatCreate(PetscObjectComm((PetscObject)As),&Au);
203:     MatSetSizes(Au,nlocs,mlocs,N,N);
204:     MatSetFromOptions(Au);
205:     MatSetUp(Au);
206:     MatGetOwnershipRange(Au,&Istart,&Iend);
207:     for (II=Istart;II<Iend;II++) {
208:       MatSetValue(Au,II,II,0.5,INSERT_VALUES);
209:     }
210:     MatAssemblyBegin(Au,MAT_FINAL_ASSEMBLY);
211:     MatAssemblyEnd(Au,MAT_FINAL_ASSEMBLY);
212:     EPSKrylovSchurUpdateSubcommMats(eps,1.0,-1.0,Au,0.0,0.0,NULL,DIFFERENT_NONZERO_PATTERN,PETSC_TRUE);
213:     MatDestroy(&Au);
214:   }

216:   EPSDestroy(&eps);
217:   MatDestroy(&A);
218:   MatDestroy(&B);
219:   SlepcFinalize();
220:   return ierr;
221: }