Actual source code: ex27.c

slepc-3.8.2 2017-12-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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[] = "Simple nonlinear eigenproblem using the NLEIGS solver.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = matrix dimension.\n"
 14:   "  -split <0/1>, to select the split form in the problem definition (enabled by default)\n";

 16: /*
 17:    Solve T(lambda)x=0 using NLEIGS solver
 18:       with T(lambda) = -D+sqrt(lambda)*I
 19:       where D is the Laplacian operator in 1 dimension
 20:       and with the interpolation interval [.01,16]
 21: */

 23: #include <slepcnep.h>

 25: /*
 26:    User-defined routines
 27: */
 28: PetscErrorCode FormFunction(NEP,PetscScalar,Mat,Mat,void*);
 29: PetscErrorCode ComputeSingularities(NEP,PetscInt*,PetscScalar*,void*);

 31: int main(int argc,char **argv)
 32: {
 33:   NEP            nep;             /* nonlinear eigensolver context */
 34:   Mat            F,A[2];
 35:   NEPType        type;
 36:   PetscInt       n=100,nev,Istart,Iend,i;
 38:   PetscBool      terse,split=PETSC_TRUE;
 39:   RG             rg;
 40:   FN             f[2];
 41:   PetscScalar    coeffs;

 43:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 44:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 45:   PetscOptionsGetBool(NULL,NULL,"-split",&split,NULL);
 46:   PetscPrintf(PETSC_COMM_WORLD,"\nSquare root eigenproblem, n=%D%s\n\n",n,split?" (in split form)":"");

 48:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 49:      Create nonlinear eigensolver context
 50:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 52:   NEPCreate(PETSC_COMM_WORLD,&nep);

 54:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 55:      Select the NLEIGS solver and set required options for it
 56:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 58:   NEPSetType(nep,NEPNLEIGS);
 59:   NEPNLEIGSSetSingularitiesFunction(nep,ComputeSingularities,NULL);
 60:   NEPGetRG(nep,&rg);
 61:   RGSetType(rg,RGINTERVAL);
 62: #if defined(PETSC_USE_COMPLEX)
 63:   RGIntervalSetEndpoints(rg,0.01,16.0,-0.001,0.001);
 64: #else
 65:   RGIntervalSetEndpoints(rg,0.01,16.0,0,0);
 66: #endif
 67:   NEPSetTarget(nep,1.1);

 69:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 70:      Define the nonlinear problem
 71:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 73:   if (split) {
 74:     /*
 75:        Create matrices for the split form
 76:     */
 77:     MatCreate(PETSC_COMM_WORLD,&A[0]);
 78:     MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n);
 79:     MatSetFromOptions(A[0]);
 80:     MatSetUp(A[0]);
 81:     MatGetOwnershipRange(A[0],&Istart,&Iend);
 82:     for (i=Istart;i<Iend;i++) {
 83:       if (i>0) { MatSetValue(A[0],i,i-1,1.0,INSERT_VALUES); }
 84:       if (i<n-1) { MatSetValue(A[0],i,i+1,1.0,INSERT_VALUES); }
 85:       MatSetValue(A[0],i,i,-2.0,INSERT_VALUES);
 86:     }
 87:     MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY);
 88:     MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY);

 90:     MatCreate(PETSC_COMM_WORLD,&A[1]);
 91:     MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,n,n);
 92:     MatSetFromOptions(A[1]);
 93:     MatSetUp(A[1]);
 94:     MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY);
 95:     MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY);
 96:     MatShift(A[1],1.0);

 98:     /*
 99:        Define funcions for the split form
100:      */
101:     FNCreate(PETSC_COMM_WORLD,&f[0]);
102:     FNSetType(f[0],FNRATIONAL);
103:     coeffs = 1.0;
104:     FNRationalSetNumerator(f[0],1,&coeffs);
105:     FNCreate(PETSC_COMM_WORLD,&f[1]);
106:     FNSetType(f[1],FNSQRT);
107:     NEPSetSplitOperator(nep,2,A,f,SUBSET_NONZERO_PATTERN);

