Actual source code: test1.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: */

 11: static char help[] = "Tests B-orthonormality of eigenvectors in a GHEP problem.\n\n";

 13: #include <slepceps.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat               A,B;        /* matrices */
 18:   EPS               eps;        /* eigenproblem solver context */
 19:   ST                st;
 20:   Vec               *X,v;
 21:   PetscReal         lev,tol=1000*PETSC_MACHINE_EPSILON;
 22:   PetscInt          N,n=45,m,Istart,Iend,II,i,j,nconv;
 23:   PetscBool         flag;
 24:   EPSPowerShiftType variant;
 25:   PetscErrorCode    ierr;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 30:   if (!flag) m=n;
 31:   N = n*m;
 32:   PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized Symmetric Eigenproblem, N=%D (%Dx%D grid)\n\n",N,n,m);

 34:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 35:      Compute the matrices that define the eigensystem, Ax=kBx
 36:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 38:   MatCreate(PETSC_COMM_WORLD,&A);
 39:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 40:   MatSetFromOptions(A);
 41:   MatSetUp(A);

 43:   MatCreate(PETSC_COMM_WORLD,&B);
 44:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
 45:   MatSetFromOptions(B);
 46:   MatSetUp(B);

 48:   MatGetOwnershipRange(A,&Istart,&Iend);
 49:   for (II=Istart;II<Iend;II++) {
 50:     i = II/n; j = II-i*n;
 51:     if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
 52:     if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
 53:     if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
 54:     if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
 55:     MatSetValue(A,II,II,4.0,INSERT_VALUES);
 56:     MatSetValue(B,II,II,2.0/PetscLogScalar(II+2),INSERT_VALUES);
 57:   }

 59:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 60:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 61:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 62:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 63:   MatCreateVecs(B,&v,NULL);

 65:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 66:                 Create the eigensolver and set various options
 67:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 69:   EPSCreate(PETSC_COMM_WORLD,&eps);
 70:   EPSSetOperators(eps,A,B);
 71:   EPSSetProblemType(eps,EPS_GHEP);
 72:   EPSSetTolerances(eps,tol,PETSC_DEFAULT);
 73:   EPSSetConvergenceTest(eps,EPS_CONV_NORM);
 74:   EPSSetFromOptions(eps);

 76:   /* illustrate how to extract parameters from specific solver types */
 77:   PetscObjectTypeCompare((PetscObject)eps,EPSPOWER,&flag);
 78:   if (flag) {
 79:     EPSGetST(eps,&st);
 80:     PetscObjectTypeCompare((PetscObject)st,STSHIFT,&flag);
 81:     if (flag) {
 82:       EPSPowerGetShiftType(eps,&variant);
 83:       PetscPrintf(PETSC_COMM_WORLD,"Type of shifts used during power iteration: %s\n",EPSPowerShiftTypes[variant]);
 84:     }
 85:   }

 87:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 88:                       Solve the eigensystem
 89:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 91:   EPSSolve(eps);

 93:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 94:                     Display solution and clean up
 95:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 97:   EPSGetTolerances(eps,&tol,NULL);
 98:   EPSErrorView(eps,EPS_ERROR_BACKWARD,NULL);
 99:   EPSGetConverged(eps,&nconv);
100:   if (nconv>1) {
101:     VecDuplicateVecs(v,nconv,&X);
102:     for (i=0;i<nconv;i++) {
103:       EPSGetEigenvector(eps,i,X[i],NULL);
104:     }
105:     VecCheckOrthogonality(X,nconv,NULL,nconv,B,NULL,&lev);
106:     if (lev<10*tol) {
107:       PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n");
108:     } else {
109:       PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)lev);
110:     }
111:     VecDestroyVecs(nconv,&X);
112:   }

114:   EPSDestroy(&eps);
115:   MatDestroy(&A);
116:   MatDestroy(&B);
117:   VecDestroy(&v);
118:   SlepcFinalize();
119:   return ierr;
120: }

