LCOV - code coverage report
Current view: top level - eps/tests - test22.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 122 142 85.9 %
Date: 2024-12-18 00:51:33 Functions: 2 3 66.7 %
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[] = "Illustrates how to obtain invariant subspaces. "
      12             :   "Based on ex9.\n"
      13             :   "The command line options are:\n"
      14             :   "  -n <n>, where <n> = block dimension of the 2x2 block matrix.\n"
      15             :   "  -L <L>, where <L> = bifurcation parameter.\n"
      16             :   "  -alpha <alpha>, -beta <beta>, -delta1 <delta1>,  -delta2 <delta2>,\n"
      17             :   "       where <alpha> <beta> <delta1> <delta2> = model parameters.\n\n";
      18             : 
      19             : #include <slepceps.h>
      20             : 
      21             : /*
      22             :    This example computes the eigenvalues with largest real part of the
      23             :    following matrix
      24             : 
      25             :         A = [ tau1*T+(beta-1)*I     alpha^2*I
      26             :                   -beta*I        tau2*T-alpha^2*I ],
      27             : 
      28             :    where
      29             : 
      30             :         T = tridiag{1,-2,1}
      31             :         h = 1/(n+1)
      32             :         tau1 = delta1/(h*L)^2
      33             :         tau2 = delta2/(h*L)^2
      34             :  */
      35             : 
      36             : /* Matrix operations */
      37             : PetscErrorCode MatMult_Brussel(Mat,Vec,Vec);
      38             : PetscErrorCode MatGetDiagonal_Brussel(Mat,Vec);
      39             : 
      40             : typedef struct {
      41             :   Mat         T;
      42             :   Vec         x1,x2,y1,y2;
      43             :   PetscScalar alpha,beta,tau1,tau2,sigma;
      44             : } CTX_BRUSSEL;
      45             : 
      46           8 : int main(int argc,char **argv)
      47             : {
      48           8 :   EPS            eps;
      49           8 :   Mat            A;
      50           8 :   Vec            *Q,v;
      51           8 :   PetscScalar    delta1,delta2,L,h,kr,ki;
      52           8 :   PetscReal      errest,tol,re,im,lev;
      53           8 :   PetscInt       N=30,n,i,j,Istart,Iend,nev,nconv;
      54           8 :   CTX_BRUSSEL    *ctx;
      55           8 :   PetscBool      errok,trueres;
      56             : 
      57           8 :   PetscFunctionBeginUser;
      58           8 :   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
      59           8 :   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
      60           8 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%" PetscInt_FMT "\n\n",N));
      61             : 
      62             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      63             :         Generate the matrix
      64             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      65             : 
      66           8 :   PetscCall(PetscNew(&ctx));
      67           8 :   ctx->alpha = 2.0;
      68           8 :   ctx->beta  = 5.45;
      69           8 :   delta1     = 0.008;
      70           8 :   delta2     = 0.004;
      71           8 :   L          = 0.51302;
      72             : 
      73           8 :   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-L",&L,NULL));
      74           8 :   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-alpha",&ctx->alpha,NULL));
      75           8 :   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-beta",&ctx->beta,NULL));
      76           8 :   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta1",&delta1,NULL));
      77           8 :   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta2",&delta2,NULL));
      78             : 
      79             :   /* Create matrix T */
      80           8 :   PetscCall(MatCreate(PETSC_COMM_WORLD,&ctx->T));
      81           8 :   PetscCall(MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N));
      82           8 :   PetscCall(MatSetFromOptions(ctx->T));
      83           8 :   PetscCall(MatGetOwnershipRange(ctx->T,&Istart,&Iend));
      84         348 :   for (i=Istart;i<Iend;i++) {
      85         340 :     if (i>0) PetscCall(MatSetValue(ctx->T,i,i-1,1.0,INSERT_VALUES));
      86         340 :     if (i<N-1) PetscCall(MatSetValue(ctx->T,i,i+1,1.0,INSERT_VALUES));
      87         340 :     PetscCall(MatSetValue(ctx->T,i,i,-2.0,INSERT_VALUES));
      88             :   }
      89           8 :   PetscCall(MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY));
      90           8 :   PetscCall(MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY));
      91           8 :   PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
      92             : 
      93             :   /* Fill the remaining information in the shell matrix context
      94             :      and create auxiliary vectors */
      95           8 :   h = 1.0 / (PetscReal)(N+1);
      96           8 :   ctx->tau1 = delta1 / ((h*L)*(h*L));
      97           8 :   ctx->tau2 = delta2 / ((h*L)*(h*L));
      98           8 :   ctx->sigma = 0.0;
      99           8 :   PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x1));
     100           8 :   PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x2));
     101           8 :   PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y1));
     102           8 :   PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y2));
     103             : 
     104             :   /* Create the shell matrix */
     105           8 :   PetscCall(MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A));
     106           8 :   PetscCall(MatShellSetOperation(A,MATOP_MULT,(void(*)(void))MatMult_Brussel));
     107           8 :   PetscCall(MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Brussel));
     108             : 
     109             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
     110             :                 Create the eigensolver and solve the problem
     111             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     112             : 
     113           8 :   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
     114           8 :   PetscCall(EPSSetOperators(eps,A,NULL));
     115           8 :   PetscCall(EPSSetProblemType(eps,EPS_NHEP));
     116           8 :   PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));
     117           8 :   PetscCall(EPSSetTrueResidual(eps,PETSC_FALSE));
     118           8 :   PetscCall(EPSSetFromOptions(eps));
     119           8 :   PetscCall(EPSSolve(eps));
     120             : 
     121           8 :   PetscCall(EPSGetTrueResidual(eps,&trueres));
     122             :   /*if (trueres) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Computing true residuals explicitly\n\n"));*/
     123             : 
     124             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
     125             :                     Display solution and clean up
     126             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     127             : 
     128           8 :   PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
     129           8 :   PetscCall(EPSGetTolerances(eps,&tol,NULL));
     130           8 :   PetscCall(EPSGetConverged(eps,&nconv));
     131           8 :   if (nconv<nev) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Problem: less than %" PetscInt_FMT " eigenvalues converged\n\n",nev));
     132             :   else {
     133             :     /* Check that all converged eigenpairs satisfy the requested tolerance
     134             :        (in this example we use the solver's error estimate instead of computing
     135             :        the residual norm explicitly) */
     136             :     errok = PETSC_TRUE;
     137          40 :     for (i=0;i<nev;i++) {
     138          32 :       PetscCall(EPSGetErrorEstimate(eps,i,&errest));
     139          32 :       PetscCall(EPSGetEigenpair(eps,i,&kr,&ki,NULL,NULL));
     140          32 :       errok = (errok && errest<5.0*SlepcAbsEigenvalue(kr,ki)*tol)? PETSC_TRUE: PETSC_FALSE;
     141             :     }
     142           8 :     if (!errok) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Problem: some of the first %" PetscInt_FMT " relative errors are higher than the tolerance\n\n",nev));
     143             :     else {
     144           8 :       PetscCall(PetscPrintf(PETSC_COMM_WORLD," All requested eigenvalues computed up to the required tolerance:"));
     145          16 :       for (i=0;i<=(nev-1)/8;i++) {
     146           8 :         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n     "));
     147          40 :         for (j=0;j<PetscMin(8,nev-8*i);j++) {
     148          32 :           PetscCall(EPSGetEigenpair(eps,8*i+j,&kr,&ki,NULL,NULL));
     149             : #if defined(PETSC_USE_COMPLEX)
     150          32 :           re = PetscRealPart(kr);
     151          32 :           im = PetscImaginaryPart(kr);
     152             : #else
     153             :           re = kr;
     154             :           im = ki;
     155             : #endif
     156          32 :           if (im!=0.0 && PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0;
     157          32 :           if (re!=0.0 && PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0;
     158          32 :           if (im!=0.0) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%.5f%+.5fi",(double)re,(double)im));
     159          20 :           else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%.5f",(double)re));
     160          32 :           if (8*i+j+1<nev) PetscCall(PetscPrintf(PETSC_COMM_WORLD,", "));
     161             :         }
     162             :       }
     163           8 :       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n\n"));
     164             :     }
     165             :   }
     166             : 
     167             :   /* Get an orthogonal basis of the invariant subspace and check it is indeed
     168             :      orthogonal (note that eigenvectors are not orthogonal in this case) */
     169           8 :   if (nconv>1) {
     170           8 :     PetscCall(MatCreateVecs(A,&v,NULL));
     171           8 :     PetscCall(VecDuplicateVecs(v,nconv,&Q));
     172           8 :     PetscCall(EPSGetInvariantSubspace(eps,Q));
     173           8 :     PetscCall(VecCheckOrthonormality(Q,nconv,NULL,nconv,NULL,NULL,&lev));
     174           8 :     if (lev<10*tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n"));
     175           0 :     else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)lev));
     176           8 :     PetscCall(VecDestroyVecs(nconv,&Q));
     177           8 :     PetscCall(VecDestroy(&v));
     178             :   }
     179             : 
     180           8 :   PetscCall(EPSDestroy(&eps));
     181           8 :   PetscCall(MatDestroy(&A));
     182           8 :   PetscCall(MatDestroy(&ctx->T));
     183           8 :   PetscCall(VecDestroy(&ctx->x1));
     184           8 :   PetscCall(VecDestroy(&ctx->x2));
     185           8 :   PetscCall(VecDestroy(&ctx->y1));
     186           8 :   PetscCall(VecDestroy(&ctx->y2));
     187           8 :   PetscCall(PetscFree(ctx));
     188           8 :   PetscCall(SlepcFinalize());
     189             :   return 0;
     190             : }
     191             : 
     192        6309 : PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
     193             : {
     194        6309 :   PetscInt          n;
     195        6309 :   const PetscScalar *px;
     196        6309 :   PetscScalar       *py;
     197        6309 :   CTX_BRUSSEL       *ctx;
     198             : 
     199        6309 :   PetscFunctionBeginUser;
     200        6309 :   PetscCall(MatShellGetContext(A,&ctx));
     201        6309 :   PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
     202        6309 :   PetscCall(VecGetArrayRead(x,&px));
     203        6309 :   PetscCall(VecGetArray(y,&py));
     204        6309 :   PetscCall(VecPlaceArray(ctx->x1,px));
     205        6309 :   PetscCall(VecPlaceArray(ctx->x2,px+n));
     206        6309 :   PetscCall(VecPlaceArray(ctx->y1,py));
     207        6309 :   PetscCall(VecPlaceArray(ctx->y2,py+n));
     208             : 
     209        6309 :   PetscCall(MatMult(ctx->T,ctx->x1,ctx->y1));
     210        6309 :   PetscCall(VecScale(ctx->y1,ctx->tau1));
     211        6309 :   PetscCall(VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1));
     212        6309 :   PetscCall(VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2));
     213             : 
     214        6309 :   PetscCall(MatMult(ctx->T,ctx->x2,ctx->y2));
     215        6309 :   PetscCall(VecScale(ctx->y2,ctx->tau2));
     216        6309 :   PetscCall(VecAXPY(ctx->y2,-ctx->beta,ctx->x1));
     217        6309 :   PetscCall(VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2));
     218             : 
     219        6309 :   PetscCall(VecRestoreArrayRead(x,&px));
     220        6309 :   PetscCall(VecRestoreArray(y,&py));
     221        6309 :   PetscCall(VecResetArray(ctx->x1));
     222        6309 :   PetscCall(VecResetArray(ctx->x2));
     223        6309 :   PetscCall(VecResetArray(ctx->y1));
     224        6309 :   PetscCall(VecResetArray(ctx->y2));
     225        6309 :   PetscFunctionReturn(PETSC_SUCCESS);
     226             : }
     227             : 
     228           0 : PetscErrorCode MatGetDiagonal_Brussel(Mat A,Vec diag)
     229             : {
     230           0 :   Vec            d1,d2;
     231           0 :   PetscInt       n;
     232           0 :   PetscScalar    *pd;
     233           0 :   MPI_Comm       comm;
     234           0 :   CTX_BRUSSEL    *ctx;
     235             : 
     236           0 :   PetscFunctionBeginUser;
     237           0 :   PetscCall(MatShellGetContext(A,&ctx));
     238           0 :   PetscCall(PetscObjectGetComm((PetscObject)A,&comm));
     239           0 :   PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
     240           0 :   PetscCall(VecGetArray(diag,&pd));
     241           0 :   PetscCall(VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd,&d1));
     242           0 :   PetscCall(VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd+n,&d2));
     243             : 
     244           0 :   PetscCall(VecSet(d1,-2.0*ctx->tau1 + ctx->beta - 1.0 + ctx->sigma));
     245           0 :   PetscCall(VecSet(d2,-2.0*ctx->tau2 - ctx->alpha*ctx->alpha + ctx->sigma));
     246             : 
     247           0 :   PetscCall(VecDestroy(&d1));
     248           0 :   PetscCall(VecDestroy(&d2));
     249           0 :   PetscCall(VecRestoreArray(diag,&pd));
     250           0 :   PetscFunctionReturn(PETSC_SUCCESS);
     251             : }
     252             : 
     253             : /*TEST
     254             : 
     255             :    test:
     256             :       suffix: 1
     257             :       args: -eps_nev 4 -eps_true_residual {{0 1}}
     258             :       requires: !single
     259             : 
     260             :    test:
     261             :       suffix: 2
     262             :       args: -eps_nev 4 -eps_true_residual -eps_balance oneside -eps_tol 1e-7
     263             :       requires: !single
     264             : 
     265             :    test:
     266             :       suffix: 3
     267             :       args: -n 50 -eps_nev 4 -eps_ncv 16 -eps_type subspace -eps_largest_magnitude -bv_orthog_block {{gs tsqr chol tsqrchol svqb}}
     268             :       requires: !single
     269             : 
     270             : TEST*/

Generated by: LCOV version 1.14