Actual source code: test17.c

slepc-3.18.1 2022-11-02
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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);
 37:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 38:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 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=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " 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) EPSKrylovSchurSetPartitions(eps,size);
101:   EPSKrylovSchurGetKSP(eps,&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 = [%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]",nev,ncv,mpd);
129:   EPSKrylovSchurSetDimensions(eps,30,60,60);
130:   EPSKrylovSchurGetDimensions(eps,&nev,&ncv,&mpd);
131:   PetscPrintf(PETSC_COMM_WORLD," ... changed to [%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]\n",nev,ncv,mpd);

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

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

150:   EPSSetFromOptions(eps);

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

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

163:   EPSSolve(eps);
164:   EPSGetDimensions(eps,&nev,NULL,NULL);
165:   EPSGetInterval(eps,&int0,&int1);
166:   PetscPrintf(PETSC_COMM_WORLD," Found %" PetscInt_FMT " eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1);

168:   if (showinertia) {
169:     EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias);
170:     PetscPrintf(PETSC_COMM_WORLD," Used %" PetscInt_FMT " shifts (inertia):\n",k);
171:     for (i=0;i<k;i++) PetscPrintf(PETSC_COMM_WORLD," .. %g (%" PetscInt_FMT ")\n",(double)shifts[i],inertias[i]);
172:     PetscFree(shifts);
173:     PetscFree(inertias);
174:   }

176:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);

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

191:     EPSKrylovSchurGetSubcommMats(eps,&As,&Bs);
192:     MatGetLocalSize(A,&nloc,NULL);
193:     MatGetLocalSize(As,&nlocs,&mlocs);
194:     PetscSynchronizedPrintf(PETSC_COMM_WORLD," Process %d owns %" PetscInt_FMT " rows of the global matrices, and %" PetscInt_FMT " rows in the subcommunicator\n",(int)rank,nloc,nlocs);
195:     PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);

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

221: /*TEST

223:    test:
224:       suffix: 1
225:       nsize: 2
226:       args: -showinertia 0 -log_exclude eps,st,rg,bv,ds
227:       requires: !single

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

235: TEST*/