Actual source code: ex50.c
slepc-3.20.2 2024-03-15
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[] = "User-defined split preconditioner when solving a quadratic eigenproblem.\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 A[3],P[3]; /* problem matrices and split preconditioner matrices */
21: PEP pep; /* polynomial eigenproblem solver context */
22: ST st;
23: PetscInt N,n=10,m,Istart,Iend,II,i,j;
24: PetscBool flag,terse;
26: PetscFunctionBeginUser;
27: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
29: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
30: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
31: if (!flag) m=n;
32: N = n*m;
33: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
35: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: Compute the matrices for (k^2*A_2+k*A_1+A_0)x=0, and their approximations
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39: /* A[0] is the 2-D Laplacian */
40: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[0]));
41: PetscCall(MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,N,N));
42: PetscCall(MatSetFromOptions(A[0]));
43: PetscCall(MatSetUp(A[0]));
44: PetscCall(MatCreate(PETSC_COMM_WORLD,&P[0]));
45: PetscCall(MatSetSizes(P[0],PETSC_DECIDE,PETSC_DECIDE,N,N));
46: PetscCall(MatSetFromOptions(P[0]));
47: PetscCall(MatSetUp(P[0]));
49: PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
50: for (II=Istart;II<Iend;II++) {
51: i = II/n; j = II-i*n;
52: if (i>0) PetscCall(MatSetValue(A[0],II,II-n,-1.0,INSERT_VALUES));
53: if (i<m-1) PetscCall(MatSetValue(A[0],II,II+n,-1.0,INSERT_VALUES));
54: if (j>0) PetscCall(MatSetValue(A[0],II,II-1,-1.0,INSERT_VALUES));
55: if (j<n-1) PetscCall(MatSetValue(A[0],II,II+1,-1.0,INSERT_VALUES));
56: PetscCall(MatSetValue(A[0],II,II,4.0,INSERT_VALUES));
57: if (j>0) PetscCall(MatSetValue(P[0],II,II-1,-1.0,INSERT_VALUES));
58: if (j<n-1) PetscCall(MatSetValue(P[0],II,II+1,-1.0,INSERT_VALUES));
59: PetscCall(MatSetValue(P[0],II,II,4.0,INSERT_VALUES));
60: }
61: PetscCall(MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY));
62: PetscCall(MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY));
63: PetscCall(MatAssemblyBegin(P[0],MAT_FINAL_ASSEMBLY));
64: PetscCall(MatAssemblyEnd(P[0],MAT_FINAL_ASSEMBLY));
66: /* A[1] is the 1-D Laplacian on horizontal lines */
67: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[1]));
68: PetscCall(MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,N,N));
69: PetscCall(MatSetFromOptions(A[1]));
70: PetscCall(MatSetUp(A[1]));
71: PetscCall(MatCreate(PETSC_COMM_WORLD,&P[1]));
72: PetscCall(MatSetSizes(P[1],PETSC_DECIDE,PETSC_DECIDE,N,N));
73: PetscCall(MatSetFromOptions(P[1]));
74: PetscCall(MatSetUp(P[1]));
76: PetscCall(MatGetOwnershipRange(A[1],&Istart,&Iend));
77: for (II=Istart;II<Iend;II++) {
78: i = II/n; j = II-i*n;
79: if (j>0) PetscCall(MatSetValue(A[1],II,II-1,-1.0,INSERT_VALUES));
80: if (j<n-1) PetscCall(MatSetValue(A[1],II,II+1,-1.0,INSERT_VALUES));
81: PetscCall(MatSetValue(A[1],II,II,2.0,INSERT_VALUES));
82: PetscCall(MatSetValue(P[1],II,II,2.0,INSERT_VALUES));
83: }
84: PetscCall(MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY));
85: PetscCall(MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY));
86: PetscCall(MatAssemblyBegin(P[1],MAT_FINAL_ASSEMBLY));
87: PetscCall(MatAssemblyEnd(P[1],MAT_FINAL_ASSEMBLY));
89: /* A[2] is a diagonal matrix */
90: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[2]));
91: PetscCall(MatSetSizes(A[2],PETSC_DECIDE,PETSC_DECIDE,N,N));
92: PetscCall(MatSetFromOptions(A[2]));
93: PetscCall(MatSetUp(A[2]));
94: PetscCall(MatCreate(PETSC_COMM_WORLD,&P[2]));
95: PetscCall(MatSetSizes(P[2],PETSC_DECIDE,PETSC_DECIDE,N,N));
96: PetscCall(MatSetFromOptions(P[2]));
97: PetscCall(MatSetUp(P[2]));
99: PetscCall(MatGetOwnershipRange(A[2],&Istart,&Iend));
100: for (II=Istart;II<Iend;II++) {
101: PetscCall(MatSetValue(A[2],II,II,(PetscReal)(II+1),INSERT_VALUES));
102: PetscCall(MatSetValue(P[2],II,II,(PetscReal)(II+1),INSERT_VALUES));
103: }
104: PetscCall(MatAssemblyBegin(A[2],MAT_FINAL_ASSEMBLY));
105: PetscCall(MatAssemblyEnd(A[2],MAT_FINAL_ASSEMBLY));
106: PetscCall(MatAssemblyBegin(P[2],MAT_FINAL_ASSEMBLY));
107: PetscCall(MatAssemblyEnd(P[2],MAT_FINAL_ASSEMBLY));
109: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: Create the eigensolver and set various options
111: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
114: PetscCall(PEPSetOperators(pep,3,A));
115: PetscCall(PEPSetProblemType(pep,PEP_HERMITIAN));
117: PetscCall(PEPGetST(pep,&st));
118: PetscCall(STSetType(st,STSINVERT));
119: PetscCall(STSetSplitPreconditioner(st,3,P,SUBSET_NONZERO_PATTERN));
121: PetscCall(PEPSetTarget(pep,-2.0));
122: PetscCall(PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE));
124: /*
125: Set solver parameters at runtime
126: */
127: PetscCall(PEPSetFromOptions(pep));
129: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130: Solve the eigensystem, display solution and clean up
131: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133: PetscCall(PEPSolve(pep));
134: /* show detailed info unless -terse option is given by user */
135: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
136: if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
137: else {
138: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
139: PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
140: PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD));
141: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
142: }
143: PetscCall(PEPDestroy(&pep));
144: PetscCall(MatDestroy(&A[0]));
145: PetscCall(MatDestroy(&A[1]));
146: PetscCall(MatDestroy(&A[2]));
147: PetscCall(MatDestroy(&P[0]));
148: PetscCall(MatDestroy(&P[1]));
149: PetscCall(MatDestroy(&P[2]));
150: PetscCall(SlepcFinalize());
151: return 0;
152: }
154: /*TEST
156: testset:
157: args: -pep_nev 4 -pep_ncv 28 -n 12 -terse
158: output_file: output/ex50_1.out
159: requires: double
160: test:
161: suffix: 1
162: args: -pep_type {{toar qarnoldi}}
163: test:
164: suffix: 1_linear
165: args: -pep_type linear -pep_general
167: testset:
168: args: -pep_all -n 12 -pep_type ciss -rg_type ellipse -rg_ellipse_center -1+1.5i -rg_ellipse_radius .3 -terse
169: output_file: output/ex50_2.out
170: requires: complex double
171: timeoutfactor: 2
172: test:
173: suffix: 2
174: test:
175: suffix: 2_par
176: nsize: 2
177: args: -pep_ciss_partitions 2
179: TEST*/