Actual source code: ex16.c

slepc-3.20.1 2023-11-27
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[] = "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;

 30:   PetscFunctionBeginUser;
 31:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

 33:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 34:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
 35:   if (!flag) m=n;
 36:   N = n*m;
 37:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " 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:   PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
 45:   PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N));
 46:   PetscCall(MatSetFromOptions(K));
 47:   PetscCall(MatSetUp(K));
 48:   PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
 49:   for (II=Istart;II<Iend;II++) {
 50:     i = II/n; j = II-i*n;
 51:     if (i>0) PetscCall(MatSetValue(K,II,II-n,-1.0,INSERT_VALUES));
 52:     if (i<m-1) PetscCall(MatSetValue(K,II,II+n,-1.0,INSERT_VALUES));
 53:     if (j>0) PetscCall(MatSetValue(K,II,II-1,-1.0,INSERT_VALUES));
 54:     if (j<n-1) PetscCall(MatSetValue(K,II,II+1,-1.0,INSERT_VALUES));
 55:     PetscCall(MatSetValue(K,II,II,4.0,INSERT_VALUES));
 56:   }
 57:   PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
 58:   PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));

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

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

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

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

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

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

110:   /*
111:      Set solver parameters at runtime
112:   */
113:   PetscCall(PEPSetFromOptions(pep));

115:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116:                       Solve the eigensystem
117:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

119:   PetscCall(PEPSolve(pep));

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

127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128:                     Display solution and clean up
129:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

170: /*TEST

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

195: TEST*/