Actual source code: test31.c
slepc-3.19.0 2023-03-31
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 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;
28: PetscFunctionBeginUser;
29: 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,"\n2-D Laplacian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
35: PetscCall(PetscOptionsGetBool(NULL,NULL,"-modify",&modify,&flag));
37: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38: Create the 2-D Laplacian
39: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
41: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
42: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
43: PetscCall(MatSetFromOptions(A));
44: PetscCall(MatSetUp(A));
45: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
46: for (II=Istart;II<Iend;II++) {
47: i = II/n; j = II-i*n;
48: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
49: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
50: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
51: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
52: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
53: }
54: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
55: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Create the eigensolver and set various options
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
62: PetscCall(EPSSetOperators(eps,A,NULL));
63: PetscCall(EPSSetProblemType(eps,EPS_HEP));
64: PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));
65: PetscCall(EPSSetWhichEigenpairs(eps,EPS_ALL));
66: PetscCall(EPSSetInterval(eps,0.5,1.3));
67: PetscCall(EPSGetST(eps,&st));
68: PetscCall(STSetType(st,STFILTER));
69: PetscCall(EPSSetFromOptions(eps));
71: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72: Solve the problem and display the solution
73: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75: PetscCall(EPSSolve(eps));
77: /* print filter information */
78: PetscCall(PetscObjectTypeCompare((PetscObject)st,STFILTER,&flag));
79: if (flag) {
80: PetscCall(STFilterGetDegree(st,°ree));
81: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Filter degree: %" PetscInt_FMT "\n",degree));
82: PetscCall(STFilterGetInterval(st,&inta,&intb));
83: PetscCall(STFilterGetRange(st,&rleft,&rright));
84: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Requested interval: [%g,%g], range: [%g,%g]\n\n",(double)inta,(double)intb,(double)rleft,(double)rright));
85: }
87: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
88: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
89: else {
90: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
91: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
92: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
93: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
94: }
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Solve the problem again after changing the matrix
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: if (modify) {
100: PetscCall(MatSetValue(A,0,0,0.3,INSERT_VALUES));
101: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
102: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
103: PetscCall(EPSSetOperators(eps,A,NULL));
104: PetscCall(EPSSolve(eps));
105: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
106: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
107: else {
108: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
109: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
110: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
111: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
112: }
113: }
115: PetscCall(EPSDestroy(&eps));
116: PetscCall(MatDestroy(&A));
117: PetscCall(SlepcFinalize());
118: return 0;
119: }
121: /*TEST
123: test:
124: suffix: 1
125: args: -terse
126: filter: sed -e "s/0.161982,7.83797/0.162007,7.83897/"
127: requires: !single
129: test:
130: suffix: 2
131: args: -modify -st_filter_range -0.5,8 -terse
132: requires: !single
134: test:
135: suffix: 3
136: args: -modify -terse
137: filter: sed -e "s/0.161982,7.83797/0.162007,7.83897/"
138: requires: !single
140: TEST*/