Actual source code: loaded_string.c

slepc-3.20.0 2023-09-29
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:    This example implements one of the problems found at
 12:        NLEVP: A Collection of Nonlinear Eigenvalue Problems,
 13:        The University of Manchester.
 14:    The details of the collection can be found at:
 15:        [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
 16:            Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.

 18:    The loaded_string problem is a rational eigenvalue problem for the
 19:    finite element model of a loaded vibrating string.
 20: */

 22: static char help[] = "Finite element model of a loaded vibrating string.\n\n"
 23:   "The command line options are:\n"
 24:   "  -n <n>, dimension of the matrices.\n"
 25:   "  -kappa <kappa>, stiffness of elastic spring.\n"
 26:   "  -mass <m>, mass of the attached load.\n\n";

 28: #include <slepcnep.h>

 30: #define NMAT 3

 32: int main(int argc,char **argv)
 33: {
 34:   Mat            A[NMAT];         /* problem matrices */
 35:   FN             f[NMAT];         /* functions to define the nonlinear operator */
 36:   NEP            nep;             /* nonlinear eigensolver context */
 37:   PetscInt       n=100,Istart,Iend,i;
 38:   PetscReal      kappa=1.0,m=1.0;
 39:   PetscScalar    sigma,numer[2],denom[2];
 40:   PetscBool      terse;

 42:   PetscFunctionBeginUser;
 43:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

 45:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 46:   PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL));
 47:   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mass",&m,NULL));
 48:   sigma = kappa/m;
 49:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Loaded vibrating string, n=%" PetscInt_FMT " kappa=%g m=%g\n\n",n,(double)kappa,(double)m));

 51:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 52:                        Build the problem matrices
 53:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 55:   /* initialize matrices */
 56:   for (i=0;i<NMAT;i++) {
 57:     PetscCall(MatCreate(PETSC_COMM_WORLD,&A[i]));
 58:     PetscCall(MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n));
 59:     PetscCall(MatSetFromOptions(A[i]));
 60:     PetscCall(MatSetUp(A[i]));
 61:   }
 62:   PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));

 64:   /* A0 */
 65:   for (i=Istart;i<Iend;i++) {
 66:     PetscCall(MatSetValue(A[0],i,i,(i==n-1)?1.0*n:2.0*n,INSERT_VALUES));
 67:     if (i>0) PetscCall(MatSetValue(A[0],i,i-1,-1.0*n,INSERT_VALUES));
 68:     if (i<n-1) PetscCall(MatSetValue(A[0],i,i+1,-1.0*n,INSERT_VALUES));
 69:   }

 71:   /* A1 */
 72:   for (i=Istart;i<Iend;i++) {
 73:     PetscCall(MatSetValue(A[1],i,i,(i==n-1)?2.0/(6.0*n):4.0/(6.0*n),INSERT_VALUES));
 74:     if (i>0) PetscCall(MatSetValue(A[1],i,i-1,1.0/(6.0*n),INSERT_VALUES));
 75:     if (i<n-1) PetscCall(MatSetValue(A[1],i,i+1,1.0/(6.0*n),INSERT_VALUES));
 76:   }

 78:   /* A2 */
 79:   if (Istart<=n-1 && n-1<Iend) PetscCall(MatSetValue(A[2],n-1,n-1,kappa,INSERT_VALUES));

 81:   /* assemble matrices */
 82:   for (i=0;i<NMAT;i++) PetscCall(MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY));
 83:   for (i=0;i<NMAT;i++) PetscCall(MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY));

 85:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 86:                        Create the problem functions
 87:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 89:   /* f1=1 */
 90:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[0]));
 91:   PetscCall(FNSetType(f[0],FNRATIONAL));
 92:   numer[0] = 1.0;
 93:   PetscCall(FNRationalSetNumerator(f[0],1,numer));

 95:   /* f2=-lambda */
 96:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[1]));
 97:   PetscCall(FNSetType(f[1],FNRATIONAL));
 98:   numer[0] = -1.0; numer[1] = 0.0;
 99:   PetscCall(FNRationalSetNumerator(f[1],2,numer));

101:   /* f3=lambda/(lambda-sigma) */
102:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[2]));
103:   PetscCall(FNSetType(f[2],FNRATIONAL));
104:   numer[0] = 1.0; numer[1] = 0.0;
105:   denom[0] = 1.0; denom[1] = -sigma;
106:   PetscCall(FNRationalSetNumerator(f[2],2,numer));
107:   PetscCall(FNRationalSetDenominator(f[2],2,denom));

109:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110:                 Create the eigensolver and solve the problem
111:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

113:   PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
114:   PetscCall(NEPSetSplitOperator(nep,3,A,f,SUBSET_NONZERO_PATTERN));
115:   PetscCall(NEPSetProblemType(nep,NEP_RATIONAL));
116:   PetscCall(NEPSetFromOptions(nep));
117:   PetscCall(NEPSolve(nep));

119:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120:                     Display solution and clean up
121:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

123:   /* show detailed info unless -terse option is given by user */
124:   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
125:   if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
126:   else {
127:     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
128:     PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD));
129:     PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
130:     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
131:   }
132:   PetscCall(NEPDestroy(&nep));
133:   for (i=0;i<NMAT;i++) {
134:     PetscCall(MatDestroy(&A[i]));
135:     PetscCall(FNDestroy(&f[i]));
136:   }
137:   PetscCall(SlepcFinalize());
138:   return 0;
139: }

