Actual source code: test12.c

slepc-3.17.2 2022-08-09
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: */

 11: static char help[] = "Test DSNEP.\n\n";

 13: #include <slepcds.h>

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

 31:   SlepcInitialize(&argc,&argv,(char*)0,help);
 32:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 33:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 34:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NEP - dimension %" PetscInt_FMT ", tau=%g.\n",n,(double)tau);
 35:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 36:   PetscOptionsGetReal(NULL,NULL,"-radius",&radius,NULL);

 38:   /* Create DS object */
 39:   DSCreate(PETSC_COMM_WORLD,&ds);
 40:   DSSetType(ds,DSNEP);
 41:   tol  = 1000*n*PETSC_MACHINE_EPSILON;
 42:   DSNEPSetRefine(ds,tol,PETSC_DECIDE);
 43:   DSSetFromOptions(ds);

 45:   /* Set functions (prior to DSAllocate) */
 46:   FNCreate(PETSC_COMM_WORLD,&f1);
 47:   FNSetType(f1,FNRATIONAL);
 48:   coeffs[0] = -1.0; coeffs[1] = 0.0;
 49:   FNRationalSetNumerator(f1,2,coeffs);

 51:   FNCreate(PETSC_COMM_WORLD,&f2);
 52:   FNSetType(f2,FNRATIONAL);
 53:   coeffs[0] = 1.0;
 54:   FNRationalSetNumerator(f2,1,coeffs);

 56:   FNCreate(PETSC_COMM_WORLD,&f3);
 57:   FNSetType(f3,FNEXP);
 58:   FNSetScale(f3,-tau,1.0);

 60:   funs[0] = f1;
 61:   funs[1] = f2;
 62:   funs[2] = f3;
 63:   DSNEPSetFN(ds,3,funs);

 65:   /* Set dimensions */
 66:   ld = n+2;  /* test leading dimension larger than n */
 67:   DSAllocate(ds,ld);
 68:   DSSetDimensions(ds,n,0,0);

 70:   /* Set region (used only in method=1) */
 71:   RGCreate(PETSC_COMM_WORLD,&rg);
 72:   RGSetType(rg,RGELLIPSE);
 73:   RGEllipseSetParameters(rg,0.0,radius,1.0);
 74:   DSNEPSetRG(ds,rg);
 75:   RGDestroy(&rg);

 77:   /* Set up viewer */
 78:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 79:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
 80:   DSView(ds,viewer);
 81:   PetscViewerPopFormat(viewer);
 82:   if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);

 84:   /* Show info about functions */
 85:   DSNEPGetNumFN(ds,&nfun);
 86:   for (i=0;i<nfun;i++) {
 87:     PetscPrintf(PETSC_COMM_WORLD,"Function %" PetscInt_FMT ":\n",i);
 88:     DSNEPGetFN(ds,i,&qfun);
 89:     FNView(qfun,NULL);
 90:   }

 92:   /* Fill matrices */
 93:   DSGetArray(ds,DS_MAT_E0,&Id);
 94:   for (i=0;i<n;i++) Id[i+i*ld]=1.0;
 95:   DSRestoreArray(ds,DS_MAT_E0,&Id);
 96:   h = PETSC_PI/(PetscReal)(n+1);
 97:   DSGetArray(ds,DS_MAT_E1,&A);
 98:   for (i=0;i<n;i++) A[i+i*ld]=-2.0/(h*h)+a;
 99:   for (i=1;i<n;i++) {
100:     A[i+(i-1)*ld]=1.0/(h*h);
101:     A[(i-1)+i*ld]=1.0/(h*h);
102:   }
103:   DSRestoreArray(ds,DS_MAT_E1,&A);
104:   DSGetArray(ds,DS_MAT_E2,&B);
105:   for (i=0;i<n;i++) {
106:     xi = (i+1)*h;
107:     B[i+i*ld] = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
108:   }
109:   DSRestoreArray(ds,DS_MAT_E2,&B);

