Actual source code: test1.c

slepc-3.21.1 2024-04-26
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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[] = "Test the solution of a PEP without calling PEPSetFromOptions (based on ex16.c).\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 14:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n"
 15:   "  -type <pep_type> = pep type to test.\n"
 16:   "  -epstype <eps_type> = eps type to test (for linear).\n\n";

 18: #include <slepcpep.h>

 20: int main(int argc,char **argv)
 21: {
 22:   Mat            M,C,K,A[3];      /* problem matrices */
 23:   PEP            pep;             /* polynomial eigenproblem solver context */
 24:   PetscInt       N,n=10,m,Istart,Iend,II,nev,i,j;
 25:   PetscReal      keep;
 26:   PetscBool      flag,isgd2,epsgiven,lock;
 27:   char           peptype[30] = "linear",epstype[30] = "";
 28:   EPS            eps;
 29:   ST             st;
 30:   KSP            ksp;
 31:   PC             pc;

 33:   PetscFunctionBeginUser;
 34:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

 36:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 37:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
 38:   if (!flag) m=n;
 39:   N = n*m;
 40:   PetscCall(PetscOptionsGetString(NULL,NULL,"-type",peptype,sizeof(peptype),NULL));
 41:   PetscCall(PetscOptionsGetString(NULL,NULL,"-epstype",epstype,sizeof(epstype),&epsgiven));
 42:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));

 44:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 45:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 46:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 48:   /* K is the 2-D Laplacian */
 49:   PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
 50:   PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N));
 51:   PetscCall(MatSetFromOptions(K));
 52:   PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
 53:   for (II=Istart;II<Iend;II++) {
 54:     i = II/n; j = II-i*n;
 55:     if (i>0) PetscCall(MatSetValue(K,II,II-n,-1.0,INSERT_VALUES));
 56:     if (i<m-1) PetscCall(MatSetValue(K,II,II+n,-1.0,INSERT_VALUES));
 57:     if (j>0) PetscCall(MatSetValue(K,II,II-1,-1.0,INSERT_VALUES));
 58:     if (j<n-1) PetscCall(MatSetValue(K,II,II+1,-1.0,INSERT_VALUES));
 59:     PetscCall(MatSetValue(K,II,II,4.0,INSERT_VALUES));
 60:   }
 61:   PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
 62:   PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));

 64:   /* C is the 1-D Laplacian on horizontal lines */
 65:   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
 66:   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N));
 67:   PetscCall(MatSetFromOptions(C));
 68:   PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
 69:   for (II=Istart;II<Iend;II++) {
 70:     i = II/n; j = II-i*n;
 71:     if (j>0) PetscCall(MatSetValue(C,II,II-1,-1.0,INSERT_VALUES));
 72:     if (j<n-1) PetscCall(MatSetValue(C,II,II+1,-1.0,INSERT_VALUES));
 73:     PetscCall(MatSetValue(C,II,II,2.0,INSERT_VALUES));
 74:   }
 75:   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
 76:   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));

 78:   /* M is a diagonal matrix */
 79:   PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
 80:   PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N));
 81:   PetscCall(MatSetFromOptions(M));
 82:   PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
 83:   for (II=Istart;II<Iend;II++) PetscCall(MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES));
 84:   PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
 85:   PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));

 87:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 88:                 Create the eigensolver and set various options
 89:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 91:   PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
 92:   A[0] = K; A[1] = C; A[2] = M;
 93:   PetscCall(PEPSetOperators(pep,3,A));
 94:   PetscCall(PEPSetProblemType(pep,PEP_GENERAL));
 95:   PetscCall(PEPSetDimensions(pep,4,20,PETSC_DEFAULT));
 96:   PetscCall(PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT));

 98:   /*
 99:      Set solver type at runtime
100:   */
101:   PetscCall(PEPSetType(pep,peptype));
102:   if (epsgiven) {
103:     PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&flag));
104:     if (flag) {
105:       PetscCall(PEPLinearGetEPS(pep,&eps));
106:       PetscCall(PetscStrcmp(epstype,"gd2",&isgd2));
107:       if (isgd2) {
108:         PetscCall(EPSSetType(eps,EPSGD));
109:         PetscCall(EPSGDSetDoubleExpansion(eps,PETSC_TRUE));
110:       } else PetscCall(EPSSetType(eps,epstype));
111:       PetscCall(EPSGetST(eps,&st));
112:       PetscCall(STGetKSP(st,&ksp));
113:       PetscCall(KSPGetPC(ksp,&pc));
114:       PetscCall(PCSetType(pc,PCJACOBI));
115:       PetscCall(PetscObjectTypeCompare((PetscObject)eps,EPSGD,&flag));
116:     }
117:     PetscCall(PEPLinearSetExplicitMatrix(pep,PETSC_TRUE));
118:   }
119:   PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPQARNOLDI,&flag));
120:   if (flag) {
121:     PetscCall(STCreate(PETSC_COMM_WORLD,&st));
122:     PetscCall(STSetTransform(st,PETSC_TRUE));
123:     PetscCall(PEPSetST(pep,st));
124:     PetscCall(STDestroy(&st));
125:     PetscCall(PEPQArnoldiGetRestart(pep,&keep));
126:     PetscCall(PEPQArnoldiGetLocking(pep,&lock));
127:     if (!lock && keep<0.6) PetscCall(PEPQArnoldiSetRestart(pep,0.6));
128:   }
129:   PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPTOAR,&flag));
130:   if (flag) {
131:     PetscCall(PEPTOARGetRestart(pep,&keep));
132:     PetscCall(PEPTOARGetLocking(pep,&lock));
133:     if (!lock && keep<0.6) PetscCall(PEPTOARSetRestart(pep,0.6));
134:   }

136:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137:                       Solve the eigensystem
138:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

140:   PetscCall(PEPSolve(pep));
141:   PetscCall(PEPGetDimensions(pep,&nev,NULL,NULL));
142:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));

144:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145:                     Display solution and clean up
146:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

148:   PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
149:   PetscCall(PEPDestroy(&pep));
150:   PetscCall(MatDestroy(&M));
151:   PetscCall(MatDestroy(&C));
152:   PetscCall(MatDestroy(&K));
153:   PetscCall(SlepcFinalize());
154:   return 0;
155: }

157: /*TEST

159:    testset:
160:       args: -m 11
161:       output_file: output/test1_1.out
162:       filter: sed -e "s/1.16403/1.16404/g" | sed -e "s/1.65362i/1.65363i/g" | sed -e "s/-1.16404-1.65363i, -1.16404+1.65363i/-1.16404+1.65363i, -1.16404-1.65363i/" | sed -e "s/-0.51784-1.31039i, -0.51784+1.31039i/-0.51784+1.31039i, -0.51784-1.31039i/"
163:       requires: !single
164:       test:
165:          suffix: 1
166:          args: -type {{toar qarnoldi linear}}
167:       test:
168:          suffix: 1_linear_gd
169:          args: -type linear -epstype gd
170:          requires: !__float128

172: TEST*/