Actual source code: wiresaw.c

slepc-3.17.1 2022-04-11
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:   SlepcInitialize(&argc,&argv,(char*)0,help);

 46:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 47:   PetscOptionsGetReal(NULL,NULL,"-v",&v,NULL);
 48:   PetscOptionsGetReal(NULL,NULL,"-eta",&eta,NULL);
 49:   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:      Compute the matrices that define the eigensystem, (k^2*M+k*D+K)x=0
 53:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

 61:   MatGetOwnershipRange(K,&Istart,&Iend);
 62:   for (j=Istart;j<Iend;j++) MatSetValue(K,j,j,(j+1)*(j+1)*PETSC_PI*PETSC_PI*(1.0-v*v),INSERT_VALUES);

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

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

 74:   MatGetOwnershipRange(D,&Istart,&Iend);
 75:   for (j=Istart;j<Iend;j++) {
 76:     for (k=0;k<n;k++) {
 77:       if ((j+k)%2) MatSetValue(D,j,k,8.0*(j+1)*(k+1)*v/((j+1)*(j+1)-(k+1)*(k+1)),INSERT_VALUES);
 78:     }
 79:   }

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

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

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

102:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103:                 Create the eigensolver and solve the problem
104:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

114:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115:                     Display solution and clean up
116:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

135: /*TEST

137:    testset:
138:       args: -pep_nev 4 -terse
139:       requires: double
140:       output_file: output/wiresaw_1.out
141:       test:
142:          suffix: 1
143:          args: -pep_type {{toar qarnoldi}}
144:       test:
145:          suffix: 1_linear_h1
146:          args: -pep_type linear -pep_linear_explicitmatrix -pep_linear_linearization 1,0 -pep_linear_st_ksp_type bcgs -pep_linear_st_pc_type kaczmarz
147:       test:
148:          suffix: 1_linear_h2
149:          args: -pep_type linear -pep_linear_explicitmatrix -pep_linear_linearization 0,1

151: TEST*/