Actual source code: test2.c

slepc-3.18.0 2022-10-01
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: */
 10: /*
 11:    Example based on spring problem in NLEVP collection [1]. See the parameters
 12:    meaning at Example 2 in [2].

 14:    [1] T. Betcke, N. J. Higham, V. Mehrmann, C. Schroder, and F. Tisseur,
 15:        NLEVP: A Collection of Nonlinear Eigenvalue Problems, MIMS EPrint
 16:        2010.98, November 2010.
 17:    [2] F. Tisseur, Backward error and condition of polynomial eigenvalue
 18:        problems, Linear Algebra and its Applications, 309 (2000), pp. 339--361,
 19:        April 2000.
 20: */

 22: static char help[] = "Test the solution of a PEP from a finite element model of "
 23:   "damped mass-spring system (problem from NLEVP collection).\n\n"
 24:   "The command line options are:\n"
 25:   "  -n <n> ... number of grid subdivisions.\n"
 26:   "  -mu <value> ... mass (default 1).\n"
 27:   "  -tau <value> ... damping constant of the dampers (default 10).\n"
 28:   "  -kappa <value> ... damping constant of the springs (default 5).\n"
 29:   "  -initv ... set an initial vector.\n\n";

 31: #include <slepcpep.h>

 33: /*
 34:    Check if computed eigenvectors have unit norm
 35: */
 36: PetscErrorCode CheckNormalizedVectors(PEP pep)
 37: {
 38:   PetscInt       i,nconv;
 39:   Mat            A;
 40:   Vec            xr,xi;
 41:   PetscReal      error=0.0,normr;
 42: #if !defined(PETSC_USE_COMPLEX)
 43:   PetscReal      normi;
 44: #endif

 47:   PEPGetConverged(pep,&nconv);
 48:   if (nconv>0) {
 49:     PEPGetOperators(pep,0,&A);
 50:     MatCreateVecs(A,&xr,&xi);
 51:     for (i=0;i<nconv;i++) {
 52:       PEPGetEigenpair(pep,i,NULL,NULL,xr,xi);
 53: #if defined(PETSC_USE_COMPLEX)
 54:       VecNorm(xr,NORM_2,&normr);
 55:       error = PetscMax(error,PetscAbsReal(normr-PetscRealConstant(1.0)));
 56: #else
 57:       VecNormBegin(xr,NORM_2,&normr);
 58:       VecNormBegin(xi,NORM_2,&normi);
 59:       VecNormEnd(xr,NORM_2,&normr);
 60:       VecNormEnd(xi,NORM_2,&normi);
 61:       error = PetscMax(error,PetscAbsReal(SlepcAbsEigenvalue(normr,normi)-PetscRealConstant(1.0)));
 62: #endif
 63:     }
 64:     VecDestroy(&xr);
 65:     VecDestroy(&xi);
 66:     if (error>100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Vectors are not normalized. Error=%g\n",(double)error);
 67:   }
 68:   return 0;
 69: }

 71: int main(int argc,char **argv)
 72: {
 73:   Mat            M,C,K,A[3];      /* problem matrices */
 74:   PEP            pep;             /* polynomial eigenproblem solver context */
 75:   PetscInt       n=30,Istart,Iend,i,nev;
 76:   PetscReal      mu=1.0,tau=10.0,kappa=5.0;
 77:   PetscBool      initv=PETSC_FALSE,skipnorm=PETSC_FALSE;
 78:   Vec            IV[2];

 81:   SlepcInitialize(&argc,&argv,(char*)0,help);

 83:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 84:   PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL);
 85:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 86:   PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL);
 87:   PetscOptionsGetBool(NULL,NULL,"-initv",&initv,NULL);
 88:   PetscOptionsGetBool(NULL,NULL,"-skipnorm",&skipnorm,NULL);

 90:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 91:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 92:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 94:   /* K is a tridiagonal */
 95:   MatCreate(PETSC_COMM_WORLD,&K);
 96:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
 97:   MatSetFromOptions(K);
 98:   MatSetUp(K);

100:   MatGetOwnershipRange(K,&Istart,&Iend);
101:   for (i=Istart;i<Iend;i++) {
102:     if (i>0) MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
103:     MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
104:     if (i<n-1) MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
105:   }

107:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
108:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

110:   /* C is a tridiagonal */
111:   MatCreate(PETSC_COMM_WORLD,&C);
112:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
113:   MatSetFromOptions(C);
114:   MatSetUp(C);

116:   MatGetOwnershipRange(C,&Istart,&Iend);
117:   for (i=Istart;i<Iend;i++) {
118:     if (i>0) MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
119:     MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
120:     if (i<n-1) MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
121:   }

123:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
124:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

126:   /* M is a diagonal matrix */
127:   MatCreate(PETSC_COMM_WORLD,&M);
128:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
129:   MatSetFromOptions(M);
130:   MatSetUp(M);
131:   MatGetOwnershipRange(M,&Istart,&Iend);
132:   for (i=Istart;i<Iend;i++) MatSetValue(M,i,i,mu,INSERT_VALUES);
133:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
134:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

136:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137:                 Create the eigensolver and set various options
138:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

140:   PEPCreate(PETSC_COMM_WORLD,&pep);
141:   A[0] = K; A[1] = C; A[2] = M;
142:   PEPSetOperators(pep,3,A);
143:   PEPSetProblemType(pep,PEP_GENERAL);
144:   PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT);
145:   if (initv) { /* initial vector */
146:     MatCreateVecs(K,&IV[0],NULL);
147:     VecSetValue(IV[0],0,-1.0,INSERT_VALUES);
148:     VecSetValue(IV[0],1,0.5,INSERT_VALUES);
149:     VecAssemblyBegin(IV[0]);
150:     VecAssemblyEnd(IV[0]);
151:     MatCreateVecs(K,&IV[1],NULL);
152:     VecSetValue(IV[1],0,4.0,INSERT_VALUES);
153:     VecSetValue(IV[1],2,1.5,INSERT_VALUES);
154:     VecAssemblyBegin(IV[1]);
155:     VecAssemblyEnd(IV[1]);
156:     PEPSetInitialSpace(pep,2,IV);
157:     VecDestroy(&IV[0]);
158:     VecDestroy(&IV[1]);
159:   }
160:   PEPSetFromOptions(pep);

162:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163:                       Solve the eigensystem
164:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

166:   PEPSolve(pep);
167:   PEPGetDimensions(pep,&nev,NULL,NULL);
168:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);

170:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171:                     Display solution and clean up
172:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

174:   PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
175:   if (!skipnorm) CheckNormalizedVectors(pep);
176:   PEPDestroy(&pep);
177:   MatDestroy(&M);
178:   MatDestroy(&C);
179:   MatDestroy(&K);
180:   SlepcFinalize();
181:   return 0;
182: }

