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 two-sided Krylov-Schur without calling EPSSetFromOptions (based on ex5.c).\n\n"
12 : "The command line options are:\n"
13 : " -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";
14 :
15 : #include <slepceps.h>
16 :
17 : /*
18 : User-defined routines
19 : */
20 : PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
21 :
22 1 : int main(int argc,char **argv)
23 : {
24 1 : Mat A; /* operator matrix */
25 1 : EPS eps; /* eigenproblem solver context */
26 1 : PetscReal tol=1000*PETSC_MACHINE_EPSILON;
27 1 : PetscInt N,m=15,nev;
28 :
29 1 : PetscFunctionBeginUser;
30 1 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
31 :
32 1 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
33 1 : N = m*(m+1)/2;
34 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n\n",N,m));
35 :
36 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37 : Compute the operator matrix that defines the eigensystem, Ax=kx
38 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39 :
40 1 : PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
41 1 : PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
42 1 : PetscCall(MatSetFromOptions(A));
43 1 : PetscCall(MatMarkovModel(m,A));
44 :
45 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46 : Create the eigensolver and set various options
47 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
48 :
49 1 : PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
50 1 : PetscCall(EPSSetOperators(eps,A,NULL));
51 1 : PetscCall(EPSSetProblemType(eps,EPS_NHEP));
52 1 : PetscCall(EPSSetTolerances(eps,tol,PETSC_CURRENT));
53 1 : PetscCall(EPSSetDimensions(eps,4,PETSC_DETERMINE,PETSC_DETERMINE));
54 1 : PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));
55 1 : PetscCall(EPSSetTwoSided(eps,PETSC_TRUE));
56 :
57 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58 : Solve the eigensystem
59 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60 :
61 1 : PetscCall(EPSSolve(eps));
62 1 : PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
63 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
64 :
65 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66 : Display solution and clean up
67 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68 :
69 1 : PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
70 1 : PetscCall(EPSDestroy(&eps));
71 1 : PetscCall(MatDestroy(&A));
72 1 : PetscCall(SlepcFinalize());
73 : return 0;
74 : }
75 :
76 1 : PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
77 : {
78 1 : const PetscReal cst = 0.5/(PetscReal)(m-1);
79 1 : PetscReal pd,pu;
80 1 : PetscInt Istart,Iend,i,j,jmax,ix=0;
81 :
82 1 : PetscFunctionBeginUser;
83 1 : PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
84 16 : for (i=1;i<=m;i++) {
85 15 : jmax = m-i+1;
86 135 : for (j=1;j<=jmax;j++) {
87 120 : ix = ix + 1;
88 120 : if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
89 120 : if (j!=jmax) {
90 105 : pd = cst*(PetscReal)(i+j-1);
91 : /* north */
92 105 : if (i==1) PetscCall(MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES));
93 91 : else PetscCall(MatSetValue(A,ix-1,ix,pd,INSERT_VALUES));
94 : /* east */
95 105 : if (j==1) PetscCall(MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES));
96 91 : else PetscCall(MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES));
97 : }
98 : /* south */
99 120 : pu = 0.5 - cst*(PetscReal)(i+j-3);
100 120 : if (j>1) PetscCall(MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES));
101 : /* west */
102 120 : if (i>1) PetscCall(MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES));
103 : }
104 : }
105 1 : PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
106 1 : PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
107 1 : PetscFunctionReturn(PETSC_SUCCESS);
108 : }
109 :
110 : /*TEST
111 :
112 : test:
113 : requires: !single
114 :
115 : TEST*/
|