Actual source code: test14.c

slepc-3.7.2 2016-07-19
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  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.

 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.

 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/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: static char help[] = "Test EPS interface functions.\n\n";

 24: #include <slepceps.h>

 28: int main(int argc,char **argv)
 29: {
 30:   Mat                A,B;         /* problem matrix */
 31:   EPS                eps;         /* eigenproblem solver context */
 32:   ST                 st;
 33:   KSP                ksp;
 34:   DS                 ds;
 35:   PetscReal          cut,tol;
 36:   PetscScalar        target;
 37:   PetscInt           n=20,i,its,nev,ncv,mpd,Istart,Iend;
 38:   PetscBool          flg;
 39:   EPSConvergedReason reason;
 40:   EPSType            type;
 41:   EPSExtraction      extr;
 42:   EPSBalance         bal;
 43:   EPSWhich           which;
 44:   EPSConv            conv;
 45:   EPSProblemType     ptype;
 46:   PetscErrorCode     ierr;
 47:   PetscViewerAndFormat *vf;

 49:   SlepcInitialize(&argc,&argv,(char*)0,help);
 50:   PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Eigenproblem, n=%D\n\n",n);

 52:   MatCreate(PETSC_COMM_WORLD,&A);
 53:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 54:   MatSetFromOptions(A);
 55:   MatSetUp(A);
 56:   MatGetOwnershipRange(A,&Istart,&Iend);
 57:   for (i=Istart;i<Iend;i++) {
 58:     MatSetValue(A,i,i,i+1,INSERT_VALUES);
 59:   }
 60:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 61:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 63:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 64:              Create eigensolver and test interface functions
 65:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 66:   EPSCreate(PETSC_COMM_WORLD,&eps);
 67:   EPSSetOperators(eps,A,NULL);
 68:   EPSGetOperators(eps,&B,NULL);
 69:   MatView(B,NULL);

 71:   EPSSetType(eps,EPSKRYLOVSCHUR);
 72:   EPSGetType(eps,&type);
 73:   PetscPrintf(PETSC_COMM_WORLD," Type set to %s\n",type);

 75:   EPSGetProblemType(eps,&ptype);
 76:   PetscPrintf(PETSC_COMM_WORLD," Problem type before changing = %d",(int)ptype);
 77:   EPSSetProblemType(eps,EPS_HEP);
 78:   EPSGetProblemType(eps,&ptype);
 79:   PetscPrintf(PETSC_COMM_WORLD," ... changed to %d.",(int)ptype);
 80:   EPSIsGeneralized(eps,&flg);
 81:   if (flg) { PetscPrintf(PETSC_COMM_WORLD," generalized"); }
 82:   EPSIsHermitian(eps,&flg);
 83:   if (flg) { PetscPrintf(PETSC_COMM_WORLD," hermitian"); }
 84:   EPSIsPositive(eps,&flg);
 85:   if (flg) { PetscPrintf(PETSC_COMM_WORLD," positive"); }

 87:   EPSGetExtraction(eps,&extr);
 88:   PetscPrintf(PETSC_COMM_WORLD,"\n Extraction before changing = %d",(int)extr);
 89:   EPSSetExtraction(eps,EPS_HARMONIC);
 90:   EPSGetExtraction(eps,&extr);
 91:   PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)extr);

 93:   EPSSetBalance(eps,EPS_BALANCE_ONESIDE,8,1e-6);
 94:   EPSGetBalance(eps,&bal,&its,&cut);
 95:   PetscPrintf(PETSC_COMM_WORLD," Balance: %s, its=%D, cutoff=%g\n",EPSBalanceTypes[bal],its,(double)cut);

 97:   EPSSetTarget(eps,4.8);
 98:   EPSGetTarget(eps,&target);
 99:   EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
100:   EPSGetWhichEigenpairs(eps,&which);
101:   PetscPrintf(PETSC_COMM_WORLD," Which = %d, target = %g\n",(int)which,(double)PetscRealPart(target));

103:   EPSSetDimensions(eps,4,PETSC_DEFAULT,PETSC_DEFAULT);
104:   EPSGetDimensions(eps,&nev,&ncv,&mpd);
105:   PetscPrintf(PETSC_COMM_WORLD," Dimensions: nev=%D, ncv=%D, mpd=%D\n",nev,ncv,mpd);

107:   EPSSetTolerances(eps,2.2e-4,200);
108:   EPSGetTolerances(eps,&tol,&its);
109:   PetscPrintf(PETSC_COMM_WORLD," Tolerance = %.5f, max_its = %D\n",(double)tol,its);

111:   EPSSetConvergenceTest(eps,EPS_CONV_ABS);
112:   EPSGetConvergenceTest(eps,&conv);
113:   PetscPrintf(PETSC_COMM_WORLD," Convergence test = %d\n",(int)conv);

115:   PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);
116:   EPSMonitorSet(eps,(PetscErrorCode (*)(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))EPSMonitorFirst,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
117:   EPSMonitorCancel(eps);

119:   EPSGetST(eps,&st);
120:   STGetKSP(st,&ksp);
121:   KSPSetTolerances(ksp,1e-8,1e-50,PETSC_DEFAULT,PETSC_DEFAULT);
122:   STView(st,NULL);
123:   EPSGetDS(eps,&ds);
124:   DSView(ds,NULL);

126:   EPSSetFromOptions(eps);
127:   EPSSolve(eps);
128:   EPSGetConvergedReason(eps,&reason);
129:   EPSGetIterationNumber(eps,&its);
130:   PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d, its=%D\n",(int)reason,its);

132:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133:                     Display solution and clean up
134:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
136:   EPSDestroy(&eps);
137:   MatDestroy(&A);
138:   SlepcFinalize();
139:   return ierr;
140: }