Actual source code: test31.c
slepc-3.15.2 2021-09-20
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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 STFILTER interface functions.\n\n"
12: "Based on ex2.\n"
13: "The command line options are:\n"
14: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
15: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
17: #include <slepceps.h>
19: int main(int argc,char **argv)
20: {
21: Mat A;
22: EPS eps;
23: ST st;
24: PetscInt N,n=10,m,Istart,Iend,II,i,j,degree;
25: PetscBool flag,modify=PETSC_FALSE,terse;
26: PetscReal inta,intb,rleft,rright;
29: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
30: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
31: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
32: if (!flag) m=n;
33: N = n*m;
34: PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%D (%Dx%D grid)\n\n",N,n,m);
35: PetscOptionsGetBool(NULL,NULL,"-modify",&modify,&flag);
37: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38: Create the 2-D Laplacian
39: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
41: MatCreate(PETSC_COMM_WORLD,&A);
42: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
43: MatSetFromOptions(A);
44: MatSetUp(A);
45: MatGetOwnershipRange(A,&Istart,&Iend);
46: for (II=Istart;II<Iend;II++) {
47: i = II/n; j = II-i*n;
48: if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
49: if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
50: if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
51: if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
52: MatSetValue(A,II,II,4.0,INSERT_VALUES);
53: }
54: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
55: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Create the eigensolver and set various options
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: EPSCreate(PETSC_COMM_WORLD,&eps);
62: EPSSetOperators(eps,A,NULL);
63: EPSSetProblemType(eps,EPS_HEP);
64: EPSSetType(eps,EPSKRYLOVSCHUR);
65: EPSSetWhichEigenpairs(eps,EPS_ALL);
66: EPSSetInterval(eps,0.5,1.3);
67: EPSGetST(eps,&st);
68: STSetType(st,STFILTER);
69: EPSSetFromOptions(eps);
71: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72: Solve the problem and display the solution
73: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75: EPSSolve(eps);
77: /* print filter information */
78: PetscObjectTypeCompare((PetscObject)st,STFILTER,&flag);
79: if (flag) {
80: STFilterGetDegree(st,°ree);
81: PetscPrintf(PETSC_COMM_WORLD," Filter degree: %D\n",degree);
82: STFilterGetInterval(st,&inta,&intb);
83: STFilterGetRange(st,&rleft,&rright);
84: PetscPrintf(PETSC_COMM_WORLD," Requested interval: [%g,%g], range: [%g,%g]\n\n",(double)inta,(double)intb,(double)rleft,(double)rright);
85: }
87: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
88: if (terse) {
89: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
90: } else {
91: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
92: EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
93: EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
94: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
95: }
97: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98: Solve the problem again after changing the matrix
99: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100: if (modify) {
101: MatSetValue(A,0,0,0.3,INSERT_VALUES);
102: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
103: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
104: EPSSetOperators(eps,A,NULL);
105: EPSSolve(eps);
106: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
107: if (terse) {
108: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
109: } else {
110: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
111: EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
112: EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
113: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
114: }
115: }
117: EPSDestroy(&eps);
118: MatDestroy(&A);
119: SlepcFinalize();
120: return ierr;
121: }
123: /*TEST
125: test:
126: suffix: 1
127: args: -terse
128: filter: sed -e "s/0.161982,7.83797/0.162007,7.83897/"
129: requires: !single
131: test:
132: suffix: 2
133: args: -modify -st_filter_range -0.5,8 -terse
134: requires: !single
136: test:
137: suffix: 3
138: args: -modify -terse
139: filter: sed -e "s/0.161982,7.83797/0.162007,7.83897/"
140: requires: !single
142: TEST*/