Actual source code: test8.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 polynomial JD.\n\n"
12: "This is based on ex16.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 <slepcpep.h>
18: int main(int argc,char **argv)
19: {
20: Mat M,C,K,A[3]; /* problem matrices */
21: PEP pep; /* polynomial eigenproblem solver context */
22: PetscInt N,n=10,m,Istart,Iend,II,i,j,midx;
23: PetscReal restart,fix;
24: PetscBool flag,reuse;
25: PEPJDProjection proj;
27: PetscFunctionBeginUser;
28: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
30: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
31: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
32: if (!flag) m=n;
33: N = n*m;
34: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
36: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
38: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40: /* K is the 2-D Laplacian */
41: PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
42: PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N));
43: PetscCall(MatSetFromOptions(K));
44: PetscCall(MatSetUp(K));
45: PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
46: for (II=Istart;II<Iend;II++) {
47: i = II/n; j = II-i*n;
48: if (i>0) PetscCall(MatSetValue(K,II,II-n,-1.0,INSERT_VALUES));
49: if (i<m-1) PetscCall(MatSetValue(K,II,II+n,-1.0,INSERT_VALUES));
50: if (j>0) PetscCall(MatSetValue(K,II,II-1,-1.0,INSERT_VALUES));
51: if (j<n-1) PetscCall(MatSetValue(K,II,II+1,-1.0,INSERT_VALUES));
52: PetscCall(MatSetValue(K,II,II,4.0,INSERT_VALUES));
53: }
54: PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
55: PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
57: /* C is the 1-D Laplacian on horizontal lines */
58: PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
59: PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N));
60: PetscCall(MatSetFromOptions(C));
61: PetscCall(MatSetUp(C));
62: PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
63: for (II=Istart;II<Iend;II++) {
64: i = II/n; j = II-i*n;
65: if (j>0) PetscCall(MatSetValue(C,II,II-1,-1.0,INSERT_VALUES));
66: if (j<n-1) PetscCall(MatSetValue(C,II,II+1,-1.0,INSERT_VALUES));
67: PetscCall(MatSetValue(C,II,II,2.0,INSERT_VALUES));
68: }
69: PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
70: PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
72: /* M is a diagonal matrix */
73: PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
74: PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N));
75: PetscCall(MatSetFromOptions(M));
76: PetscCall(MatSetUp(M));
77: PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
78: for (II=Istart;II<Iend;II++) PetscCall(MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES));
79: PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
80: PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Create the eigensolver and set various options
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
87: A[0] = K; A[1] = C; A[2] = M;
88: PetscCall(PEPSetOperators(pep,3,A));
89: PetscCall(PEPSetType(pep,PEPJD));
91: /*
92: Test interface functions of JD solver
93: */
94: PetscCall(PEPJDGetRestart(pep,&restart));
95: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Restart parameter before changing = %g",(double)restart));
96: PetscCall(PEPJDSetRestart(pep,0.3));
97: PetscCall(PEPJDGetRestart(pep,&restart));
98: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %g\n",(double)restart));
100: PetscCall(PEPJDGetFix(pep,&fix));
101: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Fix parameter before changing = %g",(double)fix));
102: PetscCall(PEPJDSetFix(pep,0.001));
103: PetscCall(PEPJDGetFix(pep,&fix));
104: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %g\n",(double)fix));
106: PetscCall(PEPJDGetReusePreconditioner(pep,&reuse));
107: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Reuse preconditioner flag before changing = %d",(int)reuse));
108: PetscCall(PEPJDSetReusePreconditioner(pep,PETSC_TRUE));
109: PetscCall(PEPJDGetReusePreconditioner(pep,&reuse));
110: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)reuse));
112: PetscCall(PEPJDGetProjection(pep,&proj));
113: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Projection type before changing = %d",(int)proj));
114: PetscCall(PEPJDSetProjection(pep,PEP_JD_PROJECTION_ORTHOGONAL));
115: PetscCall(PEPJDGetProjection(pep,&proj));
116: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)proj));
118: PetscCall(PEPJDGetMinimalityIndex(pep,&midx));
119: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Minimality index before changing = %" PetscInt_FMT,midx));
120: PetscCall(PEPJDSetMinimalityIndex(pep,2));
121: PetscCall(PEPJDGetMinimalityIndex(pep,&midx));
122: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %" PetscInt_FMT "\n",midx));
124: PetscCall(PEPSetFromOptions(pep));
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Solve the eigensystem
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: PetscCall(PEPSolve(pep));
131: PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
132: PetscCall(PEPDestroy(&pep));
133: PetscCall(MatDestroy(&M));
134: PetscCall(MatDestroy(&C));
135: PetscCall(MatDestroy(&K));
136: PetscCall(SlepcFinalize());
137: return 0;
138: }
140: /*TEST
142: test:
143: args: -n 12 -pep_nev 2 -pep_ncv 21 -pep_conv_abs
145: TEST*/