Actual source code: acoustic_wave_1d.c

slepc-3.21.0 2024-03-30
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 acoustic_wave_1d problem is a QEP from an acoustics application.
 19:    Here we solve it with the eigenvalue scaled by the imaginary unit, to be
 20:    able to use real arithmetic, so the computed eigenvalues should be scaled
 21:    back.
 22: */

 24: static char help[] = "Quadratic eigenproblem from an acoustics application (1-D).\n\n"
 25:   "The command line options are:\n"
 26:   "  -n <n>, where <n> = dimension of the matrices.\n"
 27:   "  -z <z>, where <z> = impedance (default 1.0).\n\n";

 29: #include <slepcpep.h>

 31: int main(int argc,char **argv)
 32: {
 33:   Mat            M,C,K,A[3];      /* problem matrices */
 34:   PEP            pep;             /* polynomial eigenproblem solver context */
 35:   PetscInt       n=10,Istart,Iend,i;
 36:   PetscScalar    z=1.0;
 37:   char           str[50];
 38:   PetscBool      terse;

 40:   PetscFunctionBeginUser;
 41:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

 43:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 44:   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-z",&z,NULL));
 45:   PetscCall(SlepcSNPrintfScalar(str,sizeof(str),z,PETSC_FALSE));
 46:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nAcoustic wave 1-D, n=%" PetscInt_FMT " z=%s\n\n",n,str));

 48:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 49:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 50:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 52:   /* K is a tridiagonal */
 53:   PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
 54:   PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
 55:   PetscCall(MatSetFromOptions(K));

 57:   PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
 58:   for (i=Istart;i<Iend;i++) {
 59:     if (i>0) PetscCall(MatSetValue(K,i,i-1,-1.0*n,INSERT_VALUES));
 60:     if (i<n-1) {
 61:       PetscCall(MatSetValue(K,i,i,2.0*n,INSERT_VALUES));
 62:       PetscCall(MatSetValue(K,i,i+1,-1.0*n,INSERT_VALUES));
 63:     } else PetscCall(MatSetValue(K,i,i,1.0*n,INSERT_VALUES));
 64:   }

 66:   PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
 67:   PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));

 69:   /* C is the zero matrix but one element*/
 70:   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
 71:   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n));
 72:   PetscCall(MatSetFromOptions(C));

 74:   PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
 75:   if (n-1>=Istart && n-1<Iend) PetscCall(MatSetValue(C,n-1,n-1,-2*PETSC_PI/z,INSERT_VALUES));
 76:   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
 77:   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));

 79:   /* M is a diagonal matrix */
 80:   PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
 81:   PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
 82:   PetscCall(MatSetFromOptions(M));

 84:   PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
 85:   for (i=Istart;i<Iend;i++) {
 86:     if (i<n-1) PetscCall(MatSetValue(M,i,i,4*PETSC_PI*PETSC_PI/n,INSERT_VALUES));
 87:     else PetscCall(MatSetValue(M,i,i,2*PETSC_PI*PETSC_PI/n,INSERT_VALUES));
 88:   }
 89:   PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
 90:   PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));

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

 96:   PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
 97:   A[0] = K; A[1] = C; A[2] = M;
 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:   PetscCall(MatDestroy(&M));
117:   PetscCall(MatDestroy(&C));
118:   PetscCall(MatDestroy(&K));
119:   PetscCall(SlepcFinalize());
120:   return 0;
121: }

123: /*TEST

125:    testset:
126:       args: -pep_nev 4 -pep_tol 1e-7 -n 24 -terse
127:       output_file: output/acoustic_wave_1d_1.out
128:       requires: !single
129:       test:
130:          suffix: 1
131:          args: -st_type sinvert -st_transform -pep_type {{toar qarnoldi linear}}
132:       test:
133:          suffix: 1_stoar
134:          args: -st_type sinvert -st_transform -pep_type stoar -pep_hermitian -pep_stoar_locking 0 -pep_stoar_nev 11 -pep_ncv 10
135:       test:
136:          suffix: 2
137:          args: -st_type sinvert -st_transform -pep_type toar -pep_extract {{none norm residual}}
138:       test:
139:          suffix: 3
140:          args: -st_type sinvert -pep_type linear -pep_extract {{none norm residual}}
141:       test:
142:          suffix: 4
143:          args: -pep_type jd

145: TEST*/