Actual source code: test2.c

slepc-3.17.2 2022-08-09
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:   PetscFunctionReturn(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];

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

183: /*TEST

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

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

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

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

260:    test:
261:       suffix: 5
262:       nsize: 2
263:       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
264:       output_file: output/test2_2.out

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

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

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

300:    testset:
301:       nsize: 2
302:       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
303:       output_file: output/test2_2.out
304:       test:
305:          suffix: 9_mbe
306:          args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
307:       test:
308:          suffix: 9_explicit
309:          args: -pep_refine simple -pep_refine_scheme explicit
310:       test:
311:          suffix: 9_multiple_mbe
312:          args: -pep_refine multiple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
313:          requires: !single
314:       test:
315:          suffix: 9_multiple_explicit
316:          args: -pep_refine multiple -pep_refine_scheme explicit
317:          requires: !single

319:    test:
320:       suffix: 10
321:       nsize: 4
322:       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
323:       output_file: output/test2_2.out

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

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

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

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

356: TEST*/