Actual source code: test1.c

slepc-3.17.1 2022-04-11
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:   SlepcInitialize(&argc,&argv,(char*)0,help);

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

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

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

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

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

 89:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 90:                 Create the eigensolver and set various options
 91:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

138:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139:                       Solve the eigensystem
140:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

146:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147:                     Display solution and clean up
148:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

159: /*TEST

161:    testset:
162:       args: -m 11
163:       output_file: output/test1_1.out
164:       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/"
165:       requires: !single
166:       test:
167:          suffix: 1
168:          args: -type {{toar qarnoldi linear}}
169:       test:
170:          suffix: 1_linear_gd
171:          args: -type linear -epstype gd
172:          requires: !__float128

174: TEST*/