LCOV - code coverage report
Current view: top level - sys/classes/ds/tests - test27.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 69 73 94.5 %
Date: 2024-11-23 00:39:48 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 DSHSVD with dense storage.\n\n";
      12             : 
      13             : #include <slepcds.h>
      14             : 
      15           1 : int main(int argc,char **argv)
      16             : {
      17           1 :   DS             ds;
      18           1 :   SlepcSC        sc;
      19           1 :   PetscReal      *D,sigma,rnorm,aux;
      20           1 :   PetscScalar    *A,*U,*w;
      21           1 :   PetscInt       i,j,k,n=15,m=10,m1,p=1,ld;
      22           1 :   PetscViewer    viewer;
      23           1 :   PetscBool      verbose;
      24             : 
      25           1 :   PetscFunctionBeginUser;
      26           1 :   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
      27           1 :   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
      28           1 :   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
      29           1 :   k = PetscMin(n,m);
      30           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type HSVD - dimension %" PetscInt_FMT "x%" PetscInt_FMT ".\n",n,m));
      31           1 :   PetscCall(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL));
      32           1 :   PetscCheck(p>=0 && p<=m,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Wrong value %" PetscInt_FMT ", must be 0<=p<=%" PetscInt_FMT,p,m);
      33           1 :   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
      34             : 
      35             :   /* Create DS object */
      36           1 :   PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
      37           1 :   PetscCall(DSSetType(ds,DSHSVD));
      38           1 :   PetscCall(DSSetFromOptions(ds));
      39           1 :   ld = PetscMax(n,m)+2;  /* test leading dimension larger than n */
      40           1 :   PetscCall(DSAllocate(ds,ld));
      41           1 :   PetscCall(DSSetDimensions(ds,n,0,0));
      42           1 :   PetscCall(DSHSVDSetDimensions(ds,m));
      43           1 :   PetscCall(DSHSVDGetDimensions(ds,&m1));
      44           1 :   PetscCheck(m1==m,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Inconsistent dimension value");
      45             : 
      46             :   /* Set up viewer */
      47           1 :   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
      48           1 :   PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
      49           1 :   PetscCall(DSView(ds,viewer));
      50           1 :   PetscCall(PetscViewerPopFormat(viewer));
      51             : 
      52             :   /* Fill with a rectangular Toeplitz matrix */
      53           1 :   PetscCall(DSGetArray(ds,DS_MAT_A,&A));
      54          11 :   for (i=0;i<k;i++) A[i+i*ld]=1.0;
      55           3 :   for (j=1;j<3;j++) {
      56          19 :     for (i=0;i<m-j;i++) { if ((i+j)<m) A[i+(i+j)*ld]=(PetscScalar)(j+1); }
      57             :   }
      58           7 :   for (j=1;j<n/2;j++) {
      59          75 :     for (i=0;i<n-j;i++) { if ((i+j)<n && i<m) A[(i+j)+i*ld]=-1.0; }
      60             :   }
      61           1 :   PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
      62             :   /* Fill signature matrix */
      63           1 :   PetscCall(DSGetArrayReal(ds,DS_MAT_D,&D));
      64          16 :   for (i=0;i<n;i++) D[i] = (i<p)? -1.0: 1.0;
      65           1 :   PetscCall(DSRestoreArrayReal(ds,DS_MAT_D,&D));
      66           1 :   PetscCall(DSSetState(ds,DS_STATE_RAW));
      67           1 :   if (verbose) {
      68           0 :     PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
      69           0 :     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
      70             :   }
      71           1 :   PetscCall(DSView(ds,viewer));
      72             : 
      73             :   /* Solve */
      74           1 :   PetscCall(PetscMalloc1(k,&w));
      75           1 :   PetscCall(DSGetSlepcSC(ds,&sc));
      76           1 :   sc->comparison    = SlepcCompareLargestReal;
      77           1 :   sc->comparisonctx = NULL;
      78           1 :   sc->map           = NULL;
      79           1 :   sc->mapobj        = NULL;
      80           1 :   PetscCall(DSSolve(ds,w,NULL));
      81           1 :   PetscCall(DSSort(ds,w,NULL,NULL,NULL,NULL));
      82           1 :   if (verbose) {
      83           0 :     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
      84           0 :     PetscCall(DSView(ds,viewer));
      85             :   }
      86             : 
      87             :   /* Print singular values */
      88           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n"));
      89          11 :   for (i=0;i<k;i++) {
      90          10 :     sigma = PetscRealPart(w[i]);
      91          10 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)sigma));
      92             :   }
      93             : 
      94             :   /* Singular vectors */
      95           1 :   PetscCall(DSVectors(ds,DS_MAT_U,NULL,NULL));  /* all singular vectors */
      96           1 :   j = 0;
      97           1 :   rnorm = 0.0;
      98           1 :   PetscCall(DSGetArray(ds,DS_MAT_U,&U));
      99          16 :   for (i=0;i<n;i++) {
     100          15 :     aux = PetscAbsScalar(U[i+j*ld]);
     101          15 :     rnorm += aux*aux;
     102             :   }
     103           1 :   PetscCall(DSRestoreArray(ds,DS_MAT_U,&U));
     104           1 :   rnorm = PetscSqrtReal(rnorm);
     105           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st U vector = %.3f\n",(double)rnorm));
     106             : 
     107           1 :   PetscCall(PetscFree(w));
     108           1 :   PetscCall(DSDestroy(&ds));
     109           1 :   PetscCall(SlepcFinalize());
     110             :   return 0;
     111             : }
     112             : 
     113             : /*TEST
     114             : 
     115             :    test:
     116             :       suffix: 1
     117             :       requires: !single
     118             : 
     119             : TEST*/

Generated by: LCOV version 1.14