Actual source code: test31.c
slepc-3.22.1 2024-10-28
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,NULL,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(MatGetOwnershipRange(A,&Istart,&Iend));
45: for (II=Istart;II<Iend;II++) {
46: i = II/n; j = II-i*n;
47: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
48: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
49: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
50: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
51: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
52: }
53: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
54: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Create the eigensolver and set various options
58: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
61: PetscCall(EPSSetOperators(eps,A,NULL));
62: PetscCall(EPSSetProblemType(eps,EPS_HEP));
63: PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));
64: PetscCall(EPSSetWhichEigenpairs(eps,EPS_ALL));
65: PetscCall(EPSSetInterval(eps,0.5,1.3));
66: PetscCall(EPSGetST(eps,&st));
67: PetscCall(STSetType(st,STFILTER));
68: PetscCall(EPSSetFromOptions(eps));
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Solve the problem and display the solution
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74: PetscCall(EPSSolve(eps));
76: /* print filter information */
77: PetscCall(PetscObjectTypeCompare((PetscObject)st,STFILTER,&flag));
78: if (flag) {
79: PetscCall(STFilterGetDegree(st,°ree));
80: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Filter degree: %" PetscInt_FMT "\n",degree));
81: PetscCall(STFilterGetInterval(st,&inta,&intb));
82: PetscCall(STFilterGetRange(st,&rleft,&rright));
83: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Requested interval: [%g,%g], range: [%g,%g]\n\n",(double)inta,(double)intb,(double)rleft,(double)rright));
84: }
86: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
87: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
88: else {
89: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
90: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
91: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
92: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
93: }
95: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96: Solve the problem again after changing the matrix
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: if (modify) {
99: PetscCall(MatSetValue(A,0,0,0.3,INSERT_VALUES));
100: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
101: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
102: PetscCall(EPSSetOperators(eps,A,NULL));
103: PetscCall(EPSSolve(eps));
104: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
105: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
106: else {
107: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
108: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
109: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
110: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
111: }
112: }
114: PetscCall(EPSDestroy(&eps));
115: PetscCall(MatDestroy(&A));
116: PetscCall(SlepcFinalize());
117: return 0;
118: }
120: /*TEST
122: test:
123: suffix: 1
124: args: -terse
125: filter: sed -e "s/0.161982,7.83797/0.162007,7.83897/"
126: requires: !single
128: test:
129: suffix: 2
130: args: -modify -st_filter_range -0.5,8 -terse
131: requires: !single
133: test:
134: suffix: 3
135: args: -modify -terse
136: filter: sed -e "s/0.161982,7.83797/0.162007,7.83897/"
137: requires: !single
139: TEST*/