Actual source code: test2.c

slepc-3.12.0 2019-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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: int main(int argc,char **argv)
 34: {
 35:   Mat            M,C,K,A[3];      /* problem matrices */
 36:   PEP            pep;             /* polynomial eigenproblem solver context */
 38:   PetscInt       n=30,Istart,Iend,i,nev;
 39:   PetscScalar    mu=1.0,tau=10.0,kappa=5.0;
 40:   PetscBool      initv=PETSC_FALSE;
 41:   Vec            IV[2];

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

 45:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 46:   PetscOptionsGetScalar(NULL,NULL,"-mu",&mu,NULL);
 47:   PetscOptionsGetScalar(NULL,NULL,"-tau",&tau,NULL);
 48:   PetscOptionsGetScalar(NULL,NULL,"-kappa",&kappa,NULL);
 49:   PetscOptionsGetBool(NULL,NULL,"-initv",&initv,NULL);

 51:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 52:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 53:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 55:   /* K is a tridiagonal */
 56:   MatCreate(PETSC_COMM_WORLD,&K);
 57:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
 58:   MatSetFromOptions(K);
 59:   MatSetUp(K);

 61:   MatGetOwnershipRange(K,&Istart,&Iend);
 62:   for (i=Istart;i<Iend;i++) {
 63:     if (i>0) {
 64:       MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
 65:     }
 66:     MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
 67:     if (i<n-1) {
 68:       MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
 69:     }
 70:   }

 72:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 73:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 75:   /* C is a tridiagonal */
 76:   MatCreate(PETSC_COMM_WORLD,&C);
 77:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 78:   MatSetFromOptions(C);
 79:   MatSetUp(C);

 81:   MatGetOwnershipRange(C,&Istart,&Iend);
 82:   for (i=Istart;i<Iend;i++) {
 83:     if (i>0) {
 84:       MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
 85:     }
 86:     MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
 87:     if (i<n-1) {
 88:       MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
 89:     }
 90:   }

 92:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 93:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 95:   /* M is a diagonal matrix */
 96:   MatCreate(PETSC_COMM_WORLD,&M);
 97:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
 98:   MatSetFromOptions(M);
 99:   MatSetUp(M);
100:   MatGetOwnershipRange(M,&Istart,&Iend);
101:   for (i=Istart;i<Iend;i++) {
102:     MatSetValue(M,i,i,mu,INSERT_VALUES);
103:   }
104:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
105:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

107:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108:                 Create the eigensolver and set various options
109:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

111:   PEPCreate(PETSC_COMM_WORLD,&pep);
112:   A[0] = K; A[1] = C; A[2] = M;
113:   PEPSetOperators(pep,3,A);
114:   PEPSetProblemType(pep,PEP_GENERAL);
115:   PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT);
116:   if (initv) { /* initial vector */
117:     MatCreateVecs(K,&IV[0],NULL);
118:     VecSetValue(IV[0],0,-1.0,INSERT_VALUES);
119:     VecSetValue(IV[0],1,0.5,INSERT_VALUES);
120:     VecAssemblyBegin(IV[0]);
121:     VecAssemblyEnd(IV[0]);
122:     MatCreateVecs(K,&IV[1],NULL);
123:     VecSetValue(IV[1],0,4.0,INSERT_VALUES);
124:     VecSetValue(IV[1],2,1.5,INSERT_VALUES);
125:     VecAssemblyBegin(IV[1]);
126:     VecAssemblyEnd(IV[1]);
127:     PEPSetInitialSpace(pep,2,IV);
128:     VecDestroy(&IV[0]);
129:     VecDestroy(&IV[1]);
130:   }
131:   PEPSetFromOptions(pep);

133:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134:                       Solve the eigensystem
135:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

137:   PEPSolve(pep);
138:   PEPGetDimensions(pep,&nev,NULL,NULL);
139:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

141:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142:                     Display solution and clean up
143:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

145:   PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
146:   PEPDestroy(&pep);
147:   MatDestroy(&M);
148:   MatDestroy(&C);
149:   MatDestroy(&K);
150:   SlepcFinalize();
151:   return ierr;
152: }

