Actual source code: test23.c

slepc-3.16.0 2021-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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[] = "Test interface functions of DSNEP.\n\n";

 13: #include <slepcds.h>

 15: int main(int argc,char **argv)
 16: {
 18:   DS             ds;
 19:   FN             f1,f2,f3,funs[3];
 20:   SlepcSC        sc;
 21:   PetscScalar    *Id,*A,*B,*wr,*wi,*X,*W,coeffs[2],auxr,alpha;
 22:   PetscReal      tau=0.001,h,a=20,xi,re,im,nrm,aux;
 23:   PetscInt       i,j,ii,jj,k,n=10,ld,nev,nfun,midx,ip,rits,meth,spls;
 24:   PetscViewer    viewer;
 25:   PetscBool      verbose;
 26:   RG             rg;
 27:   DSMatType      mat[3]={DS_MAT_E0,DS_MAT_E1,DS_MAT_E2};
 28: #if defined(PETSC_USE_COMPLEX)
 29:   PetscBool      flg;
 30: #else
 31:   PetscScalar    auxi;
 32: #endif

 34:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 35:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 36:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 37:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NEP - dimension %D, tau=%g.\n",n,(double)tau);
 38:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);

 40:   /* Create DS object and set options */
 41:   DSCreate(PETSC_COMM_WORLD,&ds);
 42:   DSSetType(ds,DSNEP);
 43: #if defined(PETSC_USE_COMPLEX)
 44:   DSSetMethod(ds,1);  /* contour integral */
 45: #endif
 46:   DSNEPGetRG(ds,&rg);
 47:   RGSetType(rg,RGELLIPSE);
 48:   DSNEPSetMinimality(ds,1);
 49:   DSNEPSetIntegrationPoints(ds,16);
 50:   DSNEPSetRefine(ds,PETSC_DEFAULT,2);
 51:   DSNEPSetSamplingSize(ds,25);
 52:   DSSetFromOptions(ds);

 54:   /* Print current options */
 55:   DSGetMethod(ds,&meth);
 56: #if defined(PETSC_USE_COMPLEX)
 57:   if (meth!=1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"This example requires ds_method=1");
 58:   RGIsTrivial(rg,&flg);
 59:   if (flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Must at least set the radius of the ellipse");
 60: #endif

 62:   DSNEPGetMinimality(ds,&midx);
 63:   DSNEPGetIntegrationPoints(ds,&ip);
 64:   DSNEPGetRefine(ds,NULL,&rits);
 65:   DSNEPGetSamplingSize(ds,&spls);
 66:   if (meth==1) {
 67:     PetscPrintf(PETSC_COMM_WORLD,"Contour integral method with %D integration points, minimality index %D, and sampling size %D\n",ip,midx,spls);
 68:     if (rits) {
 69:       PetscPrintf(PETSC_COMM_WORLD,"Doing %D iterations of Newton refinement\n",rits);
 70:     }
 71:   }

 73:   /* Set functions (prior to DSAllocate) */
 74:   FNCreate(PETSC_COMM_WORLD,&f1);
 75:   FNSetType(f1,FNRATIONAL);
 76:   coeffs[0] = -1.0; coeffs[1] = 0.0;
 77:   FNRationalSetNumerator(f1,2,coeffs);

 79:   FNCreate(PETSC_COMM_WORLD,&f2);
 80:   FNSetType(f2,FNRATIONAL);
 81:   coeffs[0] = 1.0;
 82:   FNRationalSetNumerator(f2,1,coeffs);

 84:   FNCreate(PETSC_COMM_WORLD,&f3);
 85:   FNSetType(f3,FNEXP);
 86:   FNSetScale(f3,-tau,1.0);

 88:   funs[0] = f1;
 89:   funs[1] = f2;
 90:   funs[2] = f3;
 91:   DSNEPSetFN(ds,3,funs);

 93:   /* Set dimensions */
 94:   ld = n;
 95:   DSAllocate(ds,ld);
 96:   DSSetDimensions(ds,n,0,0);

 98:   /* Set up viewer */
 99:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
100:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
101:   PetscViewerPopFormat(viewer);
102:   if (verbose) {
103:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
104:   }

106:   /* Fill matrices */
107:   DSGetArray(ds,DS_MAT_E0,&Id);
108:   for (i=0;i<n;i++) Id[i+i*ld]=1.0;
109:   DSRestoreArray(ds,DS_MAT_E0,&Id);
110:   h = PETSC_PI/(PetscReal)(n+1);
111:   DSGetArray(ds,DS_MAT_E1,&A);
112:   for (i=0;i<n;i++) A[i+i*ld]=-2.0/(h*h)+a;
113:   for (i=1;i<n;i++) {
114:     A[i+(i-1)*ld]=1.0/(h*h);
115:     A[(i-1)+i*ld]=1.0/(h*h);
116:   }
117:   DSRestoreArray(ds,DS_MAT_E1,&A);
118:   DSGetArray(ds,DS_MAT_E2,&B);
119:   for (i=0;i<n;i++) {
120:     xi = (i+1)*h;
121:     B[i+i*ld] = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
122:   }
123:   DSRestoreArray(ds,DS_MAT_E2,&B);

125:   if (verbose) {
126:     PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
127:     DSView(ds,viewer);
128:   }

130:   /* Solve */
131:   PetscCalloc2(n,&wr,n,&wi);
132:   DSGetSlepcSC(ds,&sc);
133:   sc->comparison    = SlepcCompareLargestMagnitude;
134:   sc->comparisonctx = NULL;
135:   sc->map           = NULL;
136:   sc->mapobj        = NULL;
137:   DSSolve(ds,wr,wi);
138:   DSSort(ds,wr,wi,NULL,NULL,NULL);

140:   if (verbose) {
141:     PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
142:     DSView(ds,viewer);
143:   }
144:   DSGetDimensions(ds,NULL,NULL,NULL,&nev);

146:   /* Print computed eigenvalues */
147:   DSNEPGetNumFN(ds,&nfun);
148:   PetscMalloc1(ld*ld,&W);
149:   DSVectors(ds,DS_MAT_X,NULL,NULL);
150:   DSGetArray(ds,DS_MAT_X,&X);
151:   PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
152:   for (i=0;i<nev;i++) {
153: #if defined(PETSC_USE_COMPLEX)
154:     re = PetscRealPart(wr[i]);
155:     im = PetscImaginaryPart(wr[i]);
156: #else
157:     re = wr[i];
158:     im = wi[i];
159: #endif
160:     /* Residual */
161:     PetscArrayzero(W,ld*ld);
162:     for (k=0;k<nfun;k++) {
163:       FNEvaluateFunction(funs[k],wr[i],&alpha);
164:       DSGetArray(ds,mat[k],&A);
165:       for (jj=0;jj<n;jj++) for (ii=0;ii<n;ii++) W[jj*ld+ii] += alpha*A[jj*ld+ii];
166:       DSRestoreArray(ds,mat[k],&A);
167:     }
168:     nrm = 0.0;
169:     for (k=0;k<n;k++) {
170:       auxr = 0.0;
171: #if !defined(PETSC_USE_COMPLEX)
172:       auxi = 0.0;
173: #endif
174:       for (j=0;j<n;j++) {
175:         auxr += W[k+j*ld]*X[i*ld+j];
176: #if !defined(PETSC_USE_COMPLEX)
177:         if (PetscAbs(wi[j])!=0.0) auxi += W[k+j*ld]*X[(i+1)*ld+j];
178: #endif
179:       }
180:       aux = SlepcAbsEigenvalue(auxr,auxi);
181:       nrm += aux*aux;
182:     }
183:     nrm = PetscSqrtReal(nrm);
184:     if (nrm>1000*n*PETSC_MACHINE_EPSILON) {
185:       PetscPrintf(PETSC_COMM_WORLD,"Warning: the residual norm of the %D-th computed eigenpair %g\n",i,(double)nrm);
186:     }
187:     if (PetscAbs(im)<1e-10) {
188:       PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)re);
189:     } else {
190:       PetscViewerASCIIPrintf(viewer,"  %.5f%+.5fi\n",(double)re,(double)im);
191:     }
192:   }
193:   PetscFree(W);
194:   DSRestoreArray(ds,DS_MAT_X,&X);
195:   DSRestoreArray(ds,DS_MAT_W,&W);
196:   PetscFree2(wr,wi);
197:   FNDestroy(&f1);
198:   FNDestroy(&f2);
199:   FNDestroy(&f3);
200:   DSDestroy(&ds);
201:   SlepcFinalize();
202:   return ierr;
203: }

205: /*TEST

207:    testset:
208:       test:
209:          suffix: 1
210:          requires: !complex
211:       test:
212:          suffix: 2
213:          args: -ds_nep_rg_ellipse_radius 10
214:          filter: sed -e "s/[+-]0\.0*i//g" | sed -e "s/37411/37410/"
215:          requires: complex

217: TEST*/