Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
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.
17 :
18 : The acoustic_wave_2d problem is a 2-D version of acoustic_wave_1d, also
19 : scaled for real arithmetic.
20 : */
21 :
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";
26 :
27 : #include <slepcpep.h>
28 :
29 6 : int main(int argc,char **argv)
30 : {
31 6 : Mat M,C,K,A[3]; /* problem matrices */
32 6 : PEP pep; /* polynomial eigenproblem solver context */
33 6 : PetscInt m=6,n,II,Istart,Iend,i,j;
34 6 : PetscScalar z=1.0;
35 6 : PetscReal h;
36 6 : char str[50];
37 6 : PetscBool terse;
38 :
39 6 : PetscFunctionBeginUser;
40 6 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
41 :
42 6 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
43 6 : PetscCheck(m>1,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"m must be at least 2");
44 6 : PetscCall(PetscOptionsGetScalar(NULL,NULL,"-z",&z,NULL));
45 6 : h = 1.0/m;
46 6 : n = m*(m-1);
47 6 : PetscCall(SlepcSNPrintfScalar(str,sizeof(str),z,PETSC_FALSE));
48 6 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nAcoustic wave 2-D, n=%" PetscInt_FMT " (m=%" PetscInt_FMT "), z=%s\n\n",n,m,str));
49 :
50 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51 : Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
52 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53 :
54 : /* K has a pattern similar to the 2D Laplacian */
55 6 : PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
56 6 : PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
57 6 : PetscCall(MatSetFromOptions(K));
58 :
59 6 : PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
60 186 : for (II=Istart;II<Iend;II++) {
61 180 : i = II/m; j = II-i*m;
62 300 : if (i>0) PetscCall(MatSetValue(K,II,II-m,(j==m-1)?-0.5:-1.0,INSERT_VALUES));
63 300 : if (i<m-2) PetscCall(MatSetValue(K,II,II+m,(j==m-1)?-0.5:-1.0,INSERT_VALUES));
64 180 : if (j>0) PetscCall(MatSetValue(K,II,II-1,-1.0,INSERT_VALUES));
65 180 : if (j<m-1) PetscCall(MatSetValue(K,II,II+1,-1.0,INSERT_VALUES));
66 330 : PetscCall(MatSetValue(K,II,II,(j==m-1)?2.0:4.0,INSERT_VALUES));
67 : }
68 :
69 6 : PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
70 6 : PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
71 :
72 : /* C is the zero matrix except for a few nonzero elements on the diagonal */
73 6 : PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
74 6 : PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n));
75 6 : PetscCall(MatSetFromOptions(C));
76 :
77 6 : PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
78 186 : for (i=Istart;i<Iend;i++) {
79 180 : if (i%m==m-1) PetscCall(MatSetValue(C,i,i,-2*PETSC_PI*h/z,INSERT_VALUES));
80 : }
81 6 : PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
82 6 : PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
83 :
84 : /* M is a diagonal matrix */
85 6 : PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
86 6 : PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
87 6 : PetscCall(MatSetFromOptions(M));
88 :
89 6 : PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
90 186 : for (i=Istart;i<Iend;i++) {
91 180 : if (i%m==m-1) PetscCall(MatSetValue(M,i,i,2*PETSC_PI*PETSC_PI*h*h,INSERT_VALUES));
92 180 : else PetscCall(MatSetValue(M,i,i,4*PETSC_PI*PETSC_PI*h*h,INSERT_VALUES));
93 : }
94 6 : PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
95 6 : PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
96 :
97 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98 : Create the eigensolver and solve the problem
99 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100 :
101 6 : PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
102 6 : A[0] = K; A[1] = C; A[2] = M;
103 6 : PetscCall(PEPSetOperators(pep,3,A));
104 6 : PetscCall(PEPSetFromOptions(pep));
105 6 : PetscCall(PEPSolve(pep));
106 :
107 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108 : Display solution and clean up
109 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110 :
111 : /* show detailed info unless -terse option is given by user */
112 6 : PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
113 6 : if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
114 : else {
115 0 : PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
116 0 : PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
117 0 : PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD));
118 0 : PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
119 : }
120 6 : PetscCall(PEPDestroy(&pep));
121 6 : PetscCall(MatDestroy(&M));
122 6 : PetscCall(MatDestroy(&C));
123 6 : PetscCall(MatDestroy(&K));
124 6 : PetscCall(SlepcFinalize());
125 : return 0;
126 : }
127 :
128 : /*TEST
129 :
130 : testset:
131 : args: -pep_nev 2 -pep_ncv 18 -terse
132 : output_file: output/acoustic_wave_2d_1.out
133 : filter: sed -e "s/2.60936i/2.60937i/g" | sed -e "s/2.60938i/2.60937i/g"
134 : test:
135 : suffix: 1
136 : args: -pep_type {{qarnoldi linear}}
137 : test:
138 : suffix: 1_toar
139 : args: -pep_type toar -pep_toar_locking 0
140 :
141 : testset:
142 : args: -pep_nev 2 -pep_ncv 18 -pep_type stoar -pep_hermitian -pep_scale scalar -st_type sinvert -terse
143 : output_file: output/acoustic_wave_2d_2.out
144 : test:
145 : suffix: 2
146 : test:
147 : suffix: 2_lin_b
148 : args: -pep_stoar_linearization 0,1
149 : test:
150 : suffix: 2_lin_ab
151 : args: -pep_stoar_linearization 0.1,0.9
152 :
153 : TEST*/
|