Subversion Repositories slepc-dev

Rev

Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1376 slepc 1
/*
2
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 3
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2575 eromero 4
   Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
6 dsic.upv.es!jroman 5
 
1672 slepc 6
   This file is part of SLEPc.
7
 
8
   SLEPc is free software: you can redistribute it and/or modify it under  the
9
   terms of version 3 of the GNU Lesser General Public License as published by
10
   the Free Software Foundation.
11
 
12
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
13
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
14
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
15
   more details.
16
 
17
   You  should have received a copy of the GNU Lesser General  Public  License
18
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
1376 slepc 19
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20
*/
21
 
161 dsic.upv.es!antodo 22
static char help[] = "Eigenvalue problem associated with a Markov model of a random walk on a triangular grid. "
1156 slepc 23
  "It is a standard nonsymmetric eigenproblem with real eigenvalues and the rightmost eigenvalue is known to be 1.\n"
6 dsic.upv.es!jroman 24
  "This example illustrates how the user can set the initial vector.\n\n"
1156 slepc 25
  "The command line options are:\n"
6 dsic.upv.es!jroman 26
  "  -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";
2342 jroman 27
 
2283 jroman 28
#include <slepceps.h>
6 dsic.upv.es!jroman 29
 
30
/*
31
   User-defined routines
32
*/
2331 jroman 33
PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
6 dsic.upv.es!jroman 34
 
35
#undef __FUNCT__
36
#define __FUNCT__ "main"
2331 jroman 37
int main(int argc,char **argv)
6 dsic.upv.es!jroman 38
{
2316 jroman 39
  Vec            v0;              /* initial vector */
40
  Mat            A;               /* operator matrix */
41
  EPS            eps;             /* eigenproblem solver context */
42
  const EPSType  type;
2503 jroman 43
  PetscReal      tol;
44
  PetscInt       N,m=15,nev,maxit,its;
2317 jroman 45
  PetscErrorCode ierr;
983 slepc 46
 
6 dsic.upv.es!jroman 47
  SlepcInitialize(&argc,&argv,(char*)0,help);
48
 
49
  ierr = PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);CHKERRQ(ierr);
50
  N = m*(m+1)/2;
2503 jroman 51
  ierr = PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n\n",N,m);CHKERRQ(ierr);
6 dsic.upv.es!jroman 52
 
53
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54
     Compute the operator matrix that defines the eigensystem, Ax=kx
55
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56
 
828 dsic.upv.es!antodo 57
  ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
58
  ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr);
6 dsic.upv.es!jroman 59
  ierr = MatSetFromOptions(A);CHKERRQ(ierr);
2331 jroman 60
  ierr = MatMarkovModel(m,A);CHKERRQ(ierr);
6 dsic.upv.es!jroman 61
 
62
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63
                Create the eigensolver and set various options
64
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65
 
66
  /*
67
     Create eigensolver context
68
  */
69
  ierr = EPSCreate(PETSC_COMM_WORLD,&eps);CHKERRQ(ierr);
70
 
71
  /*
72
     Set operators. In this case, it is a standard eigenvalue problem
73
  */
74
  ierr = EPSSetOperators(eps,A,PETSC_NULL);CHKERRQ(ierr);
408 dsic.upv.es!antodo 75
  ierr = EPSSetProblemType(eps,EPS_NHEP);CHKERRQ(ierr);
6 dsic.upv.es!jroman 76
 
77
  /*
78
     Set solver parameters at runtime
79
  */
80
  ierr = EPSSetFromOptions(eps);CHKERRQ(ierr);
81
 
82
  /*
83
     Set the initial vector. This is optional, if not done the initial
84
     vector is set to random values
85
  */
583 dsic.upv.es!jroman 86
  ierr = MatGetVecs(A,&v0,PETSC_NULL);CHKERRQ(ierr);
828 dsic.upv.es!antodo 87
  ierr = VecSet(v0,1.0);CHKERRQ(ierr);
1933 jroman 88
  ierr = EPSSetInitialSpace(eps,1,&v0);CHKERRQ(ierr);
6 dsic.upv.es!jroman 89
 
90
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91
                      Solve the eigensystem
92
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93
 
78 dsic.upv.es!antodo 94
  ierr = EPSSolve(eps);CHKERRQ(ierr);
2331 jroman 95
  ierr = EPSGetIterationNumber(eps,&its);CHKERRQ(ierr);
