Actual source code: test21.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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=1000*PETSC_MACHINE_EPSILON;
30: PetscScalar target=0.5;
31: PetscInt N,m=15,nev;
33: char str[50];
35: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
36: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
37: N = m*(m+1)/2;
38: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n",N,m);
39: PetscOptionsGetScalar(NULL,NULL,"-target",&target,NULL);
40: SlepcSNPrintfScalar(str,50,target,PETSC_FALSE);
41: PetscPrintf(PETSC_COMM_WORLD,"Searching closest eigenvalues to the right of %s.\n\n",str);
43: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44: Compute the operator matrix that defines the eigensystem, Ax=kx
45: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
47: MatCreate(PETSC_COMM_WORLD,&A);
48: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
49: MatSetFromOptions(A);
50: MatSetUp(A);
51: MatMarkovModel(m,A);
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Create a standalone spectral transformation
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: STCreate(PETSC_COMM_WORLD,&st);
58: STSetType(st,STSINVERT);
59: STSetShift(st,target);
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Create a region for filtering
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65: RGCreate(PETSC_COMM_WORLD,&rg);
66: RGSetType(rg,RGELLIPSE);
67: radius = (1.0-PetscRealPart(target))/2.0;
68: RGEllipseSetParameters(rg,target+radius,radius,0.1);
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Create the eigensolver and set various options
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74: EPSCreate(PETSC_COMM_WORLD,&eps);
75: EPSSetST(eps,st);
76: EPSSetRG(eps,rg);
77: EPSSetOperators(eps,A,NULL);
78: EPSSetProblemType(eps,EPS_NHEP);
79: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
80: EPSSetTarget(eps,target);
81: EPSSetFromOptions(eps);
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Solve the eigensystem
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: EPSSolve(eps);
88: EPSGetDimensions(eps,&nev,NULL,NULL);
89: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
91: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92: Display solution and clean up
93: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
96: EPSDestroy(&eps);
97: STDestroy(&st);
98: RGDestroy(&rg);
99: MatDestroy(&A);
100: SlepcFinalize();
101: return ierr;
102: }
104: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
105: {
106: const PetscReal cst = 0.5/(PetscReal)(m-1);
107: PetscReal pd,pu;
108: PetscInt Istart,Iend,i,j,jmax,ix=0;
109: PetscErrorCode ierr;
112: MatGetOwnershipRange(A,&Istart,&Iend);
113: for (i=1;i<=m;i++) {
114: jmax = m-i+1;
115: for (j=1;j<=jmax;j++) {
116: ix = ix + 1;
117: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
118: if (j!=jmax) {
119: pd = cst*(PetscReal)(i+j-1);
120: /* north */
121: if (i==1) {
122: MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
123: } else {
124: MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
125: }
126: /* east */
127: if (j==1) {
128: MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
129: } else {
130: MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
131: }
132: }
133: /* south */
134: pu = 0.5 - cst*(PetscReal)(i+j-3);
135: if (j>1) {
136: MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
137: }
138: /* west */
139: if (i>1) {
140: MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
141: }
142: }
143: }
144: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
145: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
146: return(0);
147: }
149: /*TEST
151: test:
152: suffix: 1
153: args: -eps_nev 4 -eps_ncv 20
154: output_file: output/test11_1.out
156: TEST*/