184: /*TEST

186:    testset:
187:       args: -pep_nev 4 -initv
188:       requires: !single
189:       output_file: output/test2_1.out
190:       test:
191:          suffix: 1
192:          args: -pep_type {{toar linear}}
193:       test:
194:          suffix: 1_toar_mgs
195:          args: -pep_type toar -bv_orthog_type mgs
196:       test:
197:          suffix: 1_qarnoldi
198:          args: -pep_type qarnoldi -bv_orthog_refine never
199:       test:
200:          suffix: 1_linear_gd
201:          args: -pep_type linear -pep_linear_eps_type gd -pep_linear_explicitmatrix

203:    testset:
204:       args: -pep_target -0.43 -pep_nev 4 -pep_ncv 20 -st_type sinvert
205:       output_file: output/test2_2.out
206:       test:
207:          suffix: 2
208:          args: -pep_type {{toar linear}}
209:       test:
210:          suffix: 2_toar_scaleboth
211:          args: -pep_type toar -pep_scale both
212:       test:
213:          suffix: 2_toar_transform
214:          args: -pep_type toar -st_transform
215:       test:
216:          suffix: 2_qarnoldi
217:          args: -pep_type qarnoldi -bv_orthog_refine always
218:       test:
219:          suffix: 2_linear_explicit
220:          args: -pep_type linear -pep_linear_explicitmatrix -pep_linear_linearization 0,1
221:       test:
222:          suffix: 2_linear_explicit_her
223:          args: -pep_type linear -pep_linear_explicitmatrix -pep_hermitian -pep_linear_linearization 0,1
224:       test:
225:          suffix: 2_stoar
226:          args: -pep_type stoar -pep_hermitian
227:       test:
228:          suffix: 2_jd
229:          args: -pep_type jd -st_type precond -pep_max_it 200 -pep_ncv 24
230:          requires: !single

232:    test:
233:       suffix: 3
234:       args: -pep_nev 12 -pep_extract {{none norm residual structured}} -pep_monitor_cancel
235:       requires: !single

237:    testset:
238:       args: -st_type sinvert -pep_target -0.43 -pep_nev 4
239:       output_file: output/test2_2.out
240:       test:
241:          suffix: 4_schur
242:          args: -pep_refine simple -pep_refine_scheme schur
243:       test:
244:          suffix: 4_mbe
245:          args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
246:       test:
247:          suffix: 4_explicit
248:          args: -pep_refine simple -pep_refine_scheme explicit
249:       test:
250:          suffix: 4_multiple_schur
251:          args: -pep_refine multiple -pep_refine_scheme schur
252:          requires: !single
253:       test:
254:          suffix: 4_multiple_mbe
255:          args: -pep_refine multiple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu -pep_refine_pc_factor_shift_type nonzero
256:       test:
257:          suffix: 4_multiple_explicit
258:          args: -pep_refine multiple -pep_refine_scheme explicit
259:          requires: !single

261:    test:
262:       suffix: 5
263:       nsize: 2
264:       args: -pep_type linear -pep_linear_explicitmatrix -pep_general -pep_target -0.43 -pep_nev 4 -pep_ncv 20 -st_type sinvert -pep_linear_st_ksp_type bcgs -pep_linear_st_pc_type bjacobi
265:       output_file: output/test2_2.out

267:    test:
268:       suffix: 6
269:       args: -pep_type linear -pep_nev 12 -pep_extract {{none norm residual}}
270:       requires: !single
271:       output_file: output/test2_3.out

273:    test:
274:       suffix: 7
275:       args: -pep_nev 12 -pep_extract {{none norm residual structured}} -pep_refine multiple
276:       requires: !single
277:       output_file: output/test2_3.out

279:    testset:
280:       args: -st_type sinvert -pep_target -0.43 -pep_nev 4 -st_transform
281:       output_file: output/test2_2.out
282:       test:
283:          suffix: 8_schur
284:          args: -pep_refine simple -pep_refine_scheme schur
285:       test:
286:          suffix: 8_mbe
287:          args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
288:       test:
289:          suffix: 8_explicit
290:          args: -pep_refine simple -pep_refine_scheme explicit
291:       test:
292:          suffix: 8_multiple_schur
293:          args: -pep_refine multiple -pep_refine_scheme schur
294:       test:
295:          suffix: 8_multiple_mbe
296:          args: -pep_refine multiple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
297:       test:
298:          suffix: 8_multiple_explicit
299:          args: -pep_refine multiple -pep_refine_scheme explicit

301:    testset:
302:       nsize: 2
303:       args: -st_type sinvert -pep_target -0.49 -pep_nev 4 -pep_refine_partitions 2 -st_ksp_type bcgs -st_pc_type bjacobi -pep_scale diagonal -pep_scale_its 4
304:       output_file: output/test2_2.out
305:       test:
306:          suffix: 9_mbe
307:          args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
308:       test:
309:          suffix: 9_explicit
310:          args: -pep_refine simple -pep_refine_scheme explicit
311:       test:
312:          suffix: 9_multiple_mbe
313:          args: -pep_refine multiple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
314:          requires: !single
315:       test:
316:          suffix: 9_multiple_explicit
317:          args: -pep_refine multiple -pep_refine_scheme explicit
318:          requires: !single

320:    test:
321:       suffix: 10
322:       nsize: 4
323:       args: -st_type sinvert -pep_target -0.43 -pep_nev 4 -pep_refine simple -pep_refine_scheme explicit -pep_refine_partitions 2 -st_ksp_type bcgs -st_pc_type bjacobi -pep_scale diagonal -pep_scale_its 4
324:       output_file: output/test2_2.out

326:    testset:
327:       args: -pep_nev 4 -initv -mat_type aijcusparse
328:       output_file: output/test2_1.out
329:       requires: cuda !single
330:       test:
331:          suffix: 11_cuda
332:          args: -pep_type {{toar linear}}
333:       test:
334:          suffix: 11_cuda_qarnoldi
335:          args: -pep_type qarnoldi -bv_orthog_refine never
336:       test:
337:          suffix: 11_cuda_linear_gd
338:          args: -pep_type linear -pep_linear_eps_type gd -pep_linear_explicitmatrix

340:    test:
341:       suffix: 12
342:       nsize: 2
343:       args: -pep_type jd -ds_parallel synchronized -pep_target -0.43 -pep_nev 4 -pep_ncv 20
344:       requires: !single

346:    test:
347:       suffix: 13
348:       args: -pep_nev 12 -pep_view_values draw -pep_monitor draw::draw_lg
349:       requires: x !single
350:       output_file: output/test2_3.out

352:    test:
353:       suffix: 14
354:       requires: complex double
355:       args: -pep_type ciss -rg_type ellipse -rg_ellipse_center -48.5 -rg_ellipse_radius 1.5 -pep_ciss_delta 1e-10

357: TEST*/