Actual source code: test12.c
slepc-3.17.0 2022-03-31
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: */
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";
18: #include <slepcpep.h>
20: int main(int argc,char **argv)
21: {
22: Mat M,C,K,A[3]; /* problem matrices */
23: PEP pep; /* polynomial eigenproblem solver context */
24: RG rg;
25: PetscInt n=30,Istart,Iend,i;
26: PetscReal mu=1.0,tau=10.0,kappa=5.0;
28: SlepcInitialize(&argc,&argv,(char*)0,help);
30: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
31: PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL);
32: PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
33: PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL);
34: 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: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
38: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40: /* K is a tridiagonal */
41: MatCreate(PETSC_COMM_WORLD,&K);
42: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
43: MatSetFromOptions(K);
44: MatSetUp(K);
46: MatGetOwnershipRange(K,&Istart,&Iend);
47: for (i=Istart;i<Iend;i++) {
48: if (i>0) MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
49: MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
50: if (i<n-1) MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
51: }
53: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
54: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
56: /* C is a tridiagonal */
57: MatCreate(PETSC_COMM_WORLD,&C);
58: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
59: MatSetFromOptions(C);
60: MatSetUp(C);
62: MatGetOwnershipRange(C,&Istart,&Iend);
63: for (i=Istart;i<Iend;i++) {
64: if (i>0) MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
65: MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
66: if (i<n-1) MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
67: }
69: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
70: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
72: /* M is a diagonal matrix */
73: MatCreate(PETSC_COMM_WORLD,&M);
74: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
75: MatSetFromOptions(M);
76: MatSetUp(M);
77: MatGetOwnershipRange(M,&Istart,&Iend);
78: for (i=Istart;i<Iend;i++) MatSetValue(M,i,i,mu,INSERT_VALUES);
79: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
80: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Create a region for filtering
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: RGCreate(PETSC_COMM_WORLD,&rg);
87: RGSetType(rg,RGINTERVAL);
88: #if defined(PETSC_USE_COMPLEX)
89: RGIntervalSetEndpoints(rg,-47.0,-35.0,-0.001,0.001);
90: #else
91: RGIntervalSetEndpoints(rg,-47.0,-35.0,-0.0,0.0);
92: #endif
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Create the eigensolver and solve the problem
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: PEPCreate(PETSC_COMM_WORLD,&pep);
99: PEPSetRG(pep,rg);
100: A[0] = K; A[1] = C; A[2] = M;
101: PEPSetOperators(pep,3,A);
102: PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT);
103: PEPSetFromOptions(pep);
104: PEPSolve(pep);
106: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107: Display solution and clean up
108: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
111: PEPDestroy(&pep);
112: MatDestroy(&M);
113: MatDestroy(&C);
114: MatDestroy(&K);
115: RGDestroy(&rg);
116: SlepcFinalize();
117: return 0;
118: }
120: /*TEST
122: test:
123: args: -pep_nev 8 -pep_type {{toar linear qarnoldi}}
124: suffix: 1
125: requires: !single
126: output_file: output/test12_1.out
128: TEST*/