Actual source code: test40.c
slepc-3.17.2 2022-08-09
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 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";
15: #include <slepceps.h>
17: /*
18: User-defined routines
19: */
20: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
22: int main(int argc,char **argv)
23: {
24: Mat A; /* operator matrix */
25: EPS eps; /* eigenproblem solver context */
26: PetscReal tol=1000*PETSC_MACHINE_EPSILON;
27: PetscInt N,m=15,nev;
29: SlepcInitialize(&argc,&argv,(char*)0,help);
31: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
32: N = m*(m+1)/2;
33: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n\n",N,m);
35: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: Compute the operator matrix that defines the eigensystem, Ax=kx
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39: MatCreate(PETSC_COMM_WORLD,&A);
40: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
41: MatSetFromOptions(A);
42: MatSetUp(A);
43: MatMarkovModel(m,A);
45: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46: Create the eigensolver and set various options
47: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
49: EPSCreate(PETSC_COMM_WORLD,&eps);
50: EPSSetOperators(eps,A,NULL);
51: EPSSetProblemType(eps,EPS_NHEP);
52: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
53: EPSSetDimensions(eps,4,PETSC_DEFAULT,PETSC_DEFAULT);
54: EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
55: EPSSetTwoSided(eps,PETSC_TRUE);
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Solve the eigensystem
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: EPSSolve(eps);
62: EPSGetDimensions(eps,&nev,NULL,NULL);
63: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
65: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66: Display solution and clean up
67: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
70: EPSDestroy(&eps);
71: MatDestroy(&A);
72: SlepcFinalize();
73: return 0;
74: }
76: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
77: {
78: const PetscReal cst = 0.5/(PetscReal)(m-1);
79: PetscReal pd,pu;
80: PetscInt Istart,Iend,i,j,jmax,ix=0;
83: MatGetOwnershipRange(A,&Istart,&Iend);
84: for (i=1;i<=m;i++) {
85: jmax = m-i+1;
86: for (j=1;j<=jmax;j++) {
87: ix = ix + 1;
88: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
89: if (j!=jmax) {
90: pd = cst*(PetscReal)(i+j-1);
91: /* north */
92: if (i==1) MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
93: else MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
94: /* east */
95: if (j==1) MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
96: else MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
97: }
98: /* south */
99: pu = 0.5 - cst*(PetscReal)(i+j-3);
100: if (j>1) MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
101: /* west */
102: if (i>1) MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
103: }
104: }
105: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
106: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
107: PetscFunctionReturn(0);
108: }
110: /*TEST
112: test:
113: requires: !single
115: TEST*/