Actual source code: spring.c
slepc-3.20.2 2024-03-15
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 spring problem is a QEP from the finite element model of a damped
19: mass-spring system. This implementation supports only scalar parameters,
20: that is all masses, dampers and springs have the same constants.
21: Furthermore, this implementation does not consider different constants
22: for dampers and springs connecting adjacent masses or masses to the ground.
23: */
25: static char help[] = "FEM model of a damped mass-spring system.\n\n"
26: "The command line options are:\n"
27: " -n <n> ... dimension of the matrices.\n"
28: " -mu <value> ... mass (default 1).\n"
29: " -tau <value> ... damping constant of the dampers (default 10).\n"
30: " -kappa <value> ... damping constant of the springs (default 5).\n\n";
32: #include <slepcpep.h>
34: int main(int argc,char **argv)
35: {
36: Mat M,C,K,A[3]; /* problem matrices */
37: PEP pep; /* polynomial eigenproblem solver context */
38: PetscInt n=5,Istart,Iend,i;
39: PetscReal mu=1.0,tau=10.0,kappa=5.0;
40: PetscBool terse;
42: PetscFunctionBeginUser;
43: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
45: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
46: PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL));
47: PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
48: PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL));
49: 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));
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: /* K is a tridiagonal */
56: PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
57: PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
58: PetscCall(MatSetFromOptions(K));
59: PetscCall(MatSetUp(K));
61: PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
62: for (i=Istart;i<Iend;i++) {
63: if (i>0) PetscCall(MatSetValue(K,i,i-1,-kappa,INSERT_VALUES));
64: PetscCall(MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES));
65: if (i<n-1) PetscCall(MatSetValue(K,i,i+1,-kappa,INSERT_VALUES));
66: }
68: PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
69: PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
71: /* C is a tridiagonal */
72: PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
73: PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n));
74: PetscCall(MatSetFromOptions(C));
75: PetscCall(MatSetUp(C));
77: PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
78: for (i=Istart;i<Iend;i++) {
79: if (i>0) PetscCall(MatSetValue(C,i,i-1,-tau,INSERT_VALUES));
80: PetscCall(MatSetValue(C,i,i,tau*3.0,INSERT_VALUES));
81: if (i<n-1) PetscCall(MatSetValue(C,i,i+1,-tau,INSERT_VALUES));
82: }
84: PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
85: PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
87: /* M is a diagonal matrix */
88: PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
89: PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
90: PetscCall(MatSetFromOptions(M));
91: PetscCall(MatSetUp(M));
92: PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
93: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(M,i,i,mu,INSERT_VALUES));
94: PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
95: PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
97: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98: Create the eigensolver and solve the problem
99: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
102: A[0] = K; A[1] = C; A[2] = M;
103: PetscCall(PEPSetOperators(pep,3,A));
104: PetscCall(PEPSetFromOptions(pep));
105: PetscCall(PEPSolve(pep));
107: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108: Display solution and clean up
109: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111: /* show detailed info unless -terse option is given by user */
112: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
113: if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
114: else {
115: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
116: PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
117: PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD));
118: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
119: }
120: PetscCall(PEPDestroy(&pep));
121: PetscCall(MatDestroy(&M));
122: PetscCall(MatDestroy(&C));
123: PetscCall(MatDestroy(&K));
124: PetscCall(SlepcFinalize());
125: return 0;
126: }
128: /*TEST
130: testset:
131: args: -pep_nev 4 -n 24 -pep_ncv 18 -pep_target -.5 -st_type sinvert -pep_scale diagonal -terse
132: output_file: output/spring_1.out
133: filter: sed -e "s/[+-]0\.0*i//g"
134: test:
135: suffix: 1
136: args: -pep_type {{toar linear}} -pep_conv_norm
137: test:
138: suffix: 1_stoar
139: args: -pep_type stoar -pep_hermitian -pep_conv_rel
140: test:
141: suffix: 1_qarnoldi
142: args: -pep_type qarnoldi -pep_conv_rel
143: test:
144: suffix: 1_cuda
145: args: -mat_type aijcusparse
146: requires: cuda
148: test:
149: suffix: 2
150: args: -pep_type jd -pep_jd_minimality_index 1 -pep_nev 4 -n 24 -pep_ncv 18 -pep_target -50 -terse
151: requires: !single
152: filter: sed -e "s/[+-]0\.0*i//g"
154: test:
155: suffix: 3
156: args: -n 300 -pep_hermitian -pep_interval -10.1,-9.5 -pep_type stoar -st_type sinvert -st_pc_type cholesky -terse
157: filter: sed -e "s/52565/52566/" | sed -e "s/90758/90759/"
158: requires: !single
160: test:
161: suffix: 4
162: args: -n 300 -pep_hyperbolic -pep_interval -9.6,-.527 -pep_type stoar -st_type sinvert -st_pc_type cholesky -terse
163: requires: !single
164: timeoutfactor: 2
166: test:
167: suffix: 5
168: args: -n 300 -pep_hyperbolic -pep_interval -.506,-.3 -pep_type stoar -st_type sinvert -st_pc_type cholesky -pep_stoar_nev 11 -terse
169: requires: !single
171: test:
172: suffix: 6
173: args: -n 24 -pep_ncv 18 -pep_target -.5 -terse -pep_type jd -pep_jd_restart .6 -pep_jd_fix .001
174: requires: !single
176: TEST*/