Actual source code: ex16.c

slepc-3.11.2 2019-07-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[] = "Simple quadratic eigenvalue problem.\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\n";

 16: #include <slepcpep.h>

 18: int main(int argc,char **argv)
 19: {
 20:   Mat            M,C,K,A[3];      /* problem matrices */
 21:   PEP            pep;             /* polynomial eigenproblem solver context */
 22:   PetscInt       N,n=10,m,Istart,Iend,II,nev,i,j,nconv;
 23:   PetscBool      flag,terse;
 24:   PetscReal      error,re,im;
 25:   PetscScalar    kr,ki;
 26:   Vec            xr,xi;

 29:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 31:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 32:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 33:   if (!flag) m=n;
 34:   N = n*m;
 35:   PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%D (%Dx%D grid)\n\n",N,n,m);

 37:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 38:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 39:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

 73:   /* M is a diagonal matrix */
 74:   MatCreate(PETSC_COMM_WORLD,&M);
 75:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
 76:   MatSetFromOptions(M);
 77:   MatSetUp(M);
 78:   MatGetOwnershipRange(M,&Istart,&Iend);
 79:   for (II=Istart;II<Iend;II++) {
 80:     MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES);
 81:   }
 82:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
 83:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

 85:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 86:                 Create the eigensolver and set various options
 87:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 89:   /*
 90:      Create eigensolver context
 91:   */
 92:   PEPCreate(PETSC_COMM_WORLD,&pep);

 94:   /*
 95:      Set matrices and problem type
 96:   */
 97:   A[0] = K; A[1] = C; A[2] = M;
 98:   PEPSetOperators(pep,3,A);
 99:   PEPSetProblemType(pep,PEP_HERMITIAN);

101:   /*
102:      Set solver parameters at runtime
103:   */
104:   PEPSetFromOptions(pep);

106:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107:                       Solve the eigensystem
108:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

110:   PEPSolve(pep);

112:   /*
113:      Optional: Get some information from the solver and display it
114:   */
115:   PEPGetDimensions(pep,&nev,NULL,NULL);
116:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

118:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119:                     Display solution and clean up
120:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

122:   /* show detailed info unless -terse option is given by user */
123:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
124:   if (terse) {
125:     PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
126:   } else {
127:     PEPGetConverged(pep,&nconv);
128:     if (nconv>0) {
129:       MatCreateVecs(M,&xr,&xi);
130:       /* display eigenvalues and relative errors */
131:       PetscPrintf(PETSC_COMM_WORLD,
132:            "\n           k          ||P(k)x||/||kx||\n"
133:            "   ----------------- ------------------\n");
134:       for (i=0;i<nconv;i++) {
135:         /* get converged eigenpairs */
136:         PEPGetEigenpair(pep,i,&kr,&ki,xr,xi);
137:         /* compute the relative error associated to each eigenpair */
138:         PEPComputeError(pep,i,PEP_ERROR_BACKWARD,&error);
139: #if defined(PETSC_USE_COMPLEX)
140:         re = PetscRealPart(kr);
141:         im = PetscImaginaryPart(kr);
142: #else
143:         re = kr;
144:         im = ki;
145: #endif
146:         if (im!=0.0) {
147:           PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi   %12g\n",(double)re,(double)im,(double)error);
148:         } else {
149:           PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g\n",(double)re,(double)error);
150:         }
151:       }
152:       PetscPrintf(PETSC_COMM_WORLD,"\n");
153:       VecDestroy(&xr);
154:       VecDestroy(&xi);
155:     }
156:   }
157:   PEPDestroy(&pep);
158:   MatDestroy(&M);
159:   MatDestroy(&C);
160:   MatDestroy(&K);
161:   SlepcFinalize();
162:   return ierr;
163: }

165: /*TEST

167:    testset:
168:       args: -pep_nev 4 -pep_ncv 21 -n 12 -terse
169:       requires: !complex
170:       output_file: output/ex16_1.out
171:       test:
172:          suffix: 1
173:          args: -pep_type {{toar qarnoldi}}
174:       test:
175:          suffix: 1_linear
176:          args: -pep_type linear -pep_linear_explicitmatrix
177:       test:
178:          suffix: 1_linear_symm
179:          args: -pep_type linear -pep_linear_explicitmatrix -pep_linear_eps_gen_indefinite
180:          requires: !single
181:       test:
182:          suffix: 1_stoar
183:          args: -pep_type stoar
184:          requires: !single
185:       test:
186:          suffix: 1_stoar_t
187:          args: -pep_type stoar -st_transform
188:          requires: !single

190: TEST*/