141: /*TEST

143:    test:
144:       suffix: 1
145:       args: -nep_type rii -nep_target 4 -terse
146:       requires: !single
147:       filter: sed -e "s/[+-]0\.0*i//g"

149:    testset:
150:       args: -nep_type interpol -rg_type interval -rg_interval_endpoints 5,700,-.1,.1 -nep_nev 7 -nep_target 5 -nep_interpol_interpolation_degree 12 -nep_refine simple -terse
151:       requires: !single
152:       output_file: output/loaded_string_2.out
153:       test:
154:          suffix: 2
155:          args: -nep_refine_scheme {{schur explicit}}
156:       test:
157:          suffix: 2_mbe
158:          args: -nep_refine_scheme mbe -nep_refine_ksp_type preonly -nep_refine_pc_type lu

160:    testset:
161:       nsize: 2
162:       args: -nep_type interpol -rg_type interval -rg_interval_endpoints 5,700,-.1,.1 -nep_nev 7 -nep_target 5 -nep_interpol_interpolation_degree 12 -nep_refine simple -nep_refine_partitions 2 -nep_interpol_st_ksp_type bcgs -nep_interpol_st_pc_type bjacobi -terse
163:       requires: !single
164:       output_file: output/loaded_string_2.out
165:       timeoutfactor: 2
166:       test:
167:          suffix: 3_explicit
168:          args: -nep_refine_scheme explicit
169:       test:
170:          suffix: 3_mbe
171:          args: -nep_refine_scheme mbe -nep_refine_ksp_type preonly -nep_refine_pc_type cholesky

173:    test:
174:       suffix: 4
175:       nsize: 4
176:       args: -nep_type interpol -rg_type interval -rg_interval_endpoints 5,700,-.1,.1 -nep_nev 7 -nep_target 5 -nep_interpol_interpolation_degree 10 -nep_refine simple -nep_refine_partitions 2 -nep_refine_scheme explicit -nep_interpol_st_ksp_type bcgs -nep_interpol_st_pc_type bjacobi -terse -log_exclude nep,pep,fn
177:       requires: !single
178:       output_file: output/loaded_string_2.out
179:       timeoutfactor: 4

181:    test:
182:       suffix: 5
183:       args: -nep_type nleigs -rg_type interval -rg_interval_endpoints 4,700,-.1,.1 -nep_nev 8 -nep_target 5 -terse
184:       filter: sed -e "s/[+-]0\.0*i//g"
185:       requires: !single

187:    test:
188:       suffix: 6
189:       args: -nep_type nleigs -rg_type interval -rg_interval_endpoints 100,700 -nep_nev 5 -nep_tol 1e-9 -nep_target 140 -nep_nleigs_interpolation_degree 15 -nep_general -terse
190:       requires: !complex !single

192:    test:
193:       suffix: 6_complex
194:       args: -nep_type nleigs -rg_type interval -rg_interval_endpoints 100,700,-.1,.1 -nep_nev 5 -nep_tol 1e-9 -nep_target 140 -nep_nleigs_interpolation_degree 15 -nep_general -terse
195:       filter: sed -e "s/[+-]0\.0*i//g"
196:       requires: complex !single
197:       output_file: output/loaded_string_6.out

199:    test:
200:       suffix: 7
201:       args: -nep_type interpol -rg_type interval -rg_interval_endpoints 5,700 -nep_nev 5 -nep_target 100 -nep_interpol_interpolation_degree 20 -nep_ncv 20 -n 20 -nep_refine simple -nep_refine_its 1 -terse
202:       requires: !complex double

204:    test:
205:       suffix: 7_complex
206:       args: -nep_type interpol -rg_type interval -rg_interval_endpoints 5,700,-.1,.1 -nep_nev 5 -nep_target 100 -nep_interpol_interpolation_degree 20 -nep_ncv 20 -n 20 -nep_refine simple -nep_refine_its 1 -terse
207:       requires: complex double
208:       output_file: output/loaded_string_7.out

210:    testset:
211:       args: -nep_target 10 -nep_nev 3 -nep_tol 5e-10 -terse
212:       requires: !single
213:       output_file: output/loaded_string_8.out
214:       filter: sed -e "s/[+-]0\.0*i//g"
215:       test:
216:          suffix: 8
217:          args: -nep_type {{rii slp narnoldi}}
218:       test:
219:          suffix: 8_rii_thres
220:          args: -nep_type rii -nep_rii_deflation_threshold 5e-10
221:       test:
222:          suffix: 8_slp_thres
223:          args: -nep_type slp -nep_slp_deflation_threshold 5e-10

225:       test:
226:          suffix: 8_slp_two_thres
227:          args: -nep_type slp -nep_slp_deflation_threshold 5e-10 -nep_two_sided

229:    test:
230:       suffix: 9
231:       args: -nep_type ciss -rg_type ellipse -rg_ellipse_center 500 -rg_ellipse_radius 500 -rg_ellipse_vscale .1 -nep_ciss_moments 4 -nep_ciss_blocksize 5 -terse
232:       requires: complex double

234: TEST*/