154: /*TEST

156:    testset:
157:       args: -pep_nev 4 -initv
158:       requires: !single
159:       output_file: output/test2_1.out
160:       test:
161:          suffix: 1
162:          args: -pep_type {{toar linear}}
163:       test:
164:          suffix: 1_toar_mgs
165:          args: -pep_type toar -bv_orthog_type mgs
166:       test:
167:          suffix: 1_qarnoldi
168:          args: -pep_type qarnoldi -bv_orthog_refine never
169:       test:
170:          suffix: 1_linear_gd
171:          args: -pep_type linear -pep_linear_eps_type gd -pep_linear_explicitmatrix

173:    testset:
174:       args: -pep_target -0.43 -pep_nev 4 -pep_ncv 20 -st_type sinvert
175:       output_file: output/test2_2.out
176:       test:
177:          suffix: 2
178:          args: -pep_type {{toar linear}}
179:       test:
180:          suffix: 2_toar_scaleboth
181:          args: -pep_type toar -pep_scale both
182:       test:
183:          suffix: 2_toar_transform
184:          args: -pep_type toar -st_transform
185:       test:
186:          suffix: 2_qarnoldi
187:          args: -pep_type qarnoldi -bv_orthog_refine always
188:       test:
189:          suffix: 2_linear_explicit
190:          args: -pep_type linear -pep_linear_explicitmatrix -pep_linear_linearization 0,1
191:       test:
192:          suffix: 2_linear_explicit_her
193:          args: -pep_type linear -pep_linear_explicitmatrix -pep_hermitian -pep_linear_linearization 0,1
194:       test:
195:          suffix: 2_stoar
196:          args: -pep_type stoar -pep_hermitian
197:          requires: !single
198:       test:
199:          suffix: 2_jd
200:          args: -pep_type jd -st_type precond -pep_max_it 200
201:          requires: !single

203:    test:
204:       suffix: 3
205:       args: -pep_nev 12 -pep_extract {{none norm residual structured}} -pep_monitor_cancel
206:       requires: !single

208:    testset:
209:       args: -st_type sinvert -pep_target -0.43 -pep_nev 4
210:       output_file: output/test2_2.out
211:       test:
212:          suffix: 4_schur
213:          args: -pep_refine simple -pep_refine_scheme schur
214:       test:
215:          suffix: 4_mbe
216:          args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
217:       test:
218:          suffix: 4_explicit
219:          args: -pep_refine simple -pep_refine_scheme explicit
220:       test:
221:          suffix: 4_multiple_schur
222:          args: -pep_refine multiple -pep_refine_scheme schur
223:       test:
224:          suffix: 4_multiple_mbe
225:          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
226:       test:
227:          suffix: 4_multiple_explicit
228:          args: -pep_refine multiple -pep_refine_scheme explicit

230:    test:
231:       suffix: 5
232:       nsize: 2
233:       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
234:       output_file: output/test2_2.out

236:    test:
237:       suffix: 6
238:       args: -pep_type linear -pep_nev 12 -pep_extract {{none norm residual}}
239:       requires: !single
240:       output_file: output/test2_3.out

242:    test:
243:       suffix: 7
244:       args: -pep_nev 12 -pep_extract {{none norm residual structured}} -pep_refine multiple
245:       requires: !single
246:       output_file: output/test2_3.out

248:    testset:
249:       args: -st_type sinvert -pep_target -0.43 -pep_nev 4 -st_transform
250:       output_file: output/test2_2.out
251:       test:
252:          suffix: 8_schur
253:          args: -pep_refine simple -pep_refine_scheme schur
254:       test:
255:          suffix: 8_mbe
256:          args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
257:       test:
258:          suffix: 8_explicit
259:          args: -pep_refine simple -pep_refine_scheme explicit
260:       test:
261:          suffix: 8_multiple_schur
262:          args: -pep_refine multiple -pep_refine_scheme schur
263:       test:
264:          suffix: 8_multiple_mbe
265:          args: -pep_refine multiple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
266:       test:
267:          suffix: 8_multiple_explicit
268:          args: -pep_refine multiple -pep_refine_scheme explicit

270:    testset:
271:       nsize: 2
272:       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
273:       output_file: output/test2_2.out
274:       requires: !single
275:       test:
276:          suffix: 9_mbe
277:          args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
278:       test:
279:          suffix: 9_explicit
280:          args: -pep_refine simple -pep_refine_scheme explicit
281:       test:
282:          suffix: 9_multiple_mbe
283:          args: -pep_refine multiple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
284:       test:
285:          suffix: 9_multiple_explicit
286:          args: -pep_refine multiple -pep_refine_scheme explicit

288:    test:
289:       suffix: 10
290:       nsize: 4
291:       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
292:       requires: double
293:       output_file: output/test2_2.out

295:    testset:
296:       args: -pep_nev 4 -initv -mat_type aijcusparse
297:       output_file: output/test2_1.out
298:       requires: cuda !single
299:       test:
300:          suffix: 11_cuda
301:          args: -pep_type {{toar linear}}
302:       test:
303:          suffix: 11_cuda_qarnoldi
304:          args: -pep_type qarnoldi -bv_orthog_refine never
305:       test:
306:          suffix: 11_cuda_linear_gd
307:          args: -pep_type linear -pep_linear_eps_type gd -pep_linear_explicitmatrix

309:    test:
310:       suffix: 12
311:       nsize: 2
312:       args: -pep_type jd -ds_parallel synchronized -pep_target -0.43 -pep_nev 4 -pep_ncv 20
313:       requires: !single

315:    test:
316:       suffix: 13
317:       args: -pep_nev 12 -pep_view_values draw -pep_monitor_lg
318:       requires: x !single
319:       output_file: output/test2_3.out

321: TEST*/