Actual source code: test2.c

slepc-3.15.0 2021-03-31
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: */
 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: {
 39:   PetscInt       i,nconv;
 40:   Mat            A;
 41:   Vec            xr,xi;
 42:   PetscReal      error=0.0,normr;
 43: #if !defined(PETSC_USE_COMPLEX)
 44:   PetscReal      normi;
 45: #endif

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

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

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

 86:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 87:   PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL);
 88:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 89:   PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL);
 90:   PetscOptionsGetBool(NULL,NULL,"-initv",&initv,NULL);
 91:   PetscOptionsGetBool(NULL,NULL,"-skipnorm",&skipnorm,NULL);

 93:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 94:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 95:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 97:   /* K is a tridiagonal */
 98:   MatCreate(PETSC_COMM_WORLD,&K);
 99:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
100:   MatSetFromOptions(K);
101:   MatSetUp(K);

103:   MatGetOwnershipRange(K,&Istart,&Iend);
104:   for (i=Istart;i<Iend;i++) {
105:     if (i>0) {
106:       MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
107:     }
108:     MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
109:     if (i<n-1) {
110:       MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
111:     }
112:   }

114:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
115:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

117:   /* C is a tridiagonal */
118:   MatCreate(PETSC_COMM_WORLD,&C);
119:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
120:   MatSetFromOptions(C);
121:   MatSetUp(C);

123:   MatGetOwnershipRange(C,&Istart,&Iend);
124:   for (i=Istart;i<Iend;i++) {
125:     if (i>0) {
126:       MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
127:     }
128:     MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
129:     if (i<n-1) {
130:       MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
131:     }
132:   }

134:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
135:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

137:   /* M is a diagonal matrix */
138:   MatCreate(PETSC_COMM_WORLD,&M);
139:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
140:   MatSetFromOptions(M);
141:   MatSetUp(M);
142:   MatGetOwnershipRange(M,&Istart,&Iend);
143:   for (i=Istart;i<Iend;i++) {
144:     MatSetValue(M,i,i,mu,INSERT_VALUES);
145:   }
146:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
147:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

149:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150:                 Create the eigensolver and set various options
151:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

153:   PEPCreate(PETSC_COMM_WORLD,&pep);
154:   A[0] = K; A[1] = C; A[2] = M;
155:   PEPSetOperators(pep,3,A);
156:   PEPSetProblemType(pep,PEP_GENERAL);
157:   PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT);
158:   if (initv) { /* initial vector */
159:     MatCreateVecs(K,&IV[0],NULL);
160:     VecSetValue(IV[0],0,-1.0,INSERT_VALUES);
161:     VecSetValue(IV[0],1,0.5,INSERT_VALUES);
162:     VecAssemblyBegin(IV[0]);
163:     VecAssemblyEnd(IV[0]);
164:     MatCreateVecs(K,&IV[1],NULL);
165:     VecSetValue(IV[1],0,4.0,INSERT_VALUES);
166:     VecSetValue(IV[1],2,1.5,INSERT_VALUES);
167:     VecAssemblyBegin(IV[1]);
168:     VecAssemblyEnd(IV[1]);
169:     PEPSetInitialSpace(pep,2,IV);
170:     VecDestroy(&IV[0]);
171:     VecDestroy(&IV[1]);
172:   }
173:   PEPSetFromOptions(pep);

175:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176:                       Solve the eigensystem
177:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

179:   PEPSolve(pep);
180:   PEPGetDimensions(pep,&nev,NULL,NULL);
181:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

183:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184:                     Display solution and clean up
185:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

187:   PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
188:   if (!skipnorm) { CheckNormalizedVectors(pep); }
189:   PEPDestroy(&pep);
190:   MatDestroy(&M);
191:   MatDestroy(&C);
192:   MatDestroy(&K);
193:   SlepcFinalize();
194:   return ierr;
195: }

