Actual source code: test3.c

slepc-3.7.0 2016-05-16
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 PEP interface functions.\n\n";

 24: #include <slepcpep.h>

 28: int main(int argc,char **argv)
 29: {
 30:   Mat                A[3],B;      /* problem matrices */
 31:   PEP                pep;         /* eigenproblem solver context */
 32:   ST                 st;
 33:   KSP                ksp;
 34:   DS                 ds;
 35:   PetscReal          tol,alpha;
 36:   PetscScalar        target;
 37:   PetscInt           n=20,i,its,nev,ncv,mpd,Istart,Iend,nmat;
 38:   PEPWhich           which;
 39:   PEPConvergedReason reason;
 40:   PEPType            type;
 41:   PEPExtract         extr;
 42:   PEPBasis           basis;
 43:   PEPScale           scale;
 44:   PEPRefine          refine;
 45:   PEPRefineScheme    rscheme;
 46:   PEPConv            conv;
 47:   PEPStop            stop;
 48:   PEPProblemType     ptype;
 49:   PetscErrorCode     ierr;
 50:   PetscViewerAndFormat *vf;

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

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

 66:   MatCreate(PETSC_COMM_WORLD,&A[1]);
 67:   MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,n,n);
 68:   MatSetFromOptions(A[1]);
 69:   MatSetUp(A[1]);
 70:   MatGetOwnershipRange(A[1],&Istart,&Iend);
 71:   for (i=Istart;i<Iend;i++) {
 72:     MatSetValue(A[1],i,i,1.0,INSERT_VALUES);
 73:   }
 74:   MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY);
 75:   MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY);

 77:   MatCreate(PETSC_COMM_WORLD,&A[2]);
 78:   MatSetSizes(A[2],PETSC_DECIDE,PETSC_DECIDE,n,n);
 79:   MatSetFromOptions(A[2]);
 80:   MatSetUp(A[2]);
 81:   MatGetOwnershipRange(A[1],&Istart,&Iend);
 82:   for (i=Istart;i<Iend;i++) {
 83:     MatSetValue(A[2],i,i,n/(PetscReal)(i+1),INSERT_VALUES);
 84:   }
 85:   MatAssemblyBegin(A[2],MAT_FINAL_ASSEMBLY);
 86:   MatAssemblyEnd(A[2],MAT_FINAL_ASSEMBLY);

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:              Create eigensolver and test interface functions
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 91:   PEPCreate(PETSC_COMM_WORLD,&pep);
 92:   PEPSetOperators(pep,3,A);
 93:   PEPGetNumMatrices(pep,&nmat);
 94:   PetscPrintf(PETSC_COMM_WORLD," Polynomial of degree %d\n",(int)nmat-1);
 95:   PEPGetOperators(pep,0,&B);
 96:   MatView(B,NULL);

 98:   PEPSetType(pep,PEPTOAR);
 99:   PEPGetType(pep,&type);
100:   PetscPrintf(PETSC_COMM_WORLD," Type set to %s\n",type);

102:   PEPGetProblemType(pep,&ptype);
103:   PetscPrintf(PETSC_COMM_WORLD," Problem type before changing = %d",(int)ptype);
104:   PEPSetProblemType(pep,PEP_HERMITIAN);
105:   PEPGetProblemType(pep,&ptype);
106:   PetscPrintf(PETSC_COMM_WORLD," ... changed to %d.",(int)ptype);

108:   PEPGetExtract(pep,&extr);
109:   PetscPrintf(PETSC_COMM_WORLD,"\n Extraction before changing = %d",(int)extr);
110:   PEPSetExtract(pep,PEP_EXTRACT_STRUCTURED);
111:   PEPGetExtract(pep,&extr);
112:   PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)extr);

114:   PEPSetScale(pep,PEP_SCALE_SCALAR,.1,NULL,NULL,5,1.0);
115:   PEPGetScale(pep,&scale,&alpha,NULL,NULL,&its,NULL);
116:   PetscPrintf(PETSC_COMM_WORLD," Scaling: %s, alpha=%g, its=%D\n",PEPScaleTypes[scale],(double)alpha,its);

118:   PEPSetBasis(pep,PEP_BASIS_CHEBYSHEV1);
119:   PEPGetBasis(pep,&basis);
120:   PetscPrintf(PETSC_COMM_WORLD," Polynomial basis: %s\n",PEPBasisTypes[basis]);

122:   PEPSetRefine(pep,PEP_REFINE_SIMPLE,1,1e-9,2,PEP_REFINE_SCHEME_SCHUR);
123:   PEPGetRefine(pep,&refine,NULL,&tol,&its,&rscheme);
124:   PetscPrintf(PETSC_COMM_WORLD," Refinement: %s, tol=%g, its=%D, scheme=%s\n",PEPRefineTypes[refine],(double)tol,its,PEPRefineSchemes[rscheme]);

126:   PEPSetTarget(pep,4.8);
127:   PEPGetTarget(pep,&target);
128:   PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE);
129:   PEPGetWhichEigenpairs(pep,&which);
130:   PetscPrintf(PETSC_COMM_WORLD," Which = %d, target = %g\n",(int)which,(double)PetscRealPart(target));

132:   PEPSetDimensions(pep,4,PETSC_DEFAULT,PETSC_DEFAULT);
133:   PEPGetDimensions(pep,&nev,&ncv,&mpd);
134:   PetscPrintf(PETSC_COMM_WORLD," Dimensions: nev=%D, ncv=%D, mpd=%D\n",nev,ncv,mpd);

136:   PEPSetTolerances(pep,2.2e-4,200);
137:   PEPGetTolerances(pep,&tol,&its);
138:   PetscPrintf(PETSC_COMM_WORLD," Tolerance = %.5f, max_its = %D\n",(double)tol,its);

140:   PEPSetConvergenceTest(pep,PEP_CONV_ABS);
141:   PEPGetConvergenceTest(pep,&conv);
142:   PEPSetStoppingTest(pep,PEP_STOP_BASIC);
143:   PEPGetStoppingTest(pep,&stop);
144:   PetscPrintf(PETSC_COMM_WORLD," Convergence test = %d, stopping test = %d\n",(int)conv,(int)stop);

146:   PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);
147:   PEPMonitorSet(pep,(PetscErrorCode (*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))PEPMonitorFirst,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
148:   PEPMonitorCancel(pep);

150:   PEPGetST(pep,&st);
151:   STGetKSP(st,&ksp);
152:   KSPSetTolerances(ksp,1e-8,1e-50,PETSC_DEFAULT,PETSC_DEFAULT);
153:   STView(st,NULL);
154:   PEPGetDS(pep,&ds);
155:   DSView(ds,NULL);

157:   PEPSetFromOptions(pep);
158:   PEPSolve(pep);
159:   PEPGetConvergedReason(pep,&reason);
160:   PEPGetIterationNumber(pep,&its);
161:   PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d, its=%D\n",(int)reason,its);

163:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164:                     Display solution and clean up
165:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166:   PEPErrorView(pep,PEP_ERROR_RELATIVE,NULL);
167:   PEPDestroy(&pep);
168:   MatDestroy(&A[0]);
169:   MatDestroy(&A[1]);
170:   MatDestroy(&A[2]);
171:   SlepcFinalize();
172:   return ierr;
173: }