Actual source code: ex16.c

slepc-3.16.1 2021-11-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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;
 27:   BV             V;
 28:   PetscRandom    rand;

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

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

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

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

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

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

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

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

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

103:   /*
104:      In complex scalars, use a real initial vector since in this example
105:      the matrices are all real, then all vectors generated by the solver
106:      will have a zero imaginary part. This is not really necessary.
107:   */
108:   PEPGetBV(pep,&V);
109:   BVGetRandomContext(V,&rand);
110:   PetscRandomSetInterval(rand,-1,1);

112:   /*
113:      Set solver parameters at runtime
114:   */
115:   PEPSetFromOptions(pep);

117:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118:                       Solve the eigensystem
119:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

121:   PEPSolve(pep);

123:   /*
124:      Optional: Get some information from the solver and display it
125:   */
126:   PEPGetDimensions(pep,&nev,NULL,NULL);
127:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

129:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130:                     Display solution and clean up
131:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

176: /*TEST

178:    testset:
179:       args: -pep_nev 4 -pep_ncv 21 -n 12 -terse
180:       output_file: output/ex16_1.out
181:       test:
182:          suffix: 1
183:          args: -pep_type {{toar qarnoldi}}
184:       test:
185:          suffix: 1_linear
186:          args: -pep_type linear -pep_linear_explicitmatrix
187:          requires: !single
188:       test:
189:          suffix: 1_linear_symm
190:          args: -pep_type linear -pep_linear_explicitmatrix -pep_linear_eps_gen_indefinite -pep_scale scalar -pep_linear_bv_definite_tol 1e-12
191:          requires: !single
192:       test:
193:          suffix: 1_stoar
194:          args: -pep_type stoar -pep_scale scalar
195:          requires: double !cuda
196:       test:
197:          suffix: 1_stoar_t
198:          args: -pep_type stoar -pep_scale scalar -st_transform
199:          requires: double !cuda

201: TEST*/