Actual source code: acoustic_wave_2d.c
slepc-3.17.1 2022-04-11
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;
39: SlepcInitialize(&argc,&argv,(char*)0,help);
41: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
43: PetscOptionsGetScalar(NULL,NULL,"-z",&z,NULL);
44: h = 1.0/m;
45: n = m*(m-1);
46: SlepcSNPrintfScalar(str,sizeof(str),z,PETSC_FALSE);
47: PetscPrintf(PETSC_COMM_WORLD,"\nAcoustic wave 2-D, n=%" PetscInt_FMT " (m=%" PetscInt_FMT "), z=%s\n\n",n,m,str);
49: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
51: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53: /* K has a pattern similar to the 2D Laplacian */
54: MatCreate(PETSC_COMM_WORLD,&K);
55: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
56: MatSetFromOptions(K);
57: MatSetUp(K);
59: MatGetOwnershipRange(K,&Istart,&Iend);
60: for (II=Istart;II<Iend;II++) {
61: i = II/m; j = II-i*m;
62: if (i>0) MatSetValue(K,II,II-m,(j==m-1)?-0.5:-1.0,INSERT_VALUES);
63: if (i<m-2) MatSetValue(K,II,II+m,(j==m-1)?-0.5:-1.0,INSERT_VALUES);
64: if (j>0) MatSetValue(K,II,II-1,-1.0,INSERT_VALUES);
65: if (j<m-1) MatSetValue(K,II,II+1,-1.0,INSERT_VALUES);
66: MatSetValue(K,II,II,(j==m-1)?2.0:4.0,INSERT_VALUES);
67: }
69: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
70: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
72: /* C is the zero matrix except for a few nonzero elements on the diagonal */
73: MatCreate(PETSC_COMM_WORLD,&C);
74: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
75: MatSetFromOptions(C);
76: MatSetUp(C);
78: MatGetOwnershipRange(C,&Istart,&Iend);
79: for (i=Istart;i<Iend;i++) {
80: if (i%m==m-1) MatSetValue(C,i,i,-2*PETSC_PI*h/z,INSERT_VALUES);
81: }
82: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
83: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
85: /* M is a diagonal matrix */
86: MatCreate(PETSC_COMM_WORLD,&M);
87: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
88: MatSetFromOptions(M);
89: MatSetUp(M);
91: MatGetOwnershipRange(M,&Istart,&Iend);
92: for (i=Istart;i<Iend;i++) {
93: if (i%m==m-1) MatSetValue(M,i,i,2*PETSC_PI*PETSC_PI*h*h,INSERT_VALUES);
94: else MatSetValue(M,i,i,4*PETSC_PI*PETSC_PI*h*h,INSERT_VALUES);
95: }
96: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
97: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
99: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: Create the eigensolver and solve the problem
101: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103: PEPCreate(PETSC_COMM_WORLD,&pep);
104: A[0] = K; A[1] = C; A[2] = M;
105: PEPSetOperators(pep,3,A);
106: PEPSetFromOptions(pep);
107: PEPSolve(pep);
109: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: Display solution and clean up
111: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113: /* show detailed info unless -terse option is given by user */
114: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
115: if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
116: else {
117: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
118: PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
119: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
120: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
121: }
122: PEPDestroy(&pep);
123: MatDestroy(&M);
124: MatDestroy(&C);
125: MatDestroy(&K);
126: SlepcFinalize();
127: return 0;
128: }
130: /*TEST
132: testset:
133: args: -pep_nev 2 -pep_ncv 18 -terse
134: output_file: output/acoustic_wave_2d_1.out
135: filter: sed -e "s/2.60936i/2.60937i/g" | sed -e "s/2.60938i/2.60937i/g"
136: test:
137: suffix: 1
138: args: -pep_type {{qarnoldi linear}}
139: test:
140: suffix: 1_toar
141: args: -pep_type toar -pep_toar_locking 0
143: testset:
144: args: -pep_nev 2 -pep_ncv 18 -pep_type stoar -pep_hermitian -pep_scale scalar -st_type sinvert -terse
145: output_file: output/acoustic_wave_2d_2.out
146: test:
147: suffix: 2
148: test:
149: suffix: 2_lin_b
150: args: -pep_stoar_linearization 0,1
151: test:
152: suffix: 2_lin_ab
153: args: -pep_stoar_linearization 0.1,0.9
155: TEST*/