Actual source code: test11.c
slepc-3.20.2 2024-03-15
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[] = "Solves the same problem as in ex5, but with a user-defined sorting criterion."
12: "It is a standard nonsymmetric eigenproblem with real eigenvalues and the rightmost eigenvalue is known to be 1.\n"
13: "This example illustrates how the user can set a custom spectrum selection.\n\n"
14: "The command line options are:\n"
15: " -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";
17: #include <slepceps.h>
19: /*
20: User-defined routines
21: */
23: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx);
24: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
26: int main(int argc,char **argv)
27: {
28: Vec v0; /* initial vector */
29: Mat A; /* operator matrix */
30: EPS eps; /* eigenproblem solver context */
31: ST st; /* spectral transformation associated */
32: PetscReal tol=PETSC_SMALL;
33: PetscScalar target=0.5;
34: PetscInt N,m=15,nev;
35: char str[50];
37: PetscFunctionBeginUser;
38: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
40: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
41: N = m*(m+1)/2;
42: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n",N,m));
43: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-target",&target,NULL));
44: PetscCall(SlepcSNPrintfScalar(str,sizeof(str),target,PETSC_FALSE));
45: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Searching closest eigenvalues to the right of %s.\n\n",str));
47: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: Compute the operator matrix that defines the eigensystem, Ax=kx
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
52: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
53: PetscCall(MatSetFromOptions(A));
54: PetscCall(MatSetUp(A));
55: PetscCall(MatMarkovModel(m,A));
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Create the eigensolver and set various options
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: /*
62: Create eigensolver context
63: */
64: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
66: /*
67: Set operators. In this case, it is a standard eigenvalue problem
68: */
69: PetscCall(EPSSetOperators(eps,A,NULL));
70: PetscCall(EPSSetProblemType(eps,EPS_NHEP));
71: PetscCall(EPSSetTolerances(eps,tol,PETSC_DEFAULT));
73: /*
74: Set the custom comparing routine in order to obtain the eigenvalues
75: closest to the target on the right only
76: */
77: PetscCall(EPSSetEigenvalueComparison(eps,MyEigenSort,&target));
79: /*
80: Set solver parameters at runtime
81: */
82: PetscCall(EPSSetFromOptions(eps));
84: /*
85: Set the preconditioner based on A - target * I
86: */
87: PetscCall(EPSGetST(eps,&st));
88: PetscCall(STSetShift(st,target));
90: /*
91: Set the initial vector. This is optional, if not done the initial
92: vector is set to random values
93: */
94: PetscCall(MatCreateVecs(A,&v0,NULL));
95: PetscCall(VecSet(v0,1.0));
96: PetscCall(EPSSetInitialSpace(eps,1,&v0));
98: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99: Solve the eigensystem
100: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102: PetscCall(EPSSolve(eps));
103: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
104: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
106: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107: Display solution and clean up
108: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
111: PetscCall(EPSDestroy(&eps));
112: PetscCall(MatDestroy(&A));
113: PetscCall(VecDestroy(&v0));
114: PetscCall(SlepcFinalize());
115: return 0;
116: }
118: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
119: {
120: const PetscReal cst = 0.5/(PetscReal)(m-1);
121: PetscReal pd,pu;
122: PetscInt Istart,Iend,i,j,jmax,ix=0;
124: PetscFunctionBeginUser;
125: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
126: for (i=1;i<=m;i++) {
127: jmax = m-i+1;
128: for (j=1;j<=jmax;j++) {
129: ix = ix + 1;
130: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
131: if (j!=jmax) {
132: pd = cst*(PetscReal)(i+j-1);
133: /* north */
134: if (i==1) PetscCall(MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES));
135: else PetscCall(MatSetValue(A,ix-1,ix,pd,INSERT_VALUES));
136: /* east */
137: if (j==1) PetscCall(MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES));
138: else PetscCall(MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES));
139: }
140: /* south */
141: pu = 0.5 - cst*(PetscReal)(i+j-3);
142: if (j>1) PetscCall(MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES));
143: /* west */
144: if (i>1) PetscCall(MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES));
145: }
146: }
147: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
148: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
149: PetscFunctionReturn(PETSC_SUCCESS);
150: }
152: /*
153: Function for user-defined eigenvalue ordering criterion.
155: Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
156: one of them as the preferred one according to the criterion.
157: In this example, the preferred value is the one closest to the target,
158: but on the right side.
159: */
160: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
161: {
162: PetscScalar target = *(PetscScalar*)ctx;
163: PetscReal da,db;
164: PetscBool aisright,bisright;
166: PetscFunctionBeginUser;
167: if (PetscRealPart(target) < PetscRealPart(ar)) aisright = PETSC_TRUE;
168: else aisright = PETSC_FALSE;
169: if (PetscRealPart(target) < PetscRealPart(br)) bisright = PETSC_TRUE;
170: else bisright = PETSC_FALSE;
171: if (aisright == bisright) {
172: /* both are on the same side of the target */
173: da = SlepcAbsEigenvalue(ar-target,ai);
174: db = SlepcAbsEigenvalue(br-target,bi);
175: if (da < db) *r = -1;
176: else if (da > db) *r = 1;
177: else *r = 0;
178: } else if (aisright && !bisright) *r = -1; /* 'a' is on the right */
179: else *r = 1; /* 'b' is on the right */
180: PetscFunctionReturn(PETSC_SUCCESS);
181: }
183: /*TEST
185: testset:
186: args: -eps_nev 4
187: requires: !single
188: output_file: output/test11_1.out
189: test:
190: suffix: 1
191: args: -eps_type {{krylovschur arnoldi lapack}} -st_type sinvert
192: test:
193: suffix: 1_ks_cayley
194: args: -st_type cayley -st_cayley_antishift 1
196: test:
197: suffix: 2
198: args: -target 0.77 -eps_type gd -eps_nev 4 -eps_tol 1e-7 -eps_gd_krylov_start -eps_gd_blocksize 3
199: requires: double
200: filter: sed -e "s/[+-]0\.0*i//g"
202: TEST*/