Actual source code: acoustic_wave_1d.c

slepc-3.18.2 2023-01-26
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;

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

 43:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 44:   PetscOptionsGetScalar(NULL,NULL,"-z",&z,NULL);
 45:   SlepcSNPrintfScalar(str,sizeof(str),z,PETSC_FALSE);
 46:   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:   MatCreate(PETSC_COMM_WORLD,&K);
 54:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
 55:   MatSetFromOptions(K);
 56:   MatSetUp(K);

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

 67:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 68:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

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

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

 81:   /* M is a diagonal matrix */
 82:   MatCreate(PETSC_COMM_WORLD,&M);
 83:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
 84:   MatSetFromOptions(M);
 85:   MatSetUp(M);

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

 95:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 96:                 Create the eigensolver and solve the problem
 97:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 99:   PEPCreate(PETSC_COMM_WORLD,&pep);
100:   A[0] = K; A[1] = C; A[2] = M;
101:   PEPSetOperators(pep,3,A);
102:   PEPSetFromOptions(pep);
103:   PEPSolve(pep);

105:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106:                     Display solution and clean up
107:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

109:   /* show detailed info unless -terse option is given by user */
110:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
111:   if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
112:   else {
113:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
114:     PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
115:     PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
116:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
117:   }
118:   PEPDestroy(&pep);
119:   MatDestroy(&M);
120:   MatDestroy(&C);
121:   MatDestroy(&K);
122:   SlepcFinalize();
123:   return 0;
124: }

126: /*TEST

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

148: TEST*/