Actual source code: test32.c

slepc-3.12.0 2019-09-30
Report Typos and Errors
  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[] = "Tests a GHEP problem with symmetric matrices.\n\n";

 13: #include <slepceps.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat            A,B;        /* matrices */
 18:   EPS            eps;        /* eigenproblem solver context */
 19:   ST             st;
 20:   KSP            ksp;
 21:   PC             pc;
 22:   PCType         pctype;
 23:   PetscInt       N,n=45,m,Istart,Iend,II,i,j;
 24:   PetscBool      flag;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 30:   if (!flag) m=n;
 31:   N = n*m;
 32:   PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized Symmetric Eigenproblem, N=%D (%Dx%D grid)\n\n",N,n,m);

 34:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 35:      Compute the matrices that define the eigensystem, Ax=kBx
 36:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 38:   MatCreate(PETSC_COMM_WORLD,&A);
 39:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 40:   MatSetFromOptions(A);
 41:   MatSetUp(A);

 43:   MatCreate(PETSC_COMM_WORLD,&B);
 44:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
 45:   MatSetFromOptions(B);
 46:   MatSetUp(B);

 48:   MatGetOwnershipRange(A,&Istart,&Iend);
 49:   for (II=Istart;II<Iend;II++) {
 50:     i = II/n; j = II-i*n;
 51:     if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
 52:     if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
 53:     if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
 54:     if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
 55:     MatSetValue(A,II,II,4.0,INSERT_VALUES);
 56:     MatSetValue(B,II,II,2.0/PetscLogScalar(II+2),INSERT_VALUES);
 57:   }
 58:   MatSetValue(B,0,1,0.4,INSERT_VALUES);
 59:   MatSetValue(B,1,0,0.4,INSERT_VALUES);

 61:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 62:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 63:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 64:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 66:   MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);
 67:   MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE);

 69:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 70:                 Create the eigensolver and solve the problem
 71:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 73:   EPSCreate(PETSC_COMM_WORLD,&eps);
 74:   EPSSetOperators(eps,A,B);
 75:   EPSSetProblemType(eps,EPS_GHEP);
 76:   EPSSetFromOptions(eps);
 77:   EPSSetUp(eps);
 78:   EPSGetST(eps,&st);
 79:   STGetKSP(st,&ksp);
 80:   KSPGetPC(ksp,&pc);
 81:   PCGetType(pc,&pctype);
 82:   PetscPrintf(PETSC_COMM_WORLD," Using %s for the PC\n",pctype);
 83:   EPSSolve(eps);
 84:   EPSErrorView(eps,EPS_ERROR_BACKWARD,NULL);

 86:   EPSDestroy(&eps);
 87:   MatDestroy(&A);
 88:   MatDestroy(&B);
 89:   SlepcFinalize();
 90:   return ierr;
 91: }

 93: /*TEST

 95:    test:
 96:       suffix: 1
 97:       args: -n 18 -eps_nev 3 -st_type sinvert -eps_target 1.02
 98:       requires: !single

100:    test:
101:       suffix: 2
102:       args: -n 18 -eps_type ciss -rg_interval_endpoints 1.0,1.2

104: TEST*/