Actual source code: acoustic_wave_2d.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */
 21: /*
 22:    This example implements one of the problems found at
 23:        NLEVP: A Collection of Nonlinear Eigenvalue Problems,
 24:        The University of Manchester.
 25:    The details of the collection can be found at:
 26:        [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
 27:            Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.

 29:    The acoustic_wave_2d problem is a 2-D version of acoustic_wave_1d, also
 30:    scaled for real arithmetic.
 31: */

 33: static char help[] = "NLEVP problem: acoustic_wave_2d.\n\n"
 34:   "The command line options are:\n"
 35:   "  -m <m>, where <m> = grid size, the matrices have dimension m*(m-1).\n"
 36:   "  -z <z>, where <z> = impedance (default 1.0).\n\n";

 38: #include <slepcpep.h>

 42: int main(int argc,char **argv)
 43: {
 44:   Mat            M,C,K,A[3];      /* problem matrices */
 45:   PEP            pep;             /* polynomial eigenproblem solver context */
 46:   PetscInt       m=6,n,II,Istart,Iend,i,j;
 47:   PetscScalar    z=1.0;
 48:   PetscReal      h;
 49:   char           str[50];

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

 54:   PetscOptionsGetInt(NULL,"-m",&m,NULL);
 55:   if (m<2) SETERRQ(PETSC_COMM_SELF,1,"m must be at least 2");
 56:   PetscOptionsGetScalar(NULL,"-z",&z,NULL);
 57:   h = 1.0/m;
 58:   n = m*(m-1);
 59:   SlepcSNPrintfScalar(str,50,z,PETSC_FALSE);
 60:   PetscPrintf(PETSC_COMM_WORLD,"\nAcoustic wave 2-D, n=%D (m=%D), z=%s\n\n",n,m,str);

 62:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 63:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 64:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 66:   /* K has a pattern similar to the 2D Laplacian */
 67:   MatCreate(PETSC_COMM_WORLD,&K);
 68:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
 69:   MatSetFromOptions(K);
 70:   MatSetUp(K);
 71:   
 72:   MatGetOwnershipRange(K,&Istart,&Iend);
 73:   for (II=Istart;II<Iend;II++) {
 74:     i = II/m; j = II-i*m;
 75:     if (i>0) { MatSetValue(K,II,II-m,(j==m-1)?-0.5:-1.0,INSERT_VALUES); }
 76:     if (i<m-2) { MatSetValue(K,II,II+m,(j==m-1)?-0.5:-1.0,INSERT_VALUES); }
 77:     if (j>0) { MatSetValue(K,II,II-1,-1.0,INSERT_VALUES); }
 78:     if (j<m-1) { MatSetValue(K,II,II+1,-1.0,INSERT_VALUES); }
 79:     MatSetValue(K,II,II,(j==m-1)?2.0:4.0,INSERT_VALUES);
 80:   }

 82:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 83:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 85:   /* C is the zero matrix except for a few nonzero elements on the diagonal */
 86:   MatCreate(PETSC_COMM_WORLD,&C);
 87:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 88:   MatSetFromOptions(C);
 89:   MatSetUp(C);

 91:   MatGetOwnershipRange(C,&Istart,&Iend);
 92:   for (i=Istart;i<Iend;i++) {
 93:     if (i%m==m-1) {
 94:       MatSetValue(C,i,i,-2*PETSC_PI*h/z,INSERT_VALUES);
 95:     }
 96:   }
 97:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 98:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 99:   
100:   /* M is a diagonal matrix */
101:   MatCreate(PETSC_COMM_WORLD,&M);
102:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
103:   MatSetFromOptions(M);
104:   MatSetUp(M);

106:   MatGetOwnershipRange(M,&Istart,&Iend);
107:   for (i=Istart;i<Iend;i++) {
108:     if (i%m==m-1) {
109:       MatSetValue(M,i,i,2*PETSC_PI*PETSC_PI*h*h,INSERT_VALUES);
110:     } else {
111:       MatSetValue(M,i,i,4*PETSC_PI*PETSC_PI*h*h,INSERT_VALUES);
112:     }
113:   }
114:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
115:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
116:   
117:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
118:                 Create the eigensolver and solve the problem
119:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

121:   PEPCreate(PETSC_COMM_WORLD,&pep);
122:   A[0] = K; A[1] = C; A[2] = M;
123:   PEPSetOperators(pep,3,A);
124:   PEPSetFromOptions(pep);
125:   PEPSolve(pep);

127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128:                     Display solution and clean up
129:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130:   
131:   PEPPrintSolution(pep,NULL);
132:   PEPDestroy(&pep);
133:   MatDestroy(&M);
134:   MatDestroy(&C);
135:   MatDestroy(&K);
136:   SlepcFinalize();
137:   return 0;
138: }