Actual source code: loaded_string.c

slepc-main 2024-12-17
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 one 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:    The loaded_string problem is a rational eigenvalue problem for the
 19:    finite element model of a loaded vibrating string.
 20:    This example solves the loaded_string problem by first transforming
 21:    it to a quadratic eigenvalue problem.
 22: */

 24: static char help[] = "Finite element model of a loaded vibrating string.\n\n"
 25:   "The command line options are:\n"
 26:   "  -n <n>, dimension of the matrices.\n"
 27:   "  -kappa <kappa>, stiffness of elastic spring.\n"
 28:   "  -mass <m>, mass of the attached load.\n\n";

 30: #include <slepcpep.h>

 32: #define NMAT 3

 34: int main(int argc,char **argv)
 35: {
 36:   Mat            A[3],M;      /* problem matrices */
 37:   PEP            pep;         /* polynomial eigenproblem solver context */
 38:   PetscInt       n=100,Istart,Iend,i;
 39:   PetscBool      terse;
 40:   PetscReal      kappa=1.0,m=1.0;
 41:   PetscScalar    sigma;

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

 46:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 47:   PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL));
 48:   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mass",&m,NULL));
 49:   sigma = kappa/m;
 50:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Loaded vibrating string (QEP), n=%" PetscInt_FMT " kappa=%g m=%g\n\n",n,(double)kappa,(double)m));

 52:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 53:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 54:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 55:   /* initialize matrices */
 56:   for (i=0;i<NMAT;i++) {
 57:     PetscCall(MatCreate(PETSC_COMM_WORLD,&A[i]));
 58:     PetscCall(MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n));
 59:     PetscCall(MatSetFromOptions(A[i]));
 60:   }
 61:   PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));

 63:   /* A0 */
 64:   for (i=Istart;i<Iend;i++) {
 65:     PetscCall(MatSetValue(A[0],i,i,(i==n-1)?1.0*n:2.0*n,INSERT_VALUES));
 66:     if (i>0) PetscCall(MatSetValue(A[0],i,i-1,-1.0*n,INSERT_VALUES));
 67:     if (i<n-1) PetscCall(MatSetValue(A[0],i,i+1,-1.0*n,INSERT_VALUES));
 68:   }

 70:   /* A1 */
 71:   for (i=Istart;i<Iend;i++) {
 72:     PetscCall(MatSetValue(A[1],i,i,(i==n-1)?2.0/(6.0*n):4.0/(6.0*n),INSERT_VALUES));
 73:     if (i>0) PetscCall(MatSetValue(A[1],i,i-1,1.0/(6.0*n),INSERT_VALUES));
 74:     if (i<n-1) PetscCall(MatSetValue(A[1],i,i+1,1.0/(6.0*n),INSERT_VALUES));
 75:   }

 77:   /* A2 */
 78:   if (Istart<=n-1 && n-1<Iend) PetscCall(MatSetValue(A[2],n-1,n-1,kappa,INSERT_VALUES));

 80:   /* assemble matrices */
 81:   for (i=0;i<NMAT;i++) PetscCall(MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY));
 82:   for (i=0;i<NMAT;i++) PetscCall(MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY));

 84:   /* build matrices for the QEP */
 85:   PetscCall(MatAXPY(A[2],1.0,A[0],DIFFERENT_NONZERO_PATTERN));
 86:   PetscCall(MatAXPY(A[2],sigma,A[1],SAME_NONZERO_PATTERN));
 87:   PetscCall(MatScale(A[2],-1.0));
 88:   PetscCall(MatScale(A[0],sigma));
 89:   M = A[1];
 90:   A[1] = A[2];
 91:   A[2] = M;

 93:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 94:                 Create the eigensolver and solve the problem
 95:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 97:   PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
 98:   PetscCall(PEPSetOperators(pep,3,A));
 99:   PetscCall(PEPSetFromOptions(pep));
100:   PetscCall(PEPSolve(pep));

102:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103:                     Display solution and clean up
104:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

106:   /* show detailed info unless -terse option is given by user */
107:   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
108:   if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
109:   else {
110:     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
111:     PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
112:     PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD));
113:     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
114:   }
115:   PetscCall(PEPDestroy(&pep));
116:   for (i=0;i<NMAT;i++) PetscCall(MatDestroy(&A[i]));
117:   PetscCall(SlepcFinalize());
118:   return 0;
119: }

121: /*TEST

123:    test:
124:       suffix: 1
125:       args: -pep_hyperbolic -pep_interval 4,900 -pep_type stoar -st_type sinvert -st_pc_type cholesky -terse
126:       requires: !single

128: TEST*/