122: /*TEST

124:    testset:
125:       args: -n 18 -eps_nev 4 -eps_max_it 1500
126:       requires: !single
127:       output_file: output/test1_1.out
128:       test:
129:          suffix: 1
130:          args: -eps_type {{krylovschur subspace arnoldi gd jd lapack}}
131:       test:
132:          suffix: 1_ks_nopurify
133:          args: -eps_purify 0
134:       test:
135:          suffix: 1_ks_trueres
136:          args: -eps_true_residual
137:       test:
138:          suffix: 1_ks_sinvert
139:          args: -st_type sinvert -eps_target 22
140:       test:
141:          suffix: 1_ks_cayley
142:          args: -st_type cayley -eps_target 22
143:       test:
144:          suffix: 1_lanczos
145:          args: -eps_type lanczos -eps_lanczos_reorthog full
146:       test:
147:          suffix: 1_gd2
148:          args: -eps_type gd -eps_gd_double_expansion
149:       test:
150:          suffix: 1_ciss
151:          args: -eps_type ciss -rg_interval_endpoints 20.8,22 -eps_largest_real
152:       test:
153:          suffix: 1_lobpcg
154:          args: -eps_type lobpcg -st_shift 22 -eps_largest_real

156:    test:
157:       suffix: 2
158:       requires: !single !complex
159:       args: -eps_interval .1,1.1 -eps_tol 1e-10 -st_type sinvert -st_ksp_type preonly -st_pc_type cholesky

161:    test:
162:       suffix: 3
163:       requires: !single
164:       args: -n 18 -eps_type power -eps_nev 3

166:    test:
167:       suffix: 4
168:       requires: !single
169:       args: -n 18 -eps_type power -eps_nev 3 -st_type sinvert -eps_target 1.149 -eps_power_shift_type {{constant rayleigh wilkinson}}

171:    testset:
172:       args: -n 18 -eps_nev 3 -eps_smallest_real -eps_max_it 500 -st_pc_type icc
173:       output_file: output/test1_5.out
174:       test:
175:          suffix: 5_rqcg
176:          args: -eps_type rqcg
177:       test:
178:          suffix: 5_lobpcg
179:          args: -eps_type lobpcg -eps_lobpcg_blocksize 3
180:          requires: !single

182:    testset:
183:       args: -n 18 -eps_nev 12 -eps_mpd 8 -eps_max_it 3000
184:       requires: !single
185:       output_file: output/test1_6.out
186:       test:
187:          suffix: 6
188:          args: -eps_type {{krylovschur subspace arnoldi gd}}
189:       test:
190:          suffix: 6_lanczos
191:          args: -eps_type lanczos -eps_lanczos_reorthog full

193:    testset:
194:       args: -n 18 -eps_nev 4 -eps_max_it 1500 -mat_type aijcusparse
195:       requires: cuda
196:       output_file: output/test1_1.out
197:       test:
198:          suffix: 7
199:          args: -eps_type {{krylovschur subspace arnoldi gd jd}}
200:       test:
201:          suffix: 7_ks_sinvert
202:          args: -st_type sinvert -eps_target 22
203:       test:
204:          suffix: 7_lanczos
205:          args: -eps_type lanczos -eps_lanczos_reorthog full
206:       test:
207:          suffix: 7_ciss
208:          args: -eps_type ciss -rg_interval_endpoints 20.8,22 -eps_largest_real

210:    testset:
211:       args: -n 18 -eps_nev 3 -eps_smallest_real -eps_max_it 500 -st_pc_type sor -mat_type aijcusparse
212:       requires: cuda
213:       output_file: output/test1_5.out
214:       test:
215:          suffix: 8_rqcg
216:          args: -eps_type rqcg
217:       test:
218:          suffix: 8_lobpcg
219:          args: -eps_type lobpcg -eps_lobpcg_blocksize 3

221:    testset:
222:       nsize: 2
223:       args: -n 18 -eps_nev 7 -eps_ncv 32 -ds_parallel synchronized
224:       requires: cuda
225:       filter: grep -v "orthogonality"
226:       output_file: output/test1_9.out
227:       test:
228:          suffix: 9_ks_ghep
229:          args: -eps_gen_hermitian -st_pc_type redundant -st_type sinvert
230:       test:
231:          suffix: 9_ks_gnhep
232:          args: -eps_gen_non_hermitian -st_pc_type redundant -st_type sinvert
233:       test:
234:          suffix: 9_ks_ghiep
235:          args: -eps_gen_indefinite -st_pc_type redundant -st_type sinvert
236:       test:
237:          suffix: 9_lobpcg_ghep
238:          args: -eps_gen_hermitian -eps_type lobpcg -eps_max_it 200 -eps_lobpcg_blocksize 6
239:       test:
240:          suffix: 9_jd_gnhep
241:          args: -eps_gen_non_hermitian -eps_type jd -eps_target 0
242:          requires: !complex

244: TEST*/