Actual source code: gun.c
slepc-3.17.0 2022-03-31
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 gun problem arises from model of a radio-frequency gun cavity, with
19: the complex nonlinear function
20: T(lambda) = K-lambda*M+i*lambda^(1/2)*W1+i*(lambda-108.8774^2)^(1/2)*W2
22: Data files can be downloaded from https://slepc.upv.es/datafiles
23: */
25: static char help[] = "Radio-frequency gun cavity.\n\n"
26: "The command line options are:\n"
27: "-K <filename1> -M <filename2> -W1 <filename3> -W2 <filename4>, where filename1,..,filename4 are files containing the matrices in PETSc binary form defining the GUN problem.\n\n";
29: #include <slepcnep.h>
31: #define NMAT 4
32: #define SIGMA 108.8774
34: int main(int argc,char **argv)
35: {
36: Mat A[NMAT]; /* problem matrices */
37: FN f[NMAT]; /* functions to define the nonlinear operator */
38: FN ff[2]; /* auxiliary functions to define the nonlinear operator */
39: NEP nep; /* nonlinear eigensolver context */
40: PetscBool terse,flg;
41: const char* string[NMAT]={"-K","-M","-W1","-W2"};
42: char filename[PETSC_MAX_PATH_LEN];
43: PetscScalar numer[2],sigma;
44: PetscInt i;
45: PetscViewer viewer;
47: SlepcInitialize(&argc,&argv,(char*)0,help);
49: PetscPrintf(PETSC_COMM_WORLD,"GUN problem\n\n");
50: #if !defined(PETSC_USE_COMPLEX)
51: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This example requires complex scalars!");
52: #endif
54: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55: Load the problem matrices
56: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58: for (i=0;i<NMAT;i++) {
59: PetscOptionsGetString(NULL,NULL,string[i],filename,sizeof(filename),&flg);
61: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
62: MatCreate(PETSC_COMM_WORLD,&A[i]);
63: MatSetFromOptions(A[i]);
64: MatLoad(A[i],viewer);
65: PetscViewerDestroy(&viewer);
66: }
68: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: Create the problem functions
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
72: /* f1=1 */
73: FNCreate(PETSC_COMM_WORLD,&f[0]);
74: FNSetType(f[0],FNRATIONAL);
75: numer[0] = 1.0;
76: FNRationalSetNumerator(f[0],1,numer);
78: /* f2=-lambda */
79: FNCreate(PETSC_COMM_WORLD,&f[1]);
80: FNSetType(f[1],FNRATIONAL);
81: numer[0] = -1.0; numer[1] = 0.0;
82: FNRationalSetNumerator(f[1],2,numer);
84: /* f3=i*sqrt(lambda) */
85: FNCreate(PETSC_COMM_WORLD,&f[2]);
86: FNSetType(f[2],FNSQRT);
87: FNSetScale(f[2],1.0,PETSC_i);
89: /* f4=i*sqrt(lambda-sigma^2) */
90: sigma = SIGMA*SIGMA;
91: FNCreate(PETSC_COMM_WORLD,&ff[0]);
92: FNSetType(ff[0],FNSQRT);
93: FNCreate(PETSC_COMM_WORLD,&ff[1]);
94: FNSetType(ff[1],FNRATIONAL);
95: numer[0] = 1.0; numer[1] = -sigma;
96: FNRationalSetNumerator(ff[1],2,numer);
97: FNCreate(PETSC_COMM_WORLD,&f[3]);
98: FNSetType(f[3],FNCOMBINE);
99: FNCombineSetChildren(f[3],FN_COMBINE_COMPOSE,ff[1],ff[0]);
100: FNSetScale(f[3],1.0,PETSC_i);
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Create the eigensolver and solve the problem
104: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106: NEPCreate(PETSC_COMM_WORLD,&nep);
107: NEPSetSplitOperator(nep,4,A,f,UNKNOWN_NONZERO_PATTERN);
108: NEPSetFromOptions(nep);
110: NEPSolve(nep);
112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: Display solution and clean up
114: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116: /* show detailed info unless -terse option is given by user */
117: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
118: if (terse) NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
119: else {
120: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
121: NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
122: NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
123: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
124: }
125: NEPDestroy(&nep);
126: for (i=0;i<NMAT;i++) {
127: MatDestroy(&A[i]);
128: FNDestroy(&f[i]);
129: }
130: for (i=0;i<2;i++) FNDestroy(&ff[i]);
131: SlepcFinalize();
132: return 0;
133: }
135: /*TEST
137: build:
138: requires: complex
140: test:
141: suffix: 1
142: args: -K ${DATAFILESPATH}/matrices/complex/gun_K.petsc -M ${DATAFILESPATH}/matrices/complex/gun_M.petsc -W1 ${DATAFILESPATH}/matrices/complex/gun_W1.petsc -W2 ${DATAFILESPATH}/matrices/complex/gun_W2.petsc -nep_type nleigs -rg_type polygon -rg_polygon_vertices 12500-1i,120500-1i,120500+30000i,70000+30000i -nep_target 65000 -nep_nev 24 -terse
143: requires: double complex datafilespath !defined(PETSC_USE_64BIT_INDICES)
144: timeoutfactor: 10
146: TEST*/