Subversion Repositories slepc-dev

Rev

Go to most recent revision | Only display areas with differences | Ignore whitespace | Details | Blame | Last modification | View Log | RSS feed

Rev 1961 Rev 2116
/*
/*
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
   SLEPc - Scalable Library for Eigenvalue Problem Computations
   SLEPc - Scalable Library for Eigenvalue Problem Computations
   Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
   Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
 
 
   This file is part of SLEPc.
   This file is part of SLEPc.
     
     
   SLEPc is free software: you can redistribute it and/or modify it under  the
   SLEPc is free software: you can redistribute it and/or modify it under  the
   terms of version 3 of the GNU Lesser General Public License as published by
   terms of version 3 of the GNU Lesser General Public License as published by
   the Free Software Foundation.
   the Free Software Foundation.
 
 
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
   more details.
   more details.
 
 
   You  should have received a copy of the GNU Lesser General  Public  License
   You  should have received a copy of the GNU Lesser General  Public  License
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
*/
*/
 
 
static char help[] = "Standard symmetric eigenproblem corresponding to the Laplacian operator in 2 dimensions.\n\n"
static char help[] = "Standard symmetric eigenproblem corresponding to the Laplacian operator in 2 dimensions.\n\n"
  "The command line options are:\n"
  "The command line options are:\n"
  "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
  "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
  "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
  "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
 
 
