Actual source code: test2.c

slepc-3.11.2 2019-07-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            v0;

 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,&v0,NULL);
118:     VecSetValue(v0,0,-1.0,INSERT_VALUES);
119:     VecSetValue(v0,1,0.5,INSERT_VALUES);
120:     VecAssemblyBegin(v0);
121:     VecAssemblyEnd(v0);
122:     PEPSetInitialSpace(pep,1,&v0);
123:     VecDestroy(&v0);
124:   }
125:   PEPSetFromOptions(pep);

127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128:                       Solve the eigensystem
129:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

131:   PEPSolve(pep);
132:   PEPGetDimensions(pep,&nev,NULL,NULL);
133:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

135:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136:                     Display solution and clean up
137:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

139:   PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
140:   PEPDestroy(&pep);
141:   MatDestroy(&M);
142:   MatDestroy(&C);
143:   MatDestroy(&K);
144:   SlepcFinalize();
145:   return ierr;
146: }

148: /*TEST

150:    testset:
151:       args: -pep_nev 4 -initv
152:       requires: !single
153:       output_file: output/test2_1.out
154:       test:
155:          suffix: 1
156:          args: -pep_type {{toar linear}}
157:       test:
158:          suffix: 1_qarnoldi
159:          args: -pep_type qarnoldi -bv_orthog_refine never
160:       test:
161:          suffix: 1_linear_gd
162:          args: -pep_type linear -pep_linear_eps_type gd -pep_linear_explicitmatrix

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

194:    test:
195:       suffix: 3
196:       args: -pep_nev 12 -pep_extract {{none norm residual structured}} -pep_monitor_cancel
197:       requires: !single

199:    testset:
200:       args: -st_type sinvert -pep_target -0.43 -pep_nev 4
201:       output_file: output/test2_2.out
202:       test:
203:          suffix: 4_schur
204:          args: -pep_refine simple -pep_refine_scheme schur
205:       test:
206:          suffix: 4_mbe
207:          args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
208:       test:
209:          suffix: 4_explicit
210:          args: -pep_refine simple -pep_refine_scheme explicit
211:       test:
212:          suffix: 4_multiple_schur
213:          args: -pep_refine multiple -pep_refine_scheme schur
214:       test:
215:          suffix: 4_multiple_mbe
216:          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
217:       test:
218:          suffix: 4_multiple_explicit
219:          args: -pep_refine multiple -pep_refine_scheme explicit

221:    test:
222:       suffix: 5
223:       nsize: 2
224:       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
225:       output_file: output/test2_2.out

227:    test:
228:       suffix: 6
229:       args: -pep_type linear -pep_nev 12 -pep_extract {{none norm residual}}
230:       requires: !single
231:       output_file: output/test2_3.out

233:    test:
234:       suffix: 7
235:       args: -pep_nev 12 -pep_extract {{none norm residual structured}} -pep_refine multiple
236:       requires: !single
237:       output_file: output/test2_3.out

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

261:    testset:
262:       nsize: 2
263:       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
264:       output_file: output/test2_2.out
265:       requires: !single
266:       test:
267:          suffix: 9_mbe
268:          args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
269:       test:
270:          suffix: 9_explicit
271:          args: -pep_refine simple -pep_refine_scheme explicit
272:       test:
273:          suffix: 9_multiple_mbe
274:          args: -pep_refine multiple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
275:       test:
276:          suffix: 9_multiple_explicit
277:          args: -pep_refine multiple -pep_refine_scheme explicit

279:    test:
280:       suffix: 10
281:       nsize: 4
282:       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
283:       requires: double
284:       output_file: output/test2_2.out

286:    testset:
287:       args: -pep_nev 4 -initv -mat_type aijcusparse
288:       output_file: output/test2_1.out
289:       requires: cuda
290:       test:
291:          suffix: 11_cuda
292:          args: -pep_type {{toar linear}}
293:       test:
294:          suffix: 11_cuda_qarnoldi
295:          args: -pep_type qarnoldi -bv_orthog_refine never
296:       test:
297:          suffix: 11_cuda_linear_gd
298:          args: -pep_type linear -pep_linear_eps_type gd -pep_linear_explicitmatrix

300:    test:
301:       suffix: 12
302:       nsize: 2
303:       args: -pep_type jd -ds_parallel synchronized -pep_target -0.43 -pep_nev 4 -pep_ncv 20
304:       requires: !single

306: TEST*/