LCOV - code coverage report
Current view: top level - pep/tutorials/nlevp - wiresaw.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 52 59 88.1 %
Date: 2024-11-23 00:34:26 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             :    This example implements two of the problems found at
      12             :        NLEVP: A Collection of Nonlinear Eigenvalue Problems,
      13             :        The University of Manchester.
      14             :    The details of the collection can be found at:
      15             :        [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
      16             :            Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.
      17             : 
      18             :    WIRESAW1 is a gyroscopic QEP from vibration analysis of a wiresaw,
      19             :    where the parameter V represents the speed of the wire. When the
      20             :    parameter eta is nonzero, then it turns into the WIRESAW2 problem
      21             :    (with added viscous damping, e.g. eta=0.8).
      22             : 
      23             :        [2] S. Wei and I. Kao, "Vibration analysis of wire and frequency
      24             :            response in the modern wiresaw manufacturing process", J. Sound
      25             :            Vib. 213(5):1383-1395, 2000.
      26             : */
      27             : 
      28             : static char help[] = "Vibration analysis of a wiresaw.\n\n"
      29             :   "The command line options are:\n"
      30             :   "  -n <n> ... dimension of the matrices (default 10).\n"
      31             :   "  -v <value> ... velocity of the wire (default 0.01).\n"
      32             :   "  -eta <value> ... viscous damping (default 0.0).\n\n";
      33             : 
      34             : #include <slepcpep.h>
      35             : 
      36           4 : int main(int argc,char **argv)
      37             : {
      38           4 :   Mat            M,D,K,A[3];      /* problem matrices */
      39           4 :   PEP            pep;             /* polynomial eigenproblem solver context */
      40           4 :   PetscInt       n=10,Istart,Iend,j,k;
      41           4 :   PetscReal      v=0.01,eta=0.0;
      42           4 :   PetscBool      terse;
      43             : 
      44           4 :   PetscFunctionBeginUser;
      45           4 :   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
      46             : 
      47           4 :   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
      48           4 :   PetscCall(PetscOptionsGetReal(NULL,NULL,"-v",&v,NULL));
      49           4 :   PetscCall(PetscOptionsGetReal(NULL,NULL,"-eta",&eta,NULL));
      50           4 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nVibration analysis of a wiresaw, n=%" PetscInt_FMT " v=%g eta=%g\n\n",n,(double)v,(double)eta));
      51             : 
      52             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      53             :      Compute the matrices that define the eigensystem, (k^2*M+k*D+K)x=0
      54             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      55             : 
      56             :   /* K is a diagonal matrix */
      57           4 :   PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
      58           4 :   PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
      59           4 :   PetscCall(MatSetFromOptions(K));
      60             : 
      61           4 :   PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
      62          44 :   for (j=Istart;j<Iend;j++) PetscCall(MatSetValue(K,j,j,(j+1)*(j+1)*PETSC_PI*PETSC_PI*(1.0-v*v),INSERT_VALUES));
      63             : 
      64           4 :   PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
      65           4 :   PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
      66           4 :   PetscCall(MatScale(K,0.5));
      67             : 
      68             :   /* D is a tridiagonal */
      69           4 :   PetscCall(MatCreate(PETSC_COMM_WORLD,&D));
      70           4 :   PetscCall(MatSetSizes(D,PETSC_DECIDE,PETSC_DECIDE,n,n));
      71           4 :   PetscCall(MatSetFromOptions(D));
      72             : 
      73           4 :   PetscCall(MatGetOwnershipRange(D,&Istart,&Iend));
      74          44 :   for (j=Istart;j<Iend;j++) {
      75         440 :     for (k=0;k<n;k++) {
      76         400 :       if ((j+k)%2) PetscCall(MatSetValue(D,j,k,8.0*(j+1)*(k+1)*v/((j+1)*(j+1)-(k+1)*(k+1)),INSERT_VALUES));
      77             :     }
      78             :   }
      79             : 
      80           4 :   PetscCall(MatAssemblyBegin(D,MAT_FINAL_ASSEMBLY));
      81           4 :   PetscCall(MatAssemblyEnd(D,MAT_FINAL_ASSEMBLY));
      82           4 :   PetscCall(MatScale(D,0.5));
      83             : 
      84             :   /* M is a diagonal matrix */
      85           4 :   PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
      86           4 :   PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
      87           4 :   PetscCall(MatSetFromOptions(M));
      88           4 :   PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
      89          44 :   for (j=Istart;j<Iend;j++) PetscCall(MatSetValue(M,j,j,1.0,INSERT_VALUES));
      90           4 :   PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
      91           4 :   PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
      92           4 :   PetscCall(MatScale(M,0.5));
      93             : 
      94             :   /* add damping */
      95           4 :   if (eta>0.0) {
      96           0 :     PetscCall(MatAXPY(K,eta,D,DIFFERENT_NONZERO_PATTERN)); /* K = K + eta*D */
      97           0 :     PetscCall(MatShift(D,eta)); /* D = D + eta*eye(n) */
      98             :   }
      99             : 
     100             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
     101             :                 Create the eigensolver and solve the problem
     102             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     103             : 
     104           4 :   PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
     105           4 :   A[0] = K; A[1] = D; A[2] = M;
     106           4 :   PetscCall(PEPSetOperators(pep,3,A));
     107           4 :   if (eta==0.0) PetscCall(PEPSetProblemType(pep,PEP_GYROSCOPIC));
     108           0 :   else PetscCall(PEPSetProblemType(pep,PEP_GENERAL));
     109           4 :   PetscCall(PEPSetFromOptions(pep));
     110           4 :   PetscCall(PEPSolve(pep));
     111             : 
     112             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
     113             :                     Display solution and clean up
     114             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     115             : 
     116             :   /* show detailed info unless -terse option is given by user */
     117           4 :   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
     118           4 :   if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
     119             :   else {
     120           0 :     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
     121           0 :     PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
     122           0 :     PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD));
     123           0 :     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
     124             :   }
     125           4 :   PetscCall(PEPDestroy(&pep));
     126           4 :   PetscCall(MatDestroy(&M));
     127           4 :   PetscCall(MatDestroy(&D));
     128           4 :   PetscCall(MatDestroy(&K));
     129           4 :   PetscCall(SlepcFinalize());
     130             :   return 0;
     131             : }
     132             : 
     133             : /*TEST
     134             : 
     135             :    testset:
     136             :       args: -pep_nev 4 -terse
     137             :       requires: double
     138             :       output_file: output/wiresaw_1.out
     139             :       test:
     140             :          suffix: 1
     141             :          args: -pep_type {{toar qarnoldi}}
     142             :       test:
     143             :          suffix: 1_linear_h1
     144             :          args: -pep_type linear -pep_linear_explicitmatrix -pep_linear_linearization 1,0 -pep_linear_st_ksp_type bcgs -pep_linear_st_pc_type kaczmarz
     145             :       test:
     146             :          suffix: 1_linear_h2
     147             :          args: -pep_type linear -pep_linear_explicitmatrix -pep_linear_linearization 0,1
     148             : 
     149             : TEST*/

Generated by: LCOV version 1.14