LCOV - code coverage report
Current view: top level - svd/tutorials/cnetwork - embedgsvd.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 63 69 91.3 %
Date: 2024-05-07 00:31:30 Functions: 2 2 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[] = "Computes a GSVD associated with a complex network analysis problem.\n"
      12             :   "The command line options are:\n"
      13             :   "  -file <filename>, text file in TSV format containing undirected graph.\n"
      14             :   "  -beta <beta>, optional decay parameter.\n\n";
      15             : 
      16             : /*
      17             :     Computes the Katz high-order proximity embedding of a graph via a partial GSVD.
      18             : 
      19             :     [1] M. Ou et al. "Asymmetric transitivity preserving graph embedding" Proc. of
      20             :         ACM SIGKDD 2016 - https://doi.org/10.1145/2939672.2939751
      21             : */
      22             : 
      23             : #include <slepcsvd.h>
      24             : #include "network.h"
      25             : 
      26           1 : PetscErrorCode SpectralRadius(Mat A,PetscReal *rho)
      27             : {
      28           1 :   EPS         eps;
      29           1 :   PetscInt    nconv;
      30           1 :   PetscScalar kr,ki;
      31             : 
      32           1 :   PetscFunctionBeginUser;
      33           1 :   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
      34           1 :   PetscCall(EPSSetOperators(eps,A,NULL));
      35           1 :   PetscCall(EPSSetProblemType(eps,EPS_NHEP));
      36           1 :   PetscCall(EPSSetDimensions(eps,1,PETSC_DEFAULT,PETSC_DEFAULT));
      37           1 :   PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE));
      38           1 :   PetscCall(EPSSetFromOptions(eps));
      39           1 :   PetscCall(EPSSolve(eps));
      40           1 :   PetscCall(EPSGetConverged(eps,&nconv));
      41           1 :   *rho = 0.0;
      42           1 :   if (nconv) {
      43           1 :     PetscCall(EPSGetEigenvalue(eps,0,&kr,&ki));
      44           1 :     *rho = SlepcAbsEigenvalue(kr,ki);
      45             :   }
      46           1 :   PetscCall(EPSDestroy(&eps));
      47           1 :   PetscFunctionReturn(PETSC_SUCCESS);
      48             : }
      49             : 
      50           1 : int main(int argc,char **argv)
      51             : {
      52           1 :   Graph          graph;           /* directed graph */
      53           1 :   Mat            A,B;             /* matrices */
      54           1 :   SVD            svd;             /* singular value problem solver context */
      55           1 :   PetscReal      beta=0.0,rho=0.0;
      56           1 :   char           filename[PETSC_MAX_PATH_LEN];
      57           1 :   PetscBool      flg,terse;
      58             : 
      59           1 :   PetscFunctionBeginUser;
      60           1 :   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
      61             : 
      62             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      63             :         Load graph and build adjacency matrix
      64             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      65             : 
      66           1 :   PetscCall(PetscOptionsGetString(NULL,NULL,"-graph",filename,sizeof(filename),&flg));
      67           1 :   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Must use -graph <name> to indicate a file containing the graph in TSV format");
      68             : 
      69           1 :   PetscCall(GraphCreate(PETSC_COMM_WORLD,&graph));
      70           1 :   PetscCall(GraphPreload(graph,filename));
      71           1 :   PetscCheck(graph->type==GRAPH_DIRECTED,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"A directed graph is required");
      72           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nGraph with %" PetscInt_FMT " vertices and %" PetscInt_FMT " edges.\n\n",graph->nvertices,graph->nedges));
      73           1 :   if (graph->weight!=GRAPH_WEIGHT_UNWEIGHTED) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"WARNING: ignoring weights in input graph\n"));
      74             : 
      75           1 :   PetscCall(GraphPreallocate(graph,filename));
      76           1 :   PetscCall(GraphLoadUnweighted(graph,filename));
      77             : 
      78           1 :   PetscCall(MatViewFromOptions(graph->adjacency,NULL,"-adjacency"));
      79             : 
      80             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      81             :       Build high-order proximity matrices (Katz): A=I-beta*G, B=beta*G
      82             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      83             : 
      84           1 :   PetscCall(PetscOptionsGetReal(NULL,NULL,"-beta",&beta,NULL));
      85           1 :   if (beta==0.0) {  /* compute decay parameter beta */
      86           1 :     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computing spectral radius...\n"));
      87           1 :     PetscCall(SpectralRadius(graph->adjacency,&rho));
      88           1 :     if (rho) {
      89           1 :       beta = 0.8*rho;
      90           1 :       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Done, rho=%g, setting decay parameter beta=%g\n\n",(double)rho,(double)beta));
      91             :     } else {
      92           0 :       beta = 0.5;
      93           0 :       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Failed, using default decay parameter beta=%g\n\n",(double)beta));
      94             :     }
      95             :   }
      96           1 :   PetscCall(MatDuplicate(graph->adjacency,MAT_COPY_VALUES,&A));
      97           1 :   PetscCall(MatScale(A,-beta));
      98           1 :   PetscCall(MatShift(A,1.0));
      99           1 :   PetscCall(MatDuplicate(graph->adjacency,MAT_COPY_VALUES,&B));
     100           1 :   PetscCall(MatScale(B,beta));
     101             : 
     102             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
     103             :            Create the singular value solver and set various options
     104             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     105             : 
     106           1 :   PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
     107           1 :   PetscCall(SVDSetOperators(svd,A,B));
     108           1 :   PetscCall(SVDSetProblemType(svd,SVD_GENERALIZED));
     109           1 :   PetscCall(SVDSetType(svd,SVDTRLANCZOS));
     110           1 :   PetscCall(SVDSetFromOptions(svd));
     111             : 
     112             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
     113             :                       Solve the problem and print solution
     114             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     115             : 
     116           1 :   PetscCall(SVDSolve(svd));
     117             : 
     118             :   /* show detailed info unless -terse option is given by user */
     119           1 :   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
     120           1 :   if (terse) PetscCall(SVDErrorView(svd,SVD_ERROR_NORM,NULL));
     121             :   else {
     122           0 :     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
     123           0 :     PetscCall(SVDConvergedReasonView(svd,PETSC_VIEWER_STDOUT_WORLD));
     124           0 :     PetscCall(SVDErrorView(svd,SVD_ERROR_NORM,PETSC_VIEWER_STDOUT_WORLD));
     125           0 :     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
     126             :   }
     127           1 :   PetscCall(SVDDestroy(&svd));
     128           1 :   PetscCall(GraphDestroy(&graph));
     129           1 :   PetscCall(MatDestroy(&A));
     130           1 :   PetscCall(MatDestroy(&B));
     131           1 :   PetscCall(SlepcFinalize());
     132             :   return 0;
     133             : }
     134             : /*TEST
     135             : 
     136             :    build:
     137             :       depends: network.o
     138             : 
     139             :    test:
     140             :       args: -graph ${SLEPC_DIR}/share/slepc/datafiles/graphs/out.moreno_taro_taro -svd_nsv 4 -terse
     141             :       filter: sed -e 's/4.38031/4.38032/' | sed -e 's/3.75089/3.7509/' | sed -e 's/3.00071/3.00072/'
     142             : 
     143             : TEST*/

Generated by: LCOV version 1.14