Actual source code: test21.c
slepc-3.18.1 2022-11-02
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[] = "Illustrates region filtering. "
12: "Based on ex5.\n"
13: "The command line options are:\n"
14: " -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";
16: #include <slepceps.h>
18: /*
19: User-defined routines
20: */
21: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
23: int main(int argc,char **argv)
24: {
25: Mat A;
26: EPS eps;
27: ST st;
28: RG rg;
29: PetscReal radius,tol=PETSC_SMALL;
30: PetscScalar target=0.5,kr,ki;
31: PetscInt N,m=15,nev,i,nconv;
32: PetscBool checkfile;
33: char filename[PETSC_MAX_PATH_LEN];
34: PetscViewer viewer;
35: char str[50];
38: SlepcInitialize(&argc,&argv,(char*)0,help);
39: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
40: N = m*(m+1)/2;
41: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n",N,m);
42: PetscOptionsGetScalar(NULL,NULL,"-target",&target,NULL);
43: SlepcSNPrintfScalar(str,sizeof(str),target,PETSC_FALSE);
44: PetscPrintf(PETSC_COMM_WORLD,"Searching closest eigenvalues to the right of %s.\n\n",str);
46: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: Compute the operator matrix that defines the eigensystem, Ax=kx
48: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50: MatCreate(PETSC_COMM_WORLD,&A);
51: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
52: MatSetFromOptions(A);
53: MatSetUp(A);
54: MatMarkovModel(m,A);
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Create a standalone spectral transformation
58: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60: STCreate(PETSC_COMM_WORLD,&st);
61: STSetType(st,STSINVERT);
62: STSetShift(st,target);
64: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65: Create a region for filtering
66: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68: RGCreate(PETSC_COMM_WORLD,&rg);
69: RGSetType(rg,RGELLIPSE);
70: radius = (1.0-PetscRealPart(target))/2.0;
71: RGEllipseSetParameters(rg,target+radius,radius,0.1);
73: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74: Create the eigensolver and set various options
75: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: EPSCreate(PETSC_COMM_WORLD,&eps);
78: EPSSetST(eps,st);
79: EPSSetRG(eps,rg);
80: EPSSetOperators(eps,A,NULL);
81: EPSSetProblemType(eps,EPS_NHEP);
82: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
83: EPSSetTarget(eps,target);
84: EPSSetFromOptions(eps);
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Solve the eigensystem and display solution
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: EPSSolve(eps);
91: EPSGetDimensions(eps,&nev,NULL,NULL);
92: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
93: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
95: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96: Check file containing the eigenvalues
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: PetscOptionsGetString(NULL,NULL,"-checkfile",filename,sizeof(filename),&checkfile);
99: if (checkfile) {
100: #if defined(PETSC_HAVE_COMPLEX)
101: PetscComplex *eigs,eval;
102: EPSGetConverged(eps,&nconv);
103: PetscMalloc1(nconv,&eigs);
104: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
105: PetscViewerBinaryRead(viewer,eigs,nconv,NULL,PETSC_COMPLEX);
106: PetscViewerDestroy(&viewer);
107: for (i=0;i<nconv;i++) {
108: EPSGetEigenpair(eps,i,&kr,&ki,NULL,NULL);
109: #if defined(PETSC_USE_COMPLEX)
110: eval = kr;
111: #else
112: eval = PetscCMPLX(kr,ki);
113: #endif
115: }
116: PetscFree(eigs);
117: #else
118: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"The -checkfile option requires C99 complex numbers");
119: #endif
120: }
122: EPSDestroy(&eps);
123: STDestroy(&st);
124: RGDestroy(&rg);
125: MatDestroy(&A);
126: SlepcFinalize();
127: return 0;
128: }
130: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
131: {
132: const PetscReal cst = 0.5/(PetscReal)(m-1);
133: PetscReal pd,pu;
134: PetscInt Istart,Iend,i,j,jmax,ix=0;
137: MatGetOwnershipRange(A,&Istart,&Iend);
138: for (i=1;i<=m;i++) {
139: jmax = m-i+1;
140: for (j=1;j<=jmax;j++) {
141: ix = ix + 1;
142: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
143: if (j!=jmax) {
144: pd = cst*(PetscReal)(i+j-1);
145: /* north */
146: if (i==1) MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
147: else MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
148: /* east */
149: if (j==1) MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
150: else MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
151: }
152: /* south */
153: pu = 0.5 - cst*(PetscReal)(i+j-3);
154: if (j>1) MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
155: /* west */
156: if (i>1) MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
157: }
158: }
159: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
160: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
161: return 0;
162: }
164: /*TEST
166: test:
167: suffix: 1
168: args: -eps_nev 4 -eps_ncv 20 -eps_view_values binary:myvalues.bin -checkfile myvalues.bin
169: output_file: output/test11_1.out
170: requires: !single c99_complex
172: TEST*/