#include "slepceps.h"
#include "slepceps.h"
 
 
#undef __FUNCT__
#undef __FUNCT__
#define __FUNCT__ "main"
#define __FUNCT__ "main"
int main( int argc, char **argv )
int main( int argc, char **argv )
{
{
  Mat            A;               /* operator matrix */
  Mat            A;               /* operator matrix */
  EPS            eps;             /* eigenproblem solver context */
  EPS            eps;             /* eigenproblem solver context */
  const EPSType  type;
  const EPSType  type;
  PetscReal      error, tol, re, im;
  PetscReal      error, tol, re, im;
  PetscScalar    kr, ki;
  PetscScalar    kr, ki;
  PetscErrorCode ierr;
  PetscErrorCode ierr;
  PetscInt       N, n=10, m, Istart, Iend, II, nev, maxit, i, j, its, nconv;
  PetscInt       N, n=10, m, Istart, Iend, II, nev, maxit, i, j, its, nconv;
  PetscTruth     flag;
  PetscTruth     flag;
 
 
  SlepcInitialize(&argc,&argv,(char*)0,help);
  SlepcInitialize(&argc,&argv,(char*)0,help);
 
 
  ierr = PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);CHKERRQ(ierr);
  ierr = PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);CHKERRQ(ierr);
  ierr = PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flag);CHKERRQ(ierr);
  ierr = PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flag);CHKERRQ(ierr);
  if( flag==PETSC_FALSE ) m=n;
  if( flag==PETSC_FALSE ) m=n;
  N = n*m;
  N = n*m;
  ierr = PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%d (%dx%d grid)\n\n",N,n,m);CHKERRQ(ierr);
  ierr = PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%d (%dx%d grid)\n\n",N,n,m);CHKERRQ(ierr);
 
 
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
     Compute the operator matrix that defines the eigensystem, Ax=kx
     Compute the operator matrix that defines the eigensystem, Ax=kx
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 
 
  ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
  ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
  ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr);
  ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr);
  ierr = MatSetFromOptions(A);CHKERRQ(ierr);
  ierr = MatSetFromOptions(A);CHKERRQ(ierr);
 
 
  ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);
  ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);
  for( II=Istart; II<Iend; II++ ) {
  for( II=Istart; II<Iend; II++ ) {
    i = II/n; j = II-i*n;  
    i = II/n; j = II-i*n;  
    if(i>0) { ierr = MatSetValue(A,II,II-n,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
    if(i>0) { ierr = MatSetValue(A,II,II-n,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
    if(i<m-1) { ierr = MatSetValue(A,II,II+n,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
    if(i<m-1) { ierr = MatSetValue(A,II,II+n,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
    if(j>0) { ierr = MatSetValue(A,II,II-1,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
    if(j>0) { ierr = MatSetValue(A,II,II-1,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
    if(j<n-1) { ierr = MatSetValue(A,II,II+1,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
    if(j<n-1) { ierr = MatSetValue(A,II,II+1,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
    ierr = MatSetValue(A,II,II,4.0,INSERT_VALUES);CHKERRQ(ierr);
    ierr = MatSetValue(A,II,II,4.0,INSERT_VALUES);CHKERRQ(ierr);
  }
  }
 
 
  ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
 
 
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
                Create the eigensolver and set various options
                Create the eigensolver and set various options
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 
 
  /*
  /*
     Create eigensolver context
     Create eigensolver context
  */
  */
  ierr = EPSCreate(PETSC_COMM_WORLD,&eps);CHKERRQ(ierr);
  ierr = EPSCreate(PETSC_COMM_WORLD,&eps);CHKERRQ(ierr);
 
 
  /*
  /*
     Set operators. In this case, it is a standard eigenvalue problem
     Set operators. In this case, it is a standard eigenvalue problem
  */
  */
  ierr = EPSSetOperators(eps,A,PETSC_NULL);CHKERRQ(ierr);
  ierr = EPSSetOperators(eps,A,PETSC_NULL);CHKERRQ(ierr);
  ierr = EPSSetProblemType(eps,EPS_HEP);CHKERRQ(ierr);
  ierr = EPSSetProblemType(eps,EPS_HEP);CHKERRQ(ierr);
 
 
  /*
  /*
     Set solver parameters at runtime
     Set solver parameters at runtime
  */
  */
  ierr = EPSSetFromOptions(eps);CHKERRQ(ierr);
  ierr = EPSSetFromOptions(eps);CHKERRQ(ierr);
 
 
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
                      Solve the eigensystem
                      Solve the eigensystem
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 
 
  ierr = EPSSolve(eps);CHKERRQ(ierr);
  ierr = EPSSolve(eps);CHKERRQ(ierr);
  ierr = EPSGetIterationNumber(eps, &its);CHKERRQ(ierr);
  ierr = EPSGetIterationNumber(eps, &its);CHKERRQ(ierr);
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);CHKERRQ(ierr);
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);CHKERRQ(ierr);
 
 
  /*
  /*
     Optional: Get some information from the solver and display it
     Optional: Get some information from the solver and display it
  */
  */
  ierr = EPSGetType(eps,&type);CHKERRQ(ierr);
  ierr = EPSGetType(eps,&type);CHKERRQ(ierr);
  ierr = PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);CHKERRQ(ierr);
  ierr = PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);CHKERRQ(ierr);
  ierr = EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
  ierr = EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);CHKERRQ(ierr);
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);CHKERRQ(ierr);
  ierr = EPSGetTolerances(eps,&tol,&maxit);CHKERRQ(ierr);
  ierr = EPSGetTolerances(eps,&tol,&maxit);CHKERRQ(ierr);
  ierr = PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);CHKERRQ(ierr);
  ierr = PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);CHKERRQ(ierr);
 
 
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
                    Display solution and clean up
                    Display solution and clean up
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 
 
  /*
  /*
     Get number of converged approximate eigenpairs
     Get number of converged approximate eigenpairs
  */
  */
  ierr = EPSGetConverged(eps,&nconv);CHKERRQ(ierr);
  ierr = EPSGetConverged(eps,&nconv);CHKERRQ(ierr);
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
         CHKERRQ(ierr);
         CHKERRQ(ierr);
 
 
  if (nconv>0) {
  if (nconv>0) {
    /*
    /*
       Display eigenvalues and relative errors
       Display eigenvalues and relative errors
    */
    */
    ierr = PetscPrintf(PETSC_COMM_WORLD,
    ierr = PetscPrintf(PETSC_COMM_WORLD,
         "           k          ||Ax-kx||/||kx||\n"
         "           k          ||Ax-kx||/||kx||\n"
         "   ----------------- ------------------\n" );CHKERRQ(ierr);
         "   ----------------- ------------------\n" );CHKERRQ(ierr);
 
 
    for( i=0; i<nconv; i++ ) {
    for( i=0; i<nconv; i++ ) {
      /*
      /*
        Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
        Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
        ki (imaginary part)
        ki (imaginary part)
      */
      */
      ierr = EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
      ierr = EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
      /*
      /*
         Compute the relative error associated to each eigenpair
         Compute the relative error associated to each eigenpair
      */
      */
      ierr = EPSComputeRelativeError(eps,i,&error);CHKERRQ(ierr);
      ierr = EPSComputeRelativeError(eps,i,&error);CHKERRQ(ierr);
 
 
#ifdef PETSC_USE_COMPLEX
#ifdef PETSC_USE_COMPLEX
      re = PetscRealPart(kr);
      re = PetscRealPart(kr);
      im = PetscImaginaryPart(kr);
      im = PetscImaginaryPart(kr);
#else
#else
      re = kr;
      re = kr;
      im = ki;
      im = ki;
#endif
#endif
      if (im!=0.0) {
      if (im!=0.0) {
        ierr = PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g\n",re,im,error);CHKERRQ(ierr);
        ierr = PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g\n",re,im,error);CHKERRQ(ierr);
      } else {
      } else {
        ierr = PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g\n",re,error);CHKERRQ(ierr);
        ierr = PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g\n",re,error);CHKERRQ(ierr);
      }
      }
    }
    }
    ierr = PetscPrintf(PETSC_COMM_WORLD,"\n" );CHKERRQ(ierr);
    ierr = PetscPrintf(PETSC_COMM_WORLD,"\n" );CHKERRQ(ierr);
  }
  }
 
 
  /*
  /*
     Free work space
     Free work space
  */
  */
  ierr = EPSDestroy(eps);CHKERRQ(ierr);
  ierr = EPSDestroy(eps);CHKERRQ(ierr);
  ierr = MatDestroy(A);CHKERRQ(ierr);
  ierr = MatDestroy(A);CHKERRQ(ierr);
  ierr = SlepcFinalize();CHKERRQ(ierr);
  ierr = SlepcFinalize();CHKERRQ(ierr);
  return 0;
  return 0;
}
}