Actual source code: test21.c
slepc-3.17.1 2022-04-11
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: PetscComplex *eigs,eval;
32: PetscInt N,m=15,nev,i,nconv;
33: PetscBool checkfile;
34: char filename[PETSC_MAX_PATH_LEN];
35: PetscViewer viewer;
36: char str[50];
38: SlepcInitialize(&argc,&argv,(char*)0,help);
39: #if defined(PETSC_HAVE_COMPLEX)
40: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
41: N = m*(m+1)/2;
42: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n",N,m);
43: PetscOptionsGetScalar(NULL,NULL,"-target",&target,NULL);
44: SlepcSNPrintfScalar(str,sizeof(str),target,PETSC_FALSE);
45: 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: MatCreate(PETSC_COMM_WORLD,&A);
52: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
53: MatSetFromOptions(A);
54: MatSetUp(A);
55: MatMarkovModel(m,A);
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Create a standalone spectral transformation
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: STCreate(PETSC_COMM_WORLD,&st);
62: STSetType(st,STSINVERT);
63: STSetShift(st,target);
65: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66: Create a region for filtering
67: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69: RGCreate(PETSC_COMM_WORLD,&rg);
70: RGSetType(rg,RGELLIPSE);
71: radius = (1.0-PetscRealPart(target))/2.0;
72: RGEllipseSetParameters(rg,target+radius,radius,0.1);
74: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75: Create the eigensolver and set various options
76: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78: EPSCreate(PETSC_COMM_WORLD,&eps);
79: EPSSetST(eps,st);
80: EPSSetRG(eps,rg);
81: EPSSetOperators(eps,A,NULL);
82: EPSSetProblemType(eps,EPS_NHEP);
83: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
84: EPSSetTarget(eps,target);
85: EPSSetFromOptions(eps);
87: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: Solve the eigensystem and display solution
89: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91: EPSSolve(eps);
92: EPSGetDimensions(eps,&nev,NULL,NULL);
93: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
94: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Check file containing the eigenvalues
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: PetscOptionsGetString(NULL,NULL,"-checkfile",filename,sizeof(filename),&checkfile);
100: if (checkfile) {
101: EPSGetConverged(eps,&nconv);
102: PetscMalloc1(nconv,&eigs);
103: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
104: PetscViewerBinaryRead(viewer,eigs,nconv,NULL,PETSC_COMPLEX);
105: PetscViewerDestroy(&viewer);
106: for (i=0;i<nconv;i++) {
107: EPSGetEigenpair(eps,i,&kr,&ki,NULL,NULL);
108: #if defined(PETSC_USE_COMPLEX)
109: eval = kr;
110: #else
111: eval = PetscCMPLX(kr,ki);
112: #endif
114: }
115: PetscFree(eigs);
116: }
118: EPSDestroy(&eps);
119: STDestroy(&st);
120: RGDestroy(&rg);
121: MatDestroy(&A);
122: #else
123: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This example requires C99 complex numbers");
124: #endif
125: SlepcFinalize();
126: return 0;
127: }
129: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
130: {
131: const PetscReal cst = 0.5/(PetscReal)(m-1);
132: PetscReal pd,pu;
133: PetscInt Istart,Iend,i,j,jmax,ix=0;
136: MatGetOwnershipRange(A,&Istart,&Iend);
137: for (i=1;i<=m;i++) {
138: jmax = m-i+1;
139: for (j=1;j<=jmax;j++) {
140: ix = ix + 1;
141: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
142: if (j!=jmax) {
143: pd = cst*(PetscReal)(i+j-1);
144: /* north */
145: if (i==1) MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
146: else MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
147: /* east */
148: if (j==1) MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
149: else MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
150: }
151: /* south */
152: pu = 0.5 - cst*(PetscReal)(i+j-3);
153: if (j>1) MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
154: /* west */
155: if (i>1) MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
156: }
157: }
158: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
159: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
160: PetscFunctionReturn(0);
161: }
163: /*TEST
165: test:
166: suffix: 1
167: args: -eps_nev 4 -eps_ncv 20 -eps_view_values binary:myvalues.bin -checkfile myvalues.bin
168: output_file: output/test11_1.out
169: requires: !single c99_complex
171: TEST*/