Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
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";
16 :
17 : #include <slepceps.h>
18 :
19 3 : int main(int argc,char **argv)
20 : {
21 3 : Mat A;
22 3 : EPS eps;
23 3 : ST st;
24 3 : PetscInt N,n=10,m,Istart,Iend,II,i,j,degree;
25 3 : PetscBool flag,modify=PETSC_FALSE,terse;
26 3 : PetscReal inta,intb,rleft,rright;
27 :
28 3 : PetscFunctionBeginUser;
29 3 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
30 3 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
31 3 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
32 3 : if (!flag) m=n;
33 3 : N = n*m;
34 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
35 3 : PetscCall(PetscOptionsGetBool(NULL,NULL,"-modify",&modify,&flag));
36 :
37 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38 : Create the 2-D Laplacian
39 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40 :
41 3 : PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
42 3 : PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
43 3 : PetscCall(MatSetFromOptions(A));
44 3 : PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
45 303 : for (II=Istart;II<Iend;II++) {
46 300 : i = II/n; j = II-i*n;
47 300 : if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
48 300 : if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
49 300 : if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
50 300 : if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
51 300 : PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
52 : }
53 3 : PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
54 3 : PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
55 :
56 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57 : Create the eigensolver and set various options
58 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59 :
60 3 : PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
61 3 : PetscCall(EPSSetOperators(eps,A,NULL));
62 3 : PetscCall(EPSSetProblemType(eps,EPS_HEP));
63 3 : PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));
64 3 : PetscCall(EPSSetWhichEigenpairs(eps,EPS_ALL));
65 3 : PetscCall(EPSSetInterval(eps,0.5,1.3));
66 3 : PetscCall(EPSGetST(eps,&st));
67 3 : PetscCall(STSetType(st,STFILTER));
68 3 : PetscCall(EPSSetFromOptions(eps));
69 :
70 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71 : Solve the problem and display the solution
72 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73 :
74 3 : PetscCall(EPSSolve(eps));
75 :
76 : /* print filter information */
77 3 : PetscCall(PetscObjectTypeCompare((PetscObject)st,STFILTER,&flag));
78 3 : if (flag) {
79 3 : PetscCall(STFilterGetDegree(st,°ree));
80 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Filter degree: %" PetscInt_FMT "\n",degree));
81 3 : PetscCall(STFilterGetInterval(st,&inta,&intb));
82 3 : PetscCall(STFilterGetRange(st,&rleft,&rright));
83 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Requested interval: [%g,%g], range: [%g,%g]\n\n",(double)inta,(double)intb,(double)rleft,(double)rright));
84 : }
85 :
86 3 : PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
87 3 : if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
88 : else {
89 0 : PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
90 0 : PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
91 0 : PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
92 0 : PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
93 : }
94 :
95 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96 : Solve the problem again after changing the matrix
97 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98 3 : if (modify) {
99 2 : PetscCall(MatSetValue(A,0,0,0.3,INSERT_VALUES));
100 2 : PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
101 2 : PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
102 2 : PetscCall(EPSSetOperators(eps,A,NULL));
103 2 : PetscCall(EPSSolve(eps));
104 2 : PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
105 2 : if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
106 : else {
107 0 : PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
108 0 : PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
109 0 : PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
110 0 : PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
111 : }
112 : }
113 :
114 3 : PetscCall(EPSDestroy(&eps));
115 3 : PetscCall(MatDestroy(&A));
116 3 : PetscCall(SlepcFinalize());
117 : return 0;
118 : }
119 :
120 : /*TEST
121 :
122 : test:
123 : suffix: 1
124 : args: -terse
125 : filter: sed -e "s/0.161982,7.83797/0.162007,7.83897/"
126 : requires: !single
127 :
128 : test:
129 : suffix: 2
130 : args: -modify -st_filter_range -0.5,8 -terse
131 : requires: !single
132 :
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
138 :
139 : TEST*/
|