Actual source code: test7.c
slepc-3.19.2 2023-09-05
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 STOAR.\n\n"
12: "This is based on ex38.c. The command line options are:\n"
13: " -n <n> ... dimension of the matrices.\n\n";
15: #include <slepcpep.h>
17: int main(int argc,char **argv)
18: {
19: Mat M,C,K,A[3]; /* problem matrices */
20: PEP pep; /* polynomial eigenproblem solver context */
21: ST st; /* spectral transformation context */
22: KSP ksp;
23: PC pc;
24: PetscBool showinertia=PETSC_TRUE,lock,detect,checket;
25: PetscInt n=100,Istart,Iend,i,*inertias,ns,nev,ncv,mpd;
26: PetscReal mu=1.0,tau=10.0,kappa=5.0,int0,int1,*shifts;
28: PetscFunctionBeginUser;
29: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
31: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
32: PetscCall(PetscOptionsGetBool(NULL,NULL,"-showinertia",&showinertia,NULL));
33: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSpectrum slicing on PEP, n=%" PetscInt_FMT "\n\n",n));
35: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39: /* K is a tridiagonal */
40: PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
41: PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
42: PetscCall(MatSetFromOptions(K));
43: PetscCall(MatSetUp(K));
45: PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
46: for (i=Istart;i<Iend;i++) {
47: if (i>0) PetscCall(MatSetValue(K,i,i-1,-kappa,INSERT_VALUES));
48: PetscCall(MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES));
49: if (i<n-1) PetscCall(MatSetValue(K,i,i+1,-kappa,INSERT_VALUES));
50: }
52: PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
53: PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
55: /* C is a tridiagonal */
56: PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
57: PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n));
58: PetscCall(MatSetFromOptions(C));
59: PetscCall(MatSetUp(C));
61: PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
62: for (i=Istart;i<Iend;i++) {
63: if (i>0) PetscCall(MatSetValue(C,i,i-1,-tau,INSERT_VALUES));
64: PetscCall(MatSetValue(C,i,i,tau*3.0,INSERT_VALUES));
65: if (i<n-1) PetscCall(MatSetValue(C,i,i+1,-tau,INSERT_VALUES));
66: }
68: PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
69: PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
71: /* M is a diagonal matrix */
72: PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
73: PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
74: PetscCall(MatSetFromOptions(M));
75: PetscCall(MatSetUp(M));
76: PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
77: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(M,i,i,mu,INSERT_VALUES));
78: PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
79: PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
81: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82: Create the eigensolver and solve the problem
83: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
86: A[0] = K; A[1] = C; A[2] = M;
87: PetscCall(PEPSetOperators(pep,3,A));
88: PetscCall(PEPSetProblemType(pep,PEP_HYPERBOLIC));
89: PetscCall(PEPSetType(pep,PEPSTOAR));
91: /*
92: Set interval and other settings for spectrum slicing
93: */
94: int0 = -11.3;
95: int1 = -9.5;
96: PetscCall(PEPSetInterval(pep,int0,int1));
97: PetscCall(PEPSetWhichEigenpairs(pep,PEP_ALL));
98: PetscCall(PEPGetST(pep,&st));
99: PetscCall(STSetType(st,STSINVERT));
100: PetscCall(STGetKSP(st,&ksp));
101: PetscCall(KSPSetType(ksp,KSPPREONLY));
102: PetscCall(KSPGetPC(ksp,&pc));
103: PetscCall(PCSetType(pc,PCCHOLESKY));
105: /*
106: Test interface functions of STOAR solver
107: */
108: PetscCall(PEPSTOARGetDetectZeros(pep,&detect));
109: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Detect zeros before changing = %d",(int)detect));
110: PetscCall(PEPSTOARSetDetectZeros(pep,PETSC_TRUE));
111: PetscCall(PEPSTOARGetDetectZeros(pep,&detect));
112: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)detect));
114: PetscCall(PEPSTOARGetLocking(pep,&lock));
115: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Locking flag before changing = %d",(int)lock));
116: PetscCall(PEPSTOARSetLocking(pep,PETSC_TRUE));
117: PetscCall(PEPSTOARGetLocking(pep,&lock));
118: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)lock));
120: PetscCall(PEPSTOARGetCheckEigenvalueType(pep,&checket));
121: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Check eigenvalue type flag before changing = %d",(int)checket));
122: PetscCall(PEPSTOARSetCheckEigenvalueType(pep,PETSC_FALSE));
123: PetscCall(PEPSTOARGetCheckEigenvalueType(pep,&checket));
124: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)checket));
126: PetscCall(PEPSTOARGetDimensions(pep,&nev,&ncv,&mpd));
127: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Sub-solve dimensions before changing = [%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]",nev,ncv,mpd));
128: PetscCall(PEPSTOARSetDimensions(pep,30,60,60));
129: PetscCall(PEPSTOARGetDimensions(pep,&nev,&ncv,&mpd));
130: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to [%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]\n",nev,ncv,mpd));
132: PetscCall(PEPSetFromOptions(pep));
134: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135: Compute all eigenvalues in interval and display info
136: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138: PetscCall(PEPSetUp(pep));
139: PetscCall(PEPSTOARGetInertias(pep,&ns,&shifts,&inertias));
140: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Inertias (after setup):\n"));
141: for (i=0;i<ns;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," .. %g (%" PetscInt_FMT ")\n",(double)shifts[i],inertias[i]));
142: PetscCall(PetscFree(shifts));
143: PetscCall(PetscFree(inertias));
145: PetscCall(PEPSolve(pep));
146: PetscCall(PEPGetDimensions(pep,&nev,NULL,NULL));
147: PetscCall(PEPGetInterval(pep,&int0,&int1));
148: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Found %" PetscInt_FMT " eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1));
150: if (showinertia) {
151: PetscCall(PEPSTOARGetInertias(pep,&ns,&shifts,&inertias));
152: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Used %" PetscInt_FMT " shifts (inertia):\n",ns));
153: for (i=0;i<ns;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," .. %g (%" PetscInt_FMT ")\n",(double)shifts[i],inertias[i]));
154: PetscCall(PetscFree(shifts));
155: PetscCall(PetscFree(inertias));
156: }
158: PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
160: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161: Clean up
162: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163: PetscCall(PEPDestroy(&pep));
164: PetscCall(MatDestroy(&M));
165: PetscCall(MatDestroy(&C));
166: PetscCall(MatDestroy(&K));
167: PetscCall(SlepcFinalize());
168: return 0;
169: }
171: /*TEST
173: test:
174: requires: !single
175: args: -showinertia 0
177: TEST*/