Actual source code: acoustic_wave_2d.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_2d problem is a 2-D version of acoustic_wave_1d, also
 19:    scaled for real arithmetic.
 20: */

 22: static char help[] = "Quadratic eigenproblem from an acoustics application (2-D).\n\n"
 23:   "The command line options are:\n"
 24:   "  -m <m>, where <m> = grid size, the matrices have dimension m*(m-1).\n"
 25:   "  -z <z>, where <z> = impedance (default 1.0).\n\n";

 27: #include <slepcpep.h>

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

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

 42:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 44:   PetscOptionsGetScalar(NULL,NULL,"-z",&z,NULL);
 45:   h = 1.0/m;
 46:   n = m*(m-1);
 47:   SlepcSNPrintfScalar(str,sizeof(str),z,PETSC_FALSE);
 48:   PetscPrintf(PETSC_COMM_WORLD,"\nAcoustic wave 2-D, n=%" PetscInt_FMT " (m=%" PetscInt_FMT "), z=%s\n\n",n,m,str);

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

 54:   /* K has a pattern similar to the 2D Laplacian */
 55:   MatCreate(PETSC_COMM_WORLD,&K);
 56:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
 57:   MatSetFromOptions(K);
 58:   MatSetUp(K);

 60:   MatGetOwnershipRange(K,&Istart,&Iend);
 61:   for (II=Istart;II<Iend;II++) {
 62:     i = II/m; j = II-i*m;
 63:     if (i>0) MatSetValue(K,II,II-m,(j==m-1)?-0.5:-1.0,INSERT_VALUES);
 64:     if (i<m-2) MatSetValue(K,II,II+m,(j==m-1)?-0.5:-1.0,INSERT_VALUES);
 65:     if (j>0) MatSetValue(K,II,II-1,-1.0,INSERT_VALUES);
 66:     if (j<m-1) MatSetValue(K,II,II+1,-1.0,INSERT_VALUES);
 67:     MatSetValue(K,II,II,(j==m-1)?2.0:4.0,INSERT_VALUES);
 68:   }

 70:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 71:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 73:   /* C is the zero matrix except for a few nonzero elements on the diagonal */
 74:   MatCreate(PETSC_COMM_WORLD,&C);
 75:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 76:   MatSetFromOptions(C);
 77:   MatSetUp(C);

 79:   MatGetOwnershipRange(C,&Istart,&Iend);
 80:   for (i=Istart;i<Iend;i++) {
 81:     if (i%m==m-1) MatSetValue(C,i,i,-2*PETSC_PI*h/z,INSERT_VALUES);
 82:   }
 83:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 84:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 86:   /* M is a diagonal matrix */
 87:   MatCreate(PETSC_COMM_WORLD,&M);
 88:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
 89:   MatSetFromOptions(M);
 90:   MatSetUp(M);

 92:   MatGetOwnershipRange(M,&Istart,&Iend);
 93:   for (i=Istart;i<Iend;i++) {
 94:     if (i%m==m-1) MatSetValue(M,i,i,2*PETSC_PI*PETSC_PI*h*h,INSERT_VALUES);
 95:     else MatSetValue(M,i,i,4*PETSC_PI*PETSC_PI*h*h,INSERT_VALUES);
 96:   }
 97:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
 98:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

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

104:   PEPCreate(PETSC_COMM_WORLD,&pep);
105:   A[0] = K; A[1] = C; A[2] = M;
106:   PEPSetOperators(pep,3,A);
107:   PEPSetFromOptions(pep);
108:   PEPSolve(pep);

110:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111:                     Display solution and clean up
112:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

131: /*TEST

133:    testset:
134:       args: -pep_nev 2 -pep_ncv 18 -terse
135:       output_file: output/acoustic_wave_2d_1.out
136:       filter: sed -e "s/2.60936i/2.60937i/g" | sed -e "s/2.60938i/2.60937i/g"
137:       test:
138:          suffix: 1
139:          args: -pep_type {{qarnoldi linear}}
140:       test:
141:          suffix: 1_toar
142:          args: -pep_type toar -pep_toar_locking 0

144:    testset:
145:       args: -pep_nev 2 -pep_ncv 18 -pep_type stoar -pep_hermitian -pep_scale scalar -st_type sinvert -terse
146:       output_file: output/acoustic_wave_2d_2.out
147:       test:
148:          suffix: 2
149:       test:
150:          suffix: 2_lin_b
151:          args: -pep_stoar_linearization 0,1
152:       test:
153:          suffix: 2_lin_ab
154:          args: -pep_stoar_linearization 0.1,0.9

156: TEST*/