2503 jroman 96
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);CHKERRQ(ierr);
6 dsic.upv.es!jroman 97
 
98
  /*
99
     Optional: Get some information from the solver and display it
100
  */
101
  ierr = EPSGetType(eps,&type);CHKERRQ(ierr);
102
  ierr = PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);CHKERRQ(ierr);
1575 slepc 103
  ierr = EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2503 jroman 104
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);CHKERRQ(ierr);
6 dsic.upv.es!jroman 105
  ierr = EPSGetTolerances(eps,&tol,&maxit);CHKERRQ(ierr);
2503 jroman 106
  ierr = PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4G, maxit=%D\n",tol,maxit);CHKERRQ(ierr);
6 dsic.upv.es!jroman 107
 
108
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109
                    Display solution and clean up
110
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111
 
2503 jroman 112
  ierr = EPSPrintSolution(eps,PETSC_NULL);CHKERRQ(ierr);
2308 jroman 113
  ierr = EPSDestroy(&eps);CHKERRQ(ierr);
2305 jroman 114
  ierr = MatDestroy(&A);CHKERRQ(ierr);
115
  ierr = VecDestroy(&v0);CHKERRQ(ierr);
6 dsic.upv.es!jroman 116
  ierr = SlepcFinalize();CHKERRQ(ierr);
117
  return 0;
118
}
119
 
120
#undef __FUNCT__
121
#define __FUNCT__ "MatMarkovModel"
122
/*
123
    Matrix generator for a Markov model of a random walk on a triangular grid.
124
 
125
    This subroutine generates a test matrix that models a random walk on a
126
    triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a
127
    FORTRAN subroutine to calculate the dominant invariant subspaces of a real
128
    matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
129
    papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
130
    (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
131
    algorithms. The transpose of the matrix  is stochastic and so it is known
132
    that one is an exact eigenvalue. One seeks the eigenvector of the transpose
133
    associated with the eigenvalue unity. The problem is to calculate the steady
134
    state probability distribution of the system, which is the eigevector
135
    associated with the eigenvalue one and scaled in such a way that the sum all
136
    the components is equal to one.
137
 
138
    Note: the code will actually compute the transpose of the stochastic matrix
139
    that contains the transition probabilities.
140
*/
2331 jroman 141
PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
6 dsic.upv.es!jroman 142
{
143
  const PetscReal cst = 0.5/(PetscReal)(m-1);
2331 jroman 144
  PetscReal       pd,pu;
145
  PetscInt        Istart,Iend,i,j,jmax,ix=0;
983 slepc 146
  PetscErrorCode  ierr;
6 dsic.upv.es!jroman 147
 
983 slepc 148
  PetscFunctionBegin;
6 dsic.upv.es!jroman 149
  ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);
2331 jroman 150
  for (i=1;i<=m;i++) {
6 dsic.upv.es!jroman 151
    jmax = m-i+1;
2331 jroman 152
    for (j=1;j<=jmax;j++) {
6 dsic.upv.es!jroman 153
      ix = ix + 1;
2331 jroman 154
      if (ix-1<Istart || ix>Iend) continue;  /* compute only owned rows */
155
      if (j!=jmax) {
6 dsic.upv.es!jroman 156
        pd = cst*(PetscReal)(i+j-1);
157
        /* north */
2331 jroman 158
        if (i==1) {
159
          ierr = MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);CHKERRQ(ierr);
2316 jroman 160
        } else {
2331 jroman 161
          ierr = MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);CHKERRQ(ierr);
6 dsic.upv.es!jroman 162
        }
163
        /* east */
2331 jroman 164
        if (j==1) {
165
          ierr = MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);CHKERRQ(ierr);
2316 jroman 166
        } else {
2331 jroman 167
          ierr = MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);CHKERRQ(ierr);
6 dsic.upv.es!jroman 168
        }
169
      }
170
      /* south */
171
      pu = 0.5 - cst*(PetscReal)(i+j-3);
2331 jroman 172
      if (j>1) {
173
        ierr = MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);CHKERRQ(ierr);
6 dsic.upv.es!jroman 174
      }
175
      /* west */
2331 jroman 176
      if (i>1) {
177
        ierr = MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);CHKERRQ(ierr);
6 dsic.upv.es!jroman 178
      }
179
    }
180
  }
2331 jroman 181
  ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
182
  ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
983 slepc 183
  PetscFunctionReturn(0);
6 dsic.upv.es!jroman 184
}
185