Actual source code: test17.c
slepc-3.17.1 2022-04-11
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: SlepcInitialize(&argc,&argv,(char*)0,help);
36: MPI_Comm_size(PETSC_COMM_WORLD,&size);
37: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
39: PetscOptionsGetBool(NULL,NULL,"-showinertia",&showinertia,NULL);
40: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
41: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
42: if (!flag) m=n;
43: N = n*m;
44: PetscPrintf(PETSC_COMM_WORLD,"\nSpectrum-slicing test, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);
46: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: Compute the matrices that define the eigensystem, Ax=kBx
48: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50: MatCreate(PETSC_COMM_WORLD,&A);
51: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
52: MatSetFromOptions(A);
53: MatSetUp(A);
55: MatCreate(PETSC_COMM_WORLD,&B);
56: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
57: MatSetFromOptions(B);
58: MatSetUp(B);
60: MatGetOwnershipRange(A,&Istart,&Iend);
61: for (II=Istart;II<Iend;II++) {
62: i = II/n; j = II-i*n;
63: if (i>0) MatSetValue(A,II,II-n,-1.0,INSERT_VALUES);
64: if (i<m-1) MatSetValue(A,II,II+n,-1.0,INSERT_VALUES);
65: if (j>0) MatSetValue(A,II,II-1,-1.0,INSERT_VALUES);
66: if (j<n-1) MatSetValue(A,II,II+1,-1.0,INSERT_VALUES);
67: MatSetValue(A,II,II,4.0,INSERT_VALUES);
68: MatSetValue(B,II,II,2.0,INSERT_VALUES);
69: }
70: if (Istart==0) {
71: MatSetValue(B,0,0,6.0,INSERT_VALUES);
72: MatSetValue(B,0,1,-1.0,INSERT_VALUES);
73: MatSetValue(B,1,0,-1.0,INSERT_VALUES);
74: MatSetValue(B,1,1,1.0,INSERT_VALUES);
75: }
77: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
78: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
79: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
80: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Create the eigensolver and set various options
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: EPSCreate(PETSC_COMM_WORLD,&eps);
87: EPSSetOperators(eps,A,B);
88: EPSSetProblemType(eps,EPS_GHEP);
89: EPSSetType(eps,EPSKRYLOVSCHUR);
91: /*
92: Set interval and other settings for spectrum slicing
93: */
94: EPSSetWhichEigenpairs(eps,EPS_ALL);
95: int0 = 1.1; int1 = 1.3;
96: EPSSetInterval(eps,int0,int1);
97: EPSGetST(eps,&st);
98: STSetType(st,STSINVERT);
99: if (size>1) EPSKrylovSchurSetPartitions(eps,size);
100: EPSKrylovSchurGetKSP(eps,&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 = [%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]",nev,ncv,mpd);
128: EPSKrylovSchurSetDimensions(eps,30,60,60);
129: EPSKrylovSchurGetDimensions(eps,&nev,&ncv,&mpd);
130: PetscPrintf(PETSC_COMM_WORLD," ... changed to [%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]\n",nev,ncv,mpd);
132: if (size>1) {
133: EPSKrylovSchurGetPartitions(eps,&npart);
134: PetscPrintf(PETSC_COMM_WORLD," Using %" PetscInt_FMT " partitions\n",npart);
136: PetscMalloc1(npart+1,&subint);
137: subint[0] = int0;
138: subint[npart] = int1;
139: for (i=1;i<npart;i++) subint[i] = int0+i*(int1-int0)/npart;
140: EPSKrylovSchurSetSubintervals(eps,subint);
141: PetscFree(subint);
142: EPSKrylovSchurGetSubintervals(eps,&subint);
143: PetscPrintf(PETSC_COMM_WORLD," Using sub-interval separations = ");
144: for (i=1;i<npart;i++) PetscPrintf(PETSC_COMM_WORLD," %g",(double)subint[i]);
145: PetscFree(subint);
146: PetscPrintf(PETSC_COMM_WORLD,"\n");
147: }
149: EPSSetFromOptions(eps);
151: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152: Compute all eigenvalues in interval and display info
153: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155: EPSSetUp(eps);
156: EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias);
157: PetscPrintf(PETSC_COMM_WORLD," Inertias after EPSSetUp:\n");
158: for (i=0;i<k;i++) PetscPrintf(PETSC_COMM_WORLD," .. %g (%" PetscInt_FMT ")\n",(double)shifts[i],inertias[i]);
159: PetscFree(shifts);
160: PetscFree(inertias);
162: EPSSolve(eps);
163: EPSGetDimensions(eps,&nev,NULL,NULL);
164: EPSGetInterval(eps,&int0,&int1);
165: PetscPrintf(PETSC_COMM_WORLD," Found %" PetscInt_FMT " eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1);
167: if (showinertia) {
168: EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias);
169: PetscPrintf(PETSC_COMM_WORLD," Used %" PetscInt_FMT " shifts (inertia):\n",k);
170: for (i=0;i<k;i++) PetscPrintf(PETSC_COMM_WORLD," .. %g (%" PetscInt_FMT ")\n",(double)shifts[i],inertias[i]);
171: PetscFree(shifts);
172: PetscFree(inertias);
173: }
175: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
177: if (size>1) {
178: EPSKrylovSchurGetSubcommInfo(eps,&k,&nval,&v);
179: PetscMalloc1(nval,&evals);
180: for (i=0;i<nval;i++) {
181: EPSKrylovSchurGetSubcommPairs(eps,i,&lambda,v);
182: evals[i] = PetscRealPart(lambda);
183: }
184: PetscFormatRealArray(vlist,sizeof(vlist),"%f",nval,evals);
185: PetscSynchronizedPrintf(PETSC_COMM_WORLD," Process %d has worked in sub-interval %" PetscInt_FMT ", containing %" PetscInt_FMT " eigenvalues: %s\n",(int)rank,k,nval,vlist);
186: PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
187: VecDestroy(&v);
188: PetscFree(evals);
190: EPSKrylovSchurGetSubcommMats(eps,&As,&Bs);
191: MatGetLocalSize(A,&nloc,NULL);
192: MatGetLocalSize(As,&nlocs,&mlocs);
193: 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);
194: PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
196: /* modify A on subcommunicators */
197: MatCreate(PetscObjectComm((PetscObject)As),&Au);
198: MatSetSizes(Au,nlocs,mlocs,N,N);
199: MatSetFromOptions(Au);
200: MatSetUp(Au);
201: MatGetOwnershipRange(Au,&Istart,&Iend);
202: for (II=Istart;II<Iend;II++) MatSetValue(Au,II,II,0.5,INSERT_VALUES);
203: MatAssemblyBegin(Au,MAT_FINAL_ASSEMBLY);
204: MatAssemblyEnd(Au,MAT_FINAL_ASSEMBLY);
205: PetscPrintf(PETSC_COMM_WORLD," Updating internal matrices\n");
206: EPSKrylovSchurUpdateSubcommMats(eps,1.1,-5.0,Au,1.0,0.0,NULL,DIFFERENT_NONZERO_PATTERN,PETSC_TRUE);
207: MatDestroy(&Au);
208: EPSSolve(eps);
209: EPSGetDimensions(eps,&nev,NULL,NULL);
210: EPSGetInterval(eps,&int0,&int1);
211: PetscPrintf(PETSC_COMM_WORLD," After update, found %" PetscInt_FMT " eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1);
212: }
213: EPSDestroy(&eps);
214: MatDestroy(&A);
215: MatDestroy(&B);
216: SlepcFinalize();
217: return 0;
218: }
220: /*TEST
222: test:
223: suffix: 1
224: nsize: 2
225: args: -showinertia 0 -log_exclude eps,st,rg,bv,ds
226: requires: !single
228: test:
229: suffix: 2
230: nsize: 1
231: args: -showinertia 0 -log_exclude eps,st,rg,bv,ds
232: requires: !single
234: TEST*/