Actual source code: ex16.c
slepc-3.17.2 2022-08-09
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[] = "Simple quadratic eigenvalue problem.\n\n"
12: "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,nev,i,j,nconv;
23: PetscBool flag,terse;
24: PetscReal error,re,im;
25: PetscScalar kr,ki;
26: Vec xr,xi;
27: BV V;
28: PetscRandom rand;
30: SlepcInitialize(&argc,&argv,(char*)0,help);
32: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
33: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
34: if (!flag) m=n;
35: N = n*m;
36: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);
38: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
39: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
40: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
42: /* K is the 2-D Laplacian */
43: MatCreate(PETSC_COMM_WORLD,&K);
44: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N);
45: MatSetFromOptions(K);
46: MatSetUp(K);
47: MatGetOwnershipRange(K,&Istart,&Iend);
48: for (II=Istart;II<Iend;II++) {
49: i = II/n; j = II-i*n;
50: if (i>0) MatSetValue(K,II,II-n,-1.0,INSERT_VALUES);
51: if (i<m-1) MatSetValue(K,II,II+n,-1.0,INSERT_VALUES);
52: if (j>0) MatSetValue(K,II,II-1,-1.0,INSERT_VALUES);
53: if (j<n-1) MatSetValue(K,II,II+1,-1.0,INSERT_VALUES);
54: MatSetValue(K,II,II,4.0,INSERT_VALUES);
55: }
56: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
57: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
59: /* C is the 1-D Laplacian on horizontal lines */
60: MatCreate(PETSC_COMM_WORLD,&C);
61: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
62: MatSetFromOptions(C);
63: MatSetUp(C);
64: MatGetOwnershipRange(C,&Istart,&Iend);
65: for (II=Istart;II<Iend;II++) {
66: i = II/n; j = II-i*n;
67: if (j>0) MatSetValue(C,II,II-1,-1.0,INSERT_VALUES);
68: if (j<n-1) MatSetValue(C,II,II+1,-1.0,INSERT_VALUES);
69: MatSetValue(C,II,II,2.0,INSERT_VALUES);
70: }
71: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
72: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
74: /* M is a diagonal matrix */
75: MatCreate(PETSC_COMM_WORLD,&M);
76: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
77: MatSetFromOptions(M);
78: MatSetUp(M);
79: MatGetOwnershipRange(M,&Istart,&Iend);
80: for (II=Istart;II<Iend;II++) MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES);
81: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
82: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
84: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: Create the eigensolver and set various options
86: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88: /*
89: Create eigensolver context
90: */
91: PEPCreate(PETSC_COMM_WORLD,&pep);
93: /*
94: Set matrices and problem type
95: */
96: A[0] = K; A[1] = C; A[2] = M;
97: PEPSetOperators(pep,3,A);
98: PEPSetProblemType(pep,PEP_HERMITIAN);
100: /*
101: In complex scalars, use a real initial vector since in this example
102: the matrices are all real, then all vectors generated by the solver
103: will have a zero imaginary part. This is not really necessary.
104: */
105: PEPGetBV(pep,&V);
106: BVGetRandomContext(V,&rand);
107: PetscRandomSetInterval(rand,-1,1);
109: /*
110: Set solver parameters at runtime
111: */
112: PEPSetFromOptions(pep);
114: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115: Solve the eigensystem
116: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118: PEPSolve(pep);
120: /*
121: Optional: Get some information from the solver and display it
122: */
123: PEPGetDimensions(pep,&nev,NULL,NULL);
124: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Display solution and clean up
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: /* show detailed info unless -terse option is given by user */
131: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
132: if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
133: else {
134: PEPGetConverged(pep,&nconv);
135: if (nconv>0) {
136: MatCreateVecs(M,&xr,&xi);
137: /* display eigenvalues and relative errors */
138: PetscCall(PetscPrintf(PETSC_COMM_WORLD,
139: "\n k ||P(k)x||/||kx||\n"
140: " ----------------- ------------------\n"));
141: for (i=0;i<nconv;i++) {
142: /* get converged eigenpairs */
143: PEPGetEigenpair(pep,i,&kr,&ki,xr,xi);
144: /* compute the relative error associated to each eigenpair */
145: PEPComputeError(pep,i,PEP_ERROR_BACKWARD,&error);
146: #if defined(PETSC_USE_COMPLEX)
147: re = PetscRealPart(kr);
148: im = PetscImaginaryPart(kr);
149: #else
150: re = kr;
151: im = ki;
152: #endif
153: if (im!=0.0) PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g\n",(double)re,(double)im,(double)error);
154: else PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",(double)re,(double)error);
155: }
156: PetscPrintf(PETSC_COMM_WORLD,"\n");
157: VecDestroy(&xr);
158: VecDestroy(&xi);
159: }
160: }
161: PEPDestroy(&pep);
162: MatDestroy(&M);
163: MatDestroy(&C);
164: MatDestroy(&K);
165: SlepcFinalize();
166: return 0;
167: }
169: /*TEST
171: testset:
172: args: -pep_nev 4 -pep_ncv 21 -n 12 -terse
173: output_file: output/ex16_1.out
174: test:
175: suffix: 1
176: args: -pep_type {{toar qarnoldi}}
177: test:
178: suffix: 1_linear
179: args: -pep_type linear -pep_linear_explicitmatrix
180: requires: !single
181: test:
182: suffix: 1_linear_symm
183: args: -pep_type linear -pep_linear_explicitmatrix -pep_linear_eps_gen_indefinite -pep_scale scalar -pep_linear_bv_definite_tol 1e-12
184: requires: !single
185: test:
186: suffix: 1_stoar
187: args: -pep_type stoar -pep_scale scalar
188: requires: double !cuda
189: test:
190: suffix: 1_stoar_t
191: args: -pep_type stoar -pep_scale scalar -st_transform
192: requires: double !cuda
194: TEST*/