Actual source code: test31.c
slepc-3.14.0 2020-09-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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,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);
36: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37: Create the 2-D Laplacian
38: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40: MatCreate(PETSC_COMM_WORLD,&A);
41: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
42: MatSetFromOptions(A);
43: MatSetUp(A);
44: MatGetOwnershipRange(A,&Istart,&Iend);
45: for (II=Istart;II<Iend;II++) {
46: i = II/n; j = II-i*n;
47: if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
48: if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
49: if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
50: if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
51: MatSetValue(A,II,II,4.0,INSERT_VALUES);
52: }
53: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
54: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Create the eigensolver and set various options
58: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60: EPSCreate(PETSC_COMM_WORLD,&eps);
61: EPSSetOperators(eps,A,NULL);
62: EPSSetProblemType(eps,EPS_HEP);
63: EPSSetType(eps,EPSKRYLOVSCHUR);
64: EPSSetWhichEigenpairs(eps,EPS_ALL);
65: EPSSetInterval(eps,0.5,1.3);
66: EPSGetST(eps,&st);
67: STSetType(st,STFILTER);
68: EPSSetFromOptions(eps);
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Solve the problem and display the solution
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74: EPSSolve(eps);
76: /* print filter information */
77: PetscObjectTypeCompare((PetscObject)st,STFILTER,&flag);
78: if (flag) {
79: STFilterGetDegree(st,°ree);
80: PetscPrintf(PETSC_COMM_WORLD," Filter degree: %D\n",degree);
81: STFilterGetInterval(st,&inta,&intb);
82: STFilterGetRange(st,&rleft,&rright);
83: PetscPrintf(PETSC_COMM_WORLD," Requested interval: [%g,%g], range: [%g,%g]\n\n",(double)inta,(double)intb,(double)rleft,(double)rright);
84: }
86: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
87: if (terse) {
88: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
89: } else {
90: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
91: EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
92: EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
93: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
94: }
96: EPSDestroy(&eps);
97: MatDestroy(&A);
98: SlepcFinalize();
99: return ierr;
100: }
102: /*TEST
104: test:
105: suffix: 1
106: args: -terse
107: requires: !single
109: TEST*/