Actual source code: test17.c

slepc-3.15.2 2021-09-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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,nloc,nlocs,mlocs;
 30:   PetscBool      flag,showinertia=PETSC_TRUE,lock,detect;
 31:   PetscReal      int0,int1,*shifts,keep,*subint,*evals;
 32:   PetscScalar    lambda;
 33:   char           vlist[4000];

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

152:   EPSSetFromOptions(eps);

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

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

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

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

182:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);

184:   if (size>1) {
185:     EPSKrylovSchurGetSubcommInfo(eps,&k,&nval,&v);
186:     PetscMalloc1(nval,&evals);
187:     for (i=0;i<nval;i++) {
188:       EPSKrylovSchurGetSubcommPairs(eps,i,&lambda,v);
189:       evals[i] = PetscRealPart(lambda);
190:     }
191:     PetscFormatRealArray(vlist,sizeof(vlist),"%f",nval,evals);
192:     PetscSynchronizedPrintf(PETSC_COMM_WORLD," Process %d has worked in sub-interval %D, containing %D eigenvalues: %s\n",(int)rank,k,nval,vlist);
193:     PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
194:     VecDestroy(&v);
195:     PetscFree(evals);

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

203:     /* modify A on subcommunicators */
204:     MatCreate(PetscObjectComm((PetscObject)As),&Au);
205:     MatSetSizes(Au,nlocs,mlocs,N,N);
206:     MatSetFromOptions(Au);
207:     MatSetUp(Au);
208:     MatGetOwnershipRange(Au,&Istart,&Iend);
209:     for (II=Istart;II<Iend;II++) {
210:       MatSetValue(Au,II,II,0.5,INSERT_VALUES);
211:     }
212:     MatAssemblyBegin(Au,MAT_FINAL_ASSEMBLY);
213:     MatAssemblyEnd(Au,MAT_FINAL_ASSEMBLY);
214:     PetscPrintf(PETSC_COMM_WORLD," Updating internal matrices\n");
215:     EPSKrylovSchurUpdateSubcommMats(eps,1.1,-5.0,Au,1.0,0.0,NULL,DIFFERENT_NONZERO_PATTERN,PETSC_TRUE);
216:     MatDestroy(&Au);
217:     EPSSolve(eps);
218:     EPSGetDimensions(eps,&nev,NULL,NULL);
219:     EPSGetInterval(eps,&int0,&int1);
220:     PetscPrintf(PETSC_COMM_WORLD," After update, found %D eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1);
221:   }
222:   EPSDestroy(&eps);
223:   MatDestroy(&A);
224:   MatDestroy(&B);
225:   SlepcFinalize();
226:   return ierr;
227: }

229: /*TEST

231:    test:
232:       suffix: 1
233:       nsize: 2
234:       args: -showinertia 0 -log_exclude eps,st,rg,bv,ds
235:       requires: !single

237:    test:
238:       suffix: 2
239:       nsize: 1
240:       args: -showinertia 0 -log_exclude eps,st,rg,bv,ds
241:       requires: !single

243: TEST*/