197: /*TEST

199:    testset:
200:       args: -pep_nev 4 -initv
201:       requires: !single
202:       output_file: output/test2_1.out
203:       test:
204:          suffix: 1
205:          args: -pep_type {{toar linear}}
206:       test:
207:          suffix: 1_toar_mgs
208:          args: -pep_type toar -bv_orthog_type mgs
209:       test:
210:          suffix: 1_qarnoldi
211:          args: -pep_type qarnoldi -bv_orthog_refine never
212:       test:
213:          suffix: 1_linear_gd
214:          args: -pep_type linear -pep_linear_eps_type gd -pep_linear_explicitmatrix

216:    testset:
217:       args: -pep_target -0.43 -pep_nev 4 -pep_ncv 20 -st_type sinvert
218:       output_file: output/test2_2.out
219:       test:
220:          suffix: 2
221:          args: -pep_type {{toar linear}}
222:       test:
223:          suffix: 2_toar_scaleboth
224:          args: -pep_type toar -pep_scale both
225:       test:
226:          suffix: 2_toar_transform
227:          args: -pep_type toar -st_transform
228:       test:
229:          suffix: 2_qarnoldi
230:          args: -pep_type qarnoldi -bv_orthog_refine always
231:       test:
232:          suffix: 2_linear_explicit
233:          args: -pep_type linear -pep_linear_explicitmatrix -pep_linear_linearization 0,1
234:       test:
235:          suffix: 2_linear_explicit_her
236:          args: -pep_type linear -pep_linear_explicitmatrix -pep_hermitian -pep_linear_linearization 0,1
237:       test:
238:          suffix: 2_stoar
239:          args: -pep_type stoar -pep_hermitian
240:       test:
241:          suffix: 2_jd
242:          args: -pep_type jd -st_type precond -pep_max_it 200 -pep_ncv 24
243:          requires: !single

245:    test:
246:       suffix: 3
247:       args: -pep_nev 12 -pep_extract {{none norm residual structured}} -pep_monitor_cancel
248:       requires: !single

250:    testset:
251:       args: -st_type sinvert -pep_target -0.43 -pep_nev 4
252:       output_file: output/test2_2.out
253:       test:
254:          suffix: 4_schur
255:          args: -pep_refine simple -pep_refine_scheme schur
256:       test:
257:          suffix: 4_mbe
258:          args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
259:       test:
260:          suffix: 4_explicit
261:          args: -pep_refine simple -pep_refine_scheme explicit
262:       test:
263:          suffix: 4_multiple_schur
264:          args: -pep_refine multiple -pep_refine_scheme schur
265:          requires: !single
266:       test:
267:          suffix: 4_multiple_mbe
268:          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
269:       test:
270:          suffix: 4_multiple_explicit
271:          args: -pep_refine multiple -pep_refine_scheme explicit
272:          requires: !single

274:    test:
275:       suffix: 5
276:       nsize: 2
277:       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
278:       output_file: output/test2_2.out

280:    test:
281:       suffix: 6
282:       args: -pep_type linear -pep_nev 12 -pep_extract {{none norm residual}}
283:       requires: !single
284:       output_file: output/test2_3.out

286:    test:
287:       suffix: 7
288:       args: -pep_nev 12 -pep_extract {{none norm residual structured}} -pep_refine multiple
289:       requires: !single
290:       output_file: output/test2_3.out

292:    testset:
293:       args: -st_type sinvert -pep_target -0.43 -pep_nev 4 -st_transform
294:       output_file: output/test2_2.out
295:       test:
296:          suffix: 8_schur
297:          args: -pep_refine simple -pep_refine_scheme schur
298:       test:
299:          suffix: 8_mbe
300:          args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
301:       test:
302:          suffix: 8_explicit
303:          args: -pep_refine simple -pep_refine_scheme explicit
304:       test:
305:          suffix: 8_multiple_schur
306:          args: -pep_refine multiple -pep_refine_scheme schur
307:       test:
308:          suffix: 8_multiple_mbe
309:          args: -pep_refine multiple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
310:       test:
311:          suffix: 8_multiple_explicit
312:          args: -pep_refine multiple -pep_refine_scheme explicit

314:    testset:
315:       nsize: 2
316:       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
317:       output_file: output/test2_2.out
318:       test:
319:          suffix: 9_mbe
320:          args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
321:       test:
322:          suffix: 9_explicit
323:          args: -pep_refine simple -pep_refine_scheme explicit
324:       test:
325:          suffix: 9_multiple_mbe
326:          args: -pep_refine multiple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
327:          requires: !single
328:       test:
329:          suffix: 9_multiple_explicit
330:          args: -pep_refine multiple -pep_refine_scheme explicit
331:          requires: !single

333:    test:
334:       suffix: 10
335:       nsize: 4
336:       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
337:       output_file: output/test2_2.out

339:    testset:
340:       args: -pep_nev 4 -initv -mat_type aijcusparse
341:       output_file: output/test2_1.out
342:       requires: cuda !single
343:       test:
344:          suffix: 11_cuda
345:          args: -pep_type {{toar linear}}
346:       test:
347:          suffix: 11_cuda_qarnoldi
348:          args: -pep_type qarnoldi -bv_orthog_refine never
349:       test:
350:          suffix: 11_cuda_linear_gd
351:          args: -pep_type linear -pep_linear_eps_type gd -pep_linear_explicitmatrix

353:    test:
354:       suffix: 12
355:       nsize: 2
356:       args: -pep_type jd -ds_parallel synchronized -pep_target -0.43 -pep_nev 4 -pep_ncv 20
357:       requires: !single

359:    test:
360:       suffix: 13
361:       args: -pep_nev 12 -pep_view_values draw -pep_monitor draw::draw_lg
362:       requires: x !single
363:       output_file: output/test2_3.out

365: TEST*/