Actual source code: test40.c
slepc-3.19.2 2023-09-05
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: PetscFunctionBeginUser;
30: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
32: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
33: N = m*(m+1)/2;
34: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n\n",N,m));
36: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37: Compute the operator matrix that defines the eigensystem, Ax=kx
38: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
41: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
42: PetscCall(MatSetFromOptions(A));
43: PetscCall(MatSetUp(A));
44: PetscCall(MatMarkovModel(m,A));
46: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: Create the eigensolver and set various options
48: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
51: PetscCall(EPSSetOperators(eps,A,NULL));
52: PetscCall(EPSSetProblemType(eps,EPS_NHEP));
53: PetscCall(EPSSetTolerances(eps,tol,PETSC_DEFAULT));
54: PetscCall(EPSSetDimensions(eps,4,PETSC_DEFAULT,PETSC_DEFAULT));
55: PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));
56: PetscCall(EPSSetTwoSided(eps,PETSC_TRUE));
58: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59: Solve the eigensystem
60: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
62: PetscCall(EPSSolve(eps));
63: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
64: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
66: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: Display solution and clean up
68: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
70: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
71: PetscCall(EPSDestroy(&eps));
72: PetscCall(MatDestroy(&A));
73: PetscCall(SlepcFinalize());
74: return 0;
75: }
77: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
78: {
79: const PetscReal cst = 0.5/(PetscReal)(m-1);
80: PetscReal pd,pu;
81: PetscInt Istart,Iend,i,j,jmax,ix=0;
83: PetscFunctionBeginUser;
84: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
85: for (i=1;i<=m;i++) {
86: jmax = m-i+1;
87: for (j=1;j<=jmax;j++) {
88: ix = ix + 1;
89: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
90: if (j!=jmax) {
91: pd = cst*(PetscReal)(i+j-1);
92: /* north */
93: if (i==1) PetscCall(MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES));
94: else PetscCall(MatSetValue(A,ix-1,ix,pd,INSERT_VALUES));
95: /* east */
96: if (j==1) PetscCall(MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES));
97: else PetscCall(MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES));
98: }
99: /* south */
100: pu = 0.5 - cst*(PetscReal)(i+j-3);
101: if (j>1) PetscCall(MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES));
102: /* west */
103: if (i>1) PetscCall(MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES));
104: }
105: }
106: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
107: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
108: PetscFunctionReturn(PETSC_SUCCESS);
109: }
111: /*TEST
113: test:
114: requires: !single
116: TEST*/