111:   if (verbose) {
112:     PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
113:     DSView(ds,viewer);
114:   }

116:   /* Solve */
117:   PetscCalloc2(n,&wr,n,&wi);
118:   DSGetSlepcSC(ds,&sc);
119:   sc->comparison    = SlepcCompareLargestMagnitude;
120:   sc->comparisonctx = NULL;
121:   sc->map           = NULL;
122:   sc->mapobj        = NULL;
123:   DSSolve(ds,wr,wi);
124:   DSSort(ds,wr,wi,NULL,NULL,NULL);

126:   if (verbose) {
127:     PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
128:     DSView(ds,viewer);
129:   }
130:   DSGetDimensions(ds,NULL,NULL,NULL,&nev);

132:   /* Print computed eigenvalues */
133:   PetscMalloc1(ld*ld,&W);
134:   DSVectors(ds,DS_MAT_X,NULL,NULL);
135:   DSGetArray(ds,DS_MAT_X,&X);
136:   PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
137:   for (i=0;i<nev;i++) {
138: #if defined(PETSC_USE_COMPLEX)
139:     re = PetscRealPart(wr[i]);
140:     im = PetscImaginaryPart(wr[i]);
141: #else
142:     re = wr[i];
143:     im = wi[i];
144: #endif
145:     /* Residual */
146:     PetscArrayzero(W,ld*ld);
147:     for (k=0;k<nfun;k++) {
148:       FNEvaluateFunction(funs[k],wr[i],&alpha);
149:       DSGetArray(ds,mat[k],&A);
150:       for (jj=0;jj<n;jj++) for (ii=0;ii<n;ii++) W[jj*ld+ii] += alpha*A[jj*ld+ii];
151:       DSRestoreArray(ds,mat[k],&A);
152:     }
153:     nrm = 0.0;
154:     for (k=0;k<n;k++) {
155:       auxr = 0.0;
156: #if !defined(PETSC_USE_COMPLEX)
157:       auxi = 0.0;
158: #endif
159:       for (j=0;j<n;j++) {
160:         auxr += W[k+j*ld]*X[i*ld+j];
161: #if !defined(PETSC_USE_COMPLEX)
162:         if (PetscAbs(wi[j])!=0.0) auxi += W[k+j*ld]*X[(i+1)*ld+j];
163: #endif
164:       }
165:       aux = SlepcAbsEigenvalue(auxr,auxi);
166:       nrm += aux*aux;
167:     }
168:     nrm = PetscSqrtReal(nrm);
169:     if (nrm/SlepcAbsEigenvalue(wr[i],wi[i])>tol) PetscPrintf(PETSC_COMM_WORLD,"Warning: the residual norm of the %" PetscInt_FMT "-th computed eigenpair %g\n",i,(double)nrm);
170:     if (PetscAbs(im)<1e-10) PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)re);
171:     else PetscViewerASCIIPrintf(viewer,"  %.5f%+.5fi\n",(double)re,(double)im);
172:   }
173:   PetscFree(W);
174:   DSRestoreArray(ds,DS_MAT_X,&X);
175:   DSRestoreArray(ds,DS_MAT_W,&W);
176:   PetscFree2(wr,wi);
177:   FNDestroy(&f1);
178:   FNDestroy(&f2);
179:   FNDestroy(&f3);
180:   DSDestroy(&ds);
181:   SlepcFinalize();
182:   return 0;
183: }

185: /*TEST

187:    testset:
188:       test:
189:          filter: grep -v "solving the problem"
190:          suffix: 1
191:       test:
192:          suffix: 2
193:          args: -ds_method 1 -radius 10 -ds_nep_refine_its 1
194:          filter: grep -v "solving the problem" | sed -e "s/[+-]0\.0*i//g" | sed -e "s/37411/37410/" | sed -e "s/tolerance [0-9]\.[0-9]*e[+-]\([0-9]*\)/tolerance removed/" | sed -e "s/tolerance [0-9]\.\([0-9]*\)/tolerance removed/"
195:          requires: complex

197: TEST*/