Actual source code: wiresaw.c

slepc-main 2024-11-22
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  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.

 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).

 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: */

 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";

 34: #include <slepcpep.h>

 36: int main(int argc,char **argv)
 37: {
 38:   Mat            M,D,K,A[3];      /* problem matrices */
 39:   PEP            pep;             /* polynomial eigenproblem solver context */
 40:   PetscInt       n=10,Istart,Iend,j,k;
 41:   PetscReal      v=0.01,eta=0.0;
 42:   PetscBool      terse;

 44:   PetscFunctionBeginUser;
 45:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));

 47:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 48:   PetscCall(PetscOptionsGetReal(NULL,NULL,"-v",&v,NULL));
 49:   PetscCall(PetscOptionsGetReal(NULL,NULL,"-eta",&eta,NULL));
 50:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nVibration analysis of a wiresaw, n=%" PetscInt_FMT " v=%g eta=%g\n\n",n,(double)v,(double)eta));

 52:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 53:      Compute the matrices that define the eigensystem, (k^2*M+k*D+K)x=0
 54:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 56:   /* K is a diagonal matrix */
 57:   PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
 58:   PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
 59:   PetscCall(MatSetFromOptions(K));

 61:   PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
 62:   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));

 64:   PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
 65:   PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
 66:   PetscCall(MatScale(K,0.5));

 68:   /* D is a tridiagonal */
 69:   PetscCall(MatCreate(PETSC_COMM_WORLD,&D));
 70:   PetscCall(MatSetSizes(D,PETSC_DECIDE,PETSC_DECIDE,n,n));
 71:   PetscCall(MatSetFromOptions(D));

 73:   PetscCall(MatGetOwnershipRange(D,&Istart,&Iend));
 74:   for (j=Istart;j<Iend;j++) {
 75:     for (k=0;k<n;k++) {
 76:       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:   }

 80:   PetscCall(MatAssemblyBegin(D,MAT_FINAL_ASSEMBLY));
 81:   PetscCall(MatAssemblyEnd(D,MAT_FINAL_ASSEMBLY));
 82:   PetscCall(MatScale(D,0.5));

 84:   /* M is a diagonal matrix */
 85:   PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
 86:   PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
 87:   PetscCall(MatSetFromOptions(M));
 88:   PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
 89:   for (j=Istart;j<Iend;j++) PetscCall(MatSetValue(M,j,j,1.0,INSERT_VALUES));
 90:   PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
 91:   PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
 92:   PetscCall(MatScale(M,0.5));

 94:   /* add damping */
 95:   if (eta>0.0) {
 96:     PetscCall(MatAXPY(K,eta,D,DIFFERENT_NONZERO_PATTERN)); /* K = K + eta*D */
 97:     PetscCall(MatShift(D,eta)); /* D = D + eta*eye(n) */
 98:   }

100:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101:                 Create the eigensolver and solve the problem
102:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

104:   PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
105:   A[0] = K; A[1] = D; A[2] = M;
106:   PetscCall(PEPSetOperators(pep,3,A));
107:   if (eta==0.0) PetscCall(PEPSetProblemType(pep,PEP_GYROSCOPIC));
108:   else PetscCall(PEPSetProblemType(pep,PEP_GENERAL));
109:   PetscCall(PEPSetFromOptions(pep));
110:   PetscCall(PEPSolve(pep));

112:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113:                     Display solution and clean up
114:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

116:   /* show detailed info unless -terse option is given by user */
117:   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
118:   if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
119:   else {
120:     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
121:     PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
122:     PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD));
123:     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
124:   }
125:   PetscCall(PEPDestroy(&pep));
126:   PetscCall(MatDestroy(&M));
127:   PetscCall(MatDestroy(&D));
128:   PetscCall(MatDestroy(&K));
129:   PetscCall(SlepcFinalize());
130:   return 0;
131: }

133: /*TEST

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

149: TEST*/