LCOV - code coverage report
Current view: top level - sys/classes/ds/tests - test12.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 115 119 96.6 %
Date: 2024-11-21 00:34:55 Functions: 1 1 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /*
       2             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       3             :    SLEPc - Scalable Library for Eigenvalue Problem Computations
       4             :    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
       5             : 
       6             :    This file is part of SLEPc.
       7             :    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
       8             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       9             : */
      10             : 
      11             : static char help[] = "Test DSNEP.\n\n";
      12             : 
      13             : #include <slepcds.h>
      14             : 
      15           1 : int main(int argc,char **argv)
      16             : {
      17           1 :   DS             ds;
      18           1 :   FN             f1,f2,f3,funs[3],qfun;
      19           1 :   SlepcSC        sc;
      20           1 :   PetscScalar    *Id,*A,*B,*wr,*wi,*X,*W,coeffs[2],auxr,alpha;
      21           1 :   PetscReal      tol,tau=0.001,radius=10,h,a=20,xi,re,im,nrm,aux;
      22           1 :   PetscInt       i,j,ii,jj,k,n=10,ld,nev,nfun;
      23           1 :   PetscViewer    viewer;
      24           1 :   PetscBool      verbose;
      25           1 :   RG             rg;
      26           1 :   DSMatType      mat[3]={DS_MAT_E0,DS_MAT_E1,DS_MAT_E2};
      27             : #if !defined(PETSC_USE_COMPLEX)
      28           1 :   PetscScalar    auxi;
      29             : #endif
      30             : 
      31           1 :   PetscFunctionBeginUser;
      32           1 :   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
      33           1 :   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
      34           1 :   PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
      35           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NEP - dimension %" PetscInt_FMT ", tau=%g.\n",n,(double)tau));
      36           1 :   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
      37           1 :   PetscCall(PetscOptionsGetReal(NULL,NULL,"-radius",&radius,NULL));
      38             : 
      39             :   /* Create DS object */
      40           1 :   PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
      41           1 :   PetscCall(DSSetType(ds,DSNEP));
      42           1 :   tol  = 1000*n*PETSC_MACHINE_EPSILON;
      43           1 :   PetscCall(DSNEPSetRefine(ds,tol,PETSC_DECIDE));
      44           1 :   PetscCall(DSSetFromOptions(ds));
      45             : 
      46             :   /* Set functions (prior to DSAllocate) */
      47           1 :   PetscCall(FNCreate(PETSC_COMM_WORLD,&f1));
      48           1 :   PetscCall(FNSetType(f1,FNRATIONAL));
      49           1 :   coeffs[0] = -1.0; coeffs[1] = 0.0;
      50           1 :   PetscCall(FNRationalSetNumerator(f1,2,coeffs));
      51             : 
      52           1 :   PetscCall(FNCreate(PETSC_COMM_WORLD,&f2));
      53           1 :   PetscCall(FNSetType(f2,FNRATIONAL));
      54           1 :   coeffs[0] = 1.0;
      55           1 :   PetscCall(FNRationalSetNumerator(f2,1,coeffs));
      56             : 
      57           1 :   PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
      58           1 :   PetscCall(FNSetType(f3,FNEXP));
      59           1 :   PetscCall(FNSetScale(f3,-tau,1.0));
      60             : 
      61           1 :   funs[0] = f1;
      62           1 :   funs[1] = f2;
      63           1 :   funs[2] = f3;
      64           1 :   PetscCall(DSNEPSetFN(ds,3,funs));
      65             : 
      66             :   /* Set dimensions */
      67           1 :   ld = n+2;  /* test leading dimension larger than n */
      68           1 :   PetscCall(DSAllocate(ds,ld));
      69           1 :   PetscCall(DSSetDimensions(ds,n,0,0));
      70             : 
      71             :   /* Set region (used only in method=1) */
      72           1 :   PetscCall(RGCreate(PETSC_COMM_WORLD,&rg));
      73           1 :   PetscCall(RGSetType(rg,RGELLIPSE));
      74           1 :   PetscCall(RGEllipseSetParameters(rg,0.0,radius,1.0));
      75           1 :   PetscCall(DSNEPSetRG(ds,rg));
      76           1 :   PetscCall(RGDestroy(&rg));
      77             : 
      78             :   /* Set up viewer */
      79           1 :   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
      80           1 :   PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
      81           1 :   PetscCall(DSView(ds,viewer));
      82           1 :   PetscCall(PetscViewerPopFormat(viewer));
      83           1 :   if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
      84             : 
      85             :   /* Show info about functions */
      86           1 :   PetscCall(DSNEPGetNumFN(ds,&nfun));
      87           4 :   for (i=0;i<nfun;i++) {
      88           3 :     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Function %" PetscInt_FMT ":\n",i));
      89           3 :     PetscCall(DSNEPGetFN(ds,i,&qfun));
      90           3 :     PetscCall(FNView(qfun,NULL));
      91             :   }
      92             : 
      93             :   /* Fill matrices */
      94           1 :   PetscCall(DSGetArray(ds,DS_MAT_E0,&Id));
      95          11 :   for (i=0;i<n;i++) Id[i+i*ld]=1.0;
      96           1 :   PetscCall(DSRestoreArray(ds,DS_MAT_E0,&Id));
      97           1 :   h = PETSC_PI/(PetscReal)(n+1);
      98           1 :   PetscCall(DSGetArray(ds,DS_MAT_E1,&A));
      99          11 :   for (i=0;i<n;i++) A[i+i*ld]=-2.0/(h*h)+a;
     100          10 :   for (i=1;i<n;i++) {
     101           9 :     A[i+(i-1)*ld]=1.0/(h*h);
     102           9 :     A[(i-1)+i*ld]=1.0/(h*h);
     103             :   }
     104           1 :   PetscCall(DSRestoreArray(ds,DS_MAT_E1,&A));
     105           1 :   PetscCall(DSGetArray(ds,DS_MAT_E2,&B));
     106          11 :   for (i=0;i<n;i++) {
     107          10 :     xi = (i+1)*h;
     108          10 :     B[i+i*ld] = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
     109             :   }
     110           1 :   PetscCall(DSRestoreArray(ds,DS_MAT_E2,&B));
     111             : 
     112           1 :   if (verbose) {
     113           0 :     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
     114           0 :     PetscCall(DSView(ds,viewer));
     115             :   }
     116             : 
     117             :   /* Solve */
     118           1 :   PetscCall(PetscCalloc2(n,&wr,n,&wi));
     119           1 :   PetscCall(DSGetSlepcSC(ds,&sc));
     120           1 :   sc->comparison    = SlepcCompareLargestMagnitude;
     121           1 :   sc->comparisonctx = NULL;
     122           1 :   sc->map           = NULL;
     123           1 :   sc->mapobj        = NULL;
     124           1 :   PetscCall(DSSolve(ds,wr,wi));
     125           1 :   PetscCall(DSSort(ds,wr,wi,NULL,NULL,NULL));
     126             : 
     127           1 :   if (verbose) {
     128           0 :     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
     129           0 :     PetscCall(DSView(ds,viewer));
     130             :   }
     131           1 :   PetscCall(DSGetDimensions(ds,NULL,NULL,NULL,&nev));
     132             : 
     133             :   /* Print computed eigenvalues */
     134           1 :   PetscCall(PetscMalloc1(ld*ld,&W));
     135           1 :   PetscCall(DSVectors(ds,DS_MAT_X,NULL,NULL));
     136           1 :   PetscCall(DSGetArray(ds,DS_MAT_X,&X));
     137           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n"));
     138           2 :   for (i=0;i<nev;i++) {
     139             : #if defined(PETSC_USE_COMPLEX)
     140             :     re = PetscRealPart(wr[i]);
     141             :     im = PetscImaginaryPart(wr[i]);
     142             : #else
     143           1 :     re = wr[i];
     144           1 :     im = wi[i];
     145             : #endif
     146             :     /* Residual */
     147           1 :     PetscCall(PetscArrayzero(W,ld*ld));
     148           4 :     for (k=0;k<nfun;k++) {
     149           3 :       PetscCall(FNEvaluateFunction(funs[k],wr[i],&alpha));
     150           3 :       PetscCall(DSGetArray(ds,mat[k],&A));
     151         333 :       for (jj=0;jj<n;jj++) for (ii=0;ii<n;ii++) W[jj*ld+ii] += alpha*A[jj*ld+ii];
     152           3 :       PetscCall(DSRestoreArray(ds,mat[k],&A));
     153             :     }
     154             :     nrm = 0.0;
     155          11 :     for (k=0;k<n;k++) {
     156             :       auxr = 0.0;
     157             : #if !defined(PETSC_USE_COMPLEX)
     158             :       auxi = 0.0;
     159             : #endif
     160         110 :       for (j=0;j<n;j++) {
     161         100 :         auxr += W[k+j*ld]*X[i*ld+j];
     162             : #if !defined(PETSC_USE_COMPLEX)
     163         100 :         if (PetscAbs(wi[j])!=0.0) auxi += W[k+j*ld]*X[(i+1)*ld+j];
     164             : #endif
     165             :       }
     166          10 :       aux = SlepcAbsEigenvalue(auxr,auxi);
     167          10 :       nrm += aux*aux;
     168             :     }
     169           1 :     nrm = PetscSqrtReal(nrm);
     170           1 :     if (nrm/SlepcAbsEigenvalue(wr[i],wi[i])>tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the residual norm of the %" PetscInt_FMT "-th computed eigenpair %g\n",i,(double)nrm));
     171           1 :     if (PetscAbs(im)<1e-10) PetscCall(PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)re));
     172           1 :     else PetscCall(PetscViewerASCIIPrintf(viewer,"  %.5f%+.5fi\n",(double)re,(double)im));
     173             :   }
     174           1 :   PetscCall(PetscFree(W));
     175           1 :   PetscCall(DSRestoreArray(ds,DS_MAT_X,&X));
     176           1 :   PetscCall(DSRestoreArray(ds,DS_MAT_W,&W));
     177           1 :   PetscCall(PetscFree2(wr,wi));
     178           1 :   PetscCall(FNDestroy(&f1));
     179           1 :   PetscCall(FNDestroy(&f2));
     180           1 :   PetscCall(FNDestroy(&f3));
     181           1 :   PetscCall(DSDestroy(&ds));
     182           1 :   PetscCall(SlepcFinalize());
     183             :   return 0;
     184             : }
     185             : 
     186             : /*TEST
     187             : 
     188             :    testset:
     189             :       test:
     190             :          filter: grep -v "solving the problem"
     191             :          suffix: 1
     192             :       test:
     193             :          suffix: 2
     194             :          args: -ds_method 1 -radius 10 -ds_nep_refine_its 1
     195             :          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/"
     196             :          requires: complex
     197             : 
     198             : TEST*/

Generated by: LCOV version 1.14