Actual source code: test17.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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);
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=%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: STGetKSP(st,&ksp);
101: KSPGetPC(ksp,&pc);
102: KSPSetType(ksp,KSPPREONLY);
103: PCSetType(pc,PCCHOLESKY);
105: /*
106: Test interface functions of Krylov-Schur solver
107: */
108: EPSKrylovSchurGetRestart(eps,&keep);
109: PetscPrintf(PETSC_COMM_WORLD," Restart parameter before changing = %g",(double)keep);
110: EPSKrylovSchurSetRestart(eps,0.4);
111: EPSKrylovSchurGetRestart(eps,&keep);
112: PetscPrintf(PETSC_COMM_WORLD," ... changed to %g\n",(double)keep);
114: EPSKrylovSchurGetDetectZeros(eps,&detect);
115: PetscPrintf(PETSC_COMM_WORLD," Detect zeros before changing = %d",(int)detect);
116: EPSKrylovSchurSetDetectZeros(eps,PETSC_TRUE);
117: EPSKrylovSchurGetDetectZeros(eps,&detect);
118: PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)detect);
120: EPSKrylovSchurGetLocking(eps,&lock);
121: PetscPrintf(PETSC_COMM_WORLD," Locking flag before changing = %d",(int)lock);
122: EPSKrylovSchurSetLocking(eps,PETSC_FALSE);
123: EPSKrylovSchurGetLocking(eps,&lock);
124: PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)lock);
126: EPSKrylovSchurGetDimensions(eps,&nev,&ncv,&mpd);
127: PetscPrintf(PETSC_COMM_WORLD," Sub-solve dimensions before changing = [%D,%D,%D]",nev,ncv,mpd);
128: EPSKrylovSchurSetDimensions(eps,30,60,60);
129: EPSKrylovSchurGetDimensions(eps,&nev,&ncv,&mpd);
130: PetscPrintf(PETSC_COMM_WORLD," ... changed to [%D,%D,%D]\n",nev,ncv,mpd);
132: if (size>1) {
133: EPSKrylovSchurSetPartitions(eps,size);
134: EPSKrylovSchurGetPartitions(eps,&npart);
135: PetscPrintf(PETSC_COMM_WORLD," Using %D 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++) {
160: PetscPrintf(PETSC_COMM_WORLD," .. %g (%D)\n",(double)shifts[i],inertias[i]);
161: }
162: PetscFree(shifts);
163: PetscFree(inertias);
165: EPSSolve(eps);
166: EPSGetDimensions(eps,&nev,NULL,NULL);
167: EPSGetInterval(eps,&int0,&int1);
168: PetscPrintf(PETSC_COMM_WORLD," Found %D eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1);
170: if (showinertia) {
171: EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias);
172: PetscPrintf(PETSC_COMM_WORLD," Used %D shifts (inertia):\n",k);
173: for (i=0;i<k;i++) {
174: PetscPrintf(PETSC_COMM_WORLD," .. %g (%D)\n",(double)shifts[i],inertias[i]);
175: }
176: PetscFree(shifts);
177: PetscFree(inertias);
178: }
180: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
182: if (size>1) {
183: EPSKrylovSchurGetSubcommInfo(eps,&k,&nval,&v);
184: PetscMalloc1(nval,&evals);
185: for (i=0;i<nval;i++) {
186: EPSKrylovSchurGetSubcommPairs(eps,i,&lambda,v);
187: evals[i] = PetscRealPart(lambda);
188: }
189: PetscFormatRealArray(vlist,sizeof(vlist),"%f",nval,evals);
190: PetscSynchronizedPrintf(PETSC_COMM_WORLD," Process %d has worked in sub-interval %D, containing %D eigenvalues: %s\n",(int)rank,k,nval,vlist);
191: PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
192: VecDestroy(&v);
193: PetscFree(evals);
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: }
223: /*TEST
225: test:
226: suffix: 1
227: nsize: 2
228: args: -showinertia 0 -info_exclude eps,st,rg,bv,ds -log_exclude eps,st,rg,bv,ds
229: requires: !single
231: TEST*/