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 : static char help[] = "Illustrates region filtering in PEP (based on spring).\n"
12 : "The command line options are:\n"
13 : " -n <n> ... number of grid subdivisions.\n"
14 : " -mu <value> ... mass (default 1).\n"
15 : " -tau <value> ... damping constant of the dampers (default 10).\n"
16 : " -kappa <value> ... damping constant of the springs (default 5).\n\n";
17 :
18 : #include <slepcpep.h>
19 :
20 3 : int main(int argc,char **argv)
21 : {
22 3 : Mat M,C,K,A[3]; /* problem matrices */
23 3 : PEP pep; /* polynomial eigenproblem solver context */
24 3 : RG rg;
25 3 : PetscInt n=30,Istart,Iend,i;
26 3 : PetscReal mu=1.0,tau=10.0,kappa=5.0;
27 :
28 3 : PetscFunctionBeginUser;
29 3 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
30 :
31 3 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
32 3 : PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL));
33 3 : PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
34 3 : PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL));
35 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nDamped mass-spring system, n=%" PetscInt_FMT " mu=%g tau=%g kappa=%g\n\n",n,(double)mu,(double)tau,(double)kappa));
36 :
37 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38 : Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
39 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40 :
41 : /* K is a tridiagonal */
42 3 : PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
43 3 : PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
44 3 : PetscCall(MatSetFromOptions(K));
45 :
46 3 : PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
47 93 : for (i=Istart;i<Iend;i++) {
48 90 : if (i>0) PetscCall(MatSetValue(K,i,i-1,-kappa,INSERT_VALUES));
49 90 : PetscCall(MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES));
50 90 : if (i<n-1) PetscCall(MatSetValue(K,i,i+1,-kappa,INSERT_VALUES));
51 : }
52 :
53 3 : PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
54 3 : PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
55 :
56 : /* C is a tridiagonal */
57 3 : PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
58 3 : PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n));
59 3 : PetscCall(MatSetFromOptions(C));
60 :
61 3 : PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
62 93 : for (i=Istart;i<Iend;i++) {
63 90 : if (i>0) PetscCall(MatSetValue(C,i,i-1,-tau,INSERT_VALUES));
64 90 : PetscCall(MatSetValue(C,i,i,tau*3.0,INSERT_VALUES));
65 90 : if (i<n-1) PetscCall(MatSetValue(C,i,i+1,-tau,INSERT_VALUES));
66 : }
67 :
68 3 : PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
69 3 : PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
70 :
71 : /* M is a diagonal matrix */
72 3 : PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
73 3 : PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
74 3 : PetscCall(MatSetFromOptions(M));
75 3 : PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
76 93 : for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(M,i,i,mu,INSERT_VALUES));
77 3 : PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
78 3 : PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
79 :
80 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81 : Create a region for filtering
82 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83 :
84 3 : PetscCall(RGCreate(PETSC_COMM_WORLD,&rg));
85 3 : PetscCall(RGSetType(rg,RGINTERVAL));
86 : #if defined(PETSC_USE_COMPLEX)
87 : PetscCall(RGIntervalSetEndpoints(rg,-47.0,-35.0,-0.001,0.001));
88 : #else
89 3 : PetscCall(RGIntervalSetEndpoints(rg,-47.0,-35.0,-0.0,0.0));
90 : #endif
91 :
92 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93 : Create the eigensolver and solve the problem
94 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95 :
96 3 : PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
97 3 : PetscCall(PEPSetRG(pep,rg));
98 3 : A[0] = K; A[1] = C; A[2] = M;
99 3 : PetscCall(PEPSetOperators(pep,3,A));
100 3 : PetscCall(PEPSetTolerances(pep,PETSC_SMALL,PETSC_CURRENT));
101 3 : PetscCall(PEPSetFromOptions(pep));
102 3 : PetscCall(PEPSolve(pep));
103 :
104 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105 : Display solution and clean up
106 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107 :
108 3 : PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
109 3 : PetscCall(PEPDestroy(&pep));
110 3 : PetscCall(MatDestroy(&M));
111 3 : PetscCall(MatDestroy(&C));
112 3 : PetscCall(MatDestroy(&K));
113 3 : PetscCall(RGDestroy(&rg));
114 3 : PetscCall(SlepcFinalize());
115 : return 0;
116 : }
117 :
118 : /*TEST
119 :
120 : test:
121 : args: -pep_nev 8 -pep_type {{toar linear qarnoldi}}
122 : suffix: 1
123 : requires: !single
124 : output_file: output/test12_1.out
125 :
126 : TEST*/
|