Actual source code: gun.c

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */
 21: /*
 22:    This example implements one of the problems found at
 23:        NLEVP: A Collection of Nonlinear Eigenvalue Problems,
 24:        The University of Manchester.
 25:    The details of the collection can be found at:
 26:        [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
 27:            Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.

 29:    The gun problem arises from model of a radio-frequency gun cavity, with
 30:    the complex nonlinear function
 31:    T(lambda) = K-lambda*M+i*lambda^(1/2)*W1+i*(lambda-108.8774^2)^(1/2)*W2

 33:    Data files can be downloaded from http://slepc.upv.es/datafiles
 34: */

 36: static char help[] = "Radio-frequency gun cavity.\n\n"
 37:   "The command line options are:\n"
 38:   "-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";

 40: #include <slepcnep.h>

 42: #define NMAT 4
 43: #define SIGMA 108.8774

 45: PetscErrorCode ComputeSingularities(NEP,PetscInt*,PetscScalar*,void*);

 49: int main(int argc,char **argv)
 50: {
 52:   Mat            A[NMAT];         /* problem matrices */
 53:   FN             f[NMAT];         /* functions to define the nonlinear operator */
 54:   FN             ff[2];           /* auxiliary functions to define the nonlinear operator */
 55:   NEP            nep;             /* nonlinear eigensolver context */
 56:   PetscBool      terse,flg;
 57:   const char*    string[NMAT]={"-K","-M","-W1","-W2"};
 58:   char           filename[PETSC_MAX_PATH_LEN];
 59:   PetscScalar    numer[2],sigma;
 60:   PetscInt       i;
 61:   PetscViewer    viewer;

 63:   SlepcInitialize(&argc,&argv,(char*)0,help);

 65:   PetscPrintf(PETSC_COMM_WORLD,"GUN problem\n\n");
 66: #if !defined(PETSC_USE_COMPLEX)
 67:   SETERRQ(PETSC_COMM_SELF,1,"This example requires complex scalars!");
 68: #endif

 70:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 71:                        Load the problem matrices 
 72:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 74:   for (i=0;i<NMAT;i++) {
 75:     PetscOptionsGetString(NULL,NULL,string[i],filename,PETSC_MAX_PATH_LEN,&flg);
 76:     if (!flg) SETERRQ1(PETSC_COMM_WORLD,1,"Must indicate a filename with the %s option",string[i]);
 77:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
 78:     MatCreate(PETSC_COMM_WORLD,&A[i]);
 79:     MatSetFromOptions(A[i]);
 80:     MatLoad(A[i],viewer);
 81:     PetscViewerDestroy(&viewer);
 82:   }

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

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

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

100:   /* f3=i*sqrt(lambda) */
101:   FNCreate(PETSC_COMM_WORLD,&f[2]);
102:   FNSetType(f[2],FNSQRT);
103:   FNSetScale(f[2],1.0,PETSC_i);

105:   /* f4=i*sqrt(lambda-sigma^2) */
106:   sigma = SIGMA*SIGMA;
107:   FNCreate(PETSC_COMM_WORLD,&ff[0]);
108:   FNSetType(ff[0],FNSQRT);
109:   FNCreate(PETSC_COMM_WORLD,&ff[1]);
110:   FNSetType(ff[1],FNRATIONAL);
111:   numer[0] = 1.0; numer[1] = -sigma;
112:   FNRationalSetNumerator(ff[1],2,numer);
113:   FNCreate(PETSC_COMM_WORLD,&f[3]);
114:   FNSetType(f[3],FNCOMBINE);
115:   FNCombineSetChildren(f[3],FN_COMBINE_COMPOSE,ff[1],ff[0]);
116:   FNSetScale(f[3],1.0,PETSC_i);

118:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
119:                 Create the eigensolver and solve the problem
120:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

122:   NEPCreate(PETSC_COMM_WORLD,&nep);
123:   NEPSetSplitOperator(nep,4,A,f,DIFFERENT_NONZERO_PATTERN);
124:   NEPSetFromOptions(nep);

126:   PetscObjectTypeCompare((PetscObject)nep,NEPNLEIGS,&flg);
127:   if (flg) {
128:     NEPNLEIGSSetSingularitiesFunction(nep,ComputeSingularities,NULL);
129:   }

131:   NEPSolve(nep);

133:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134:                     Display solution and clean up
135:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136:   
137:   /* show detailed info unless -terse option is given by user */
138:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
139:   if (terse) {
140:     NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
141:   } else {
142:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
143:     NEPReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
144:     NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
145:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
146:   }
147:   NEPDestroy(&nep);
148:   for (i=0;i<NMAT;i++) {
149:     MatDestroy(&A[i]);
150:     FNDestroy(&f[i]);
151:   }
152:   for (i=0;i<2;i++) {
153:     FNDestroy(&ff[i]);
154:   }
155:   SlepcFinalize();
156:   return ierr;
157: }

161: /*
162:    ComputeSingularities - Computes maxnp points (at most) in the complex plane where
163:    the function T(.) is not analytic.

165:    In this case, we discretize the singularity region (-inf,108.8774^2)~(-10e+12,-10e-12+108.8774^2) 
166: */
167: PetscErrorCode ComputeSingularities(NEP nep,PetscInt *maxnp,PetscScalar *xi,void *pt)
168: {
169:   PetscReal h;
170:   PetscInt  i;
171:   PetscReal   sigma,end;

174:   sigma = SIGMA*SIGMA;
175:   end = PetscLogReal(sigma);  
176:   h = (12.0+end)/(*maxnp-1);
177:   xi[0] = sigma;
178:   for (i=1;i<*maxnp;i++) xi[i] = -PetscPowReal(10,h*i)+sigma;
179:   return(0);
180: }