109:   } else {
110:     /*
111:        Callback form: create matrix and set Function evaluation routine
112:      */
113:     MatCreate(PETSC_COMM_WORLD,&F);
114:     MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n);
115:     MatSetFromOptions(F);
116:     MatSeqAIJSetPreallocation(F,3,NULL);
117:     MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
118:     MatSetUp(F);
119:     NEPSetFunction(nep,F,F,FormFunction,NULL);
120:   }

122:   /*
123:      Set solver parameters at runtime
124:   */
125:   NEPSetFromOptions(nep);

127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128:                       Solve the eigensystem
129:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130:   NEPSolve(nep);
131:   NEPGetType(nep,&type);
132:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n",type);
133:   NEPGetDimensions(nep,&nev,NULL,NULL);
134:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

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

140:   /* show detailed info unless -terse option is given by user */
141:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
142:   if (terse) {
143:     NEPErrorView(nep,NEP_ERROR_BACKWARD,NULL);
144:   } else {
145:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
146:     NEPReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
147:     NEPErrorView(nep,NEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
148:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
149:   }
150:   NEPDestroy(&nep);
151:   if (split) {
152:     MatDestroy(&A[0]);
153:     MatDestroy(&A[1]);
154:     FNDestroy(&f[0]);
155:     FNDestroy(&f[1]);
156:   } else {
157:     MatDestroy(&F);
158:   }
159:   SlepcFinalize();
160:   return ierr;
161: }

163: /* ------------------------------------------------------------------- */
164: /*
165:    FormFunction - Computes Function matrix  T(lambda)
166: */
167: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
168: {
170:   PetscInt       i,n,col[3],Istart,Iend;
171:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
172:   PetscScalar    value[3],t;

175:   /*
176:      Compute Function entries and insert into matrix
177:   */
178:   t = PetscSqrtScalar(lambda);
179:   MatGetSize(fun,&n,NULL);
180:   MatGetOwnershipRange(fun,&Istart,&Iend);
181:   if (Istart==0) FirstBlock=PETSC_TRUE;
182:   if (Iend==n) LastBlock=PETSC_TRUE;
183:   value[0]=1.0; value[1]=t-2.0; value[2]=1.0;
184:   for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
185:     col[0]=i-1; col[1]=i; col[2]=i+1;
186:     MatSetValues(fun,1,&i,3,col,value,INSERT_VALUES);
187:   }
188:   if (LastBlock) {
189:     i=n-1; col[0]=n-2; col[1]=n-1;
190:     MatSetValues(fun,1,&i,2,col,value,INSERT_VALUES);
191:   }
192:   if (FirstBlock) {
193:     i=0; col[0]=0; col[1]=1; value[0]=t-2.0; value[1]=1.0;
194:     MatSetValues(fun,1,&i,2,col,value,INSERT_VALUES);
195:   }

197:   /*
198:      Assemble matrix
199:   */
200:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
201:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
202:   if (fun != B) {
203:     MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
204:     MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
205:   }
206:   return(0);
207: }

209: /*
210:    ComputeSingularities - Computes maxnp points (at most) in the complex plane where
211:    the function T(.) is not analytic.

213:    In this case, we discretize the singularity region (-inf,0)~(-10e+6,-10e-6)
214: */
215: PetscErrorCode ComputeSingularities(NEP nep,PetscInt *maxnp,PetscScalar *xi,void *pt)
216: {
217:   PetscReal h;
218:   PetscInt  i;

221:   h = 11.0/(*maxnp-1);
222:   xi[0] = -1e-5; xi[*maxnp-1] = -1e+6;
223:   for (i=1;i<*maxnp-1;i++) xi[i] = -PetscPowReal(10,-5+h*i);
224:   return(0);
225: }