Actual source code: test17.c

slepc-3.22.1 2024-10-28
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];

 35:   PetscFunctionBeginUser;
 36:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 37:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
 38:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));

 40:   PetscCall(PetscOptionsGetBool(NULL,NULL,"-showinertia",&showinertia,NULL));
 41:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 42:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
 43:   if (!flag) m=n;
 44:   N = n*m;
 45:   PetscCall(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:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 52:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
 53:   PetscCall(MatSetFromOptions(A));

 55:   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
 56:   PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N));
 57:   PetscCall(MatSetFromOptions(B));

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

 76:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 77:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
 78:   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
 79:   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));

 81:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 82:                 Create the eigensolver and set various options
 83:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 85:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 86:   PetscCall(EPSSetOperators(eps,A,B));
 87:   PetscCall(EPSSetProblemType(eps,EPS_GHEP));
 88:   PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));

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

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

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

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

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

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

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

148:   PetscCall(EPSSetFromOptions(eps));

150:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151:            Compute all eigenvalues in interval and display info
152:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

174:   PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));

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

189:     PetscCall(EPSKrylovSchurGetSubcommMats(eps,&As,&Bs));
190:     PetscCall(MatGetLocalSize(A,&nloc,NULL));
191:     PetscCall(MatGetLocalSize(As,&nlocs,&mlocs));
192:     PetscCall(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));
193:     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));

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

218: /*TEST

220:    test:
221:       suffix: 1
222:       nsize: 2
223:       args: -showinertia 0 -log_exclude eps,st,rg,bv,ds
224:       requires: !single

226:    test:
227:       suffix: 2
228:       nsize: 1
229:       args: -showinertia 0 -log_exclude eps,st,rg,bv,ds
230:       requires: !single

232: TEST*/