Actual source code: spring.c
slepc-main 2024-11-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,NULL,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));
60: PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
61: for (i=Istart;i<Iend;i++) {
62: if (i>0) PetscCall(MatSetValue(K,i,i-1,-kappa,INSERT_VALUES));
63: PetscCall(MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES));
64: if (i<n-1) PetscCall(MatSetValue(K,i,i+1,-kappa,INSERT_VALUES));
65: }
67: PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
68: PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
70: /* C is a tridiagonal */
71: PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
72: PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n));
73: PetscCall(MatSetFromOptions(C));
75: PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
76: for (i=Istart;i<Iend;i++) {
77: if (i>0) PetscCall(MatSetValue(C,i,i-1,-tau,INSERT_VALUES));
78: PetscCall(MatSetValue(C,i,i,tau*3.0,INSERT_VALUES));
79: if (i<n-1) PetscCall(MatSetValue(C,i,i+1,-tau,INSERT_VALUES));
80: }
82: PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
83: PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
85: /* M is a diagonal matrix */
86: PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
87: PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
88: PetscCall(MatSetFromOptions(M));
89: PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
90: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(M,i,i,mu,INSERT_VALUES));
91: PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
92: PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Create the eigensolver and solve the problem
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
99: A[0] = K; A[1] = C; A[2] = M;
100: PetscCall(PEPSetOperators(pep,3,A));
101: PetscCall(PEPSetFromOptions(pep));
102: PetscCall(PEPSolve(pep));
104: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: Display solution and clean up
106: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: /* show detailed info unless -terse option is given by user */
109: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
110: if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
111: else {
112: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
113: PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
114: PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD));
115: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
116: }
117: PetscCall(PEPDestroy(&pep));
118: PetscCall(MatDestroy(&M));
119: PetscCall(MatDestroy(&C));
120: PetscCall(MatDestroy(&K));
121: PetscCall(SlepcFinalize());
122: return 0;
123: }
125: /*TEST
127: testset:
128: args: -pep_nev 4 -n 24 -pep_ncv 18 -pep_target -.5 -st_type sinvert -pep_scale diagonal -terse
129: output_file: output/spring_1.out
130: filter: sed -e "s/[+-]0\.0*i//g"
131: test:
132: suffix: 1
133: args: -pep_type {{toar linear}} -pep_conv_norm
134: test:
135: suffix: 1_stoar
136: args: -pep_type stoar -pep_hermitian -pep_conv_rel
137: test:
138: suffix: 1_qarnoldi
139: args: -pep_type qarnoldi -pep_conv_rel
140: test:
141: suffix: 1_cuda
142: args: -mat_type aijcusparse
143: requires: cuda
144: test:
145: suffix: 1_hip
146: args: -mat_type aijhipsparse
147: requires: hip
149: test:
150: suffix: 2
151: args: -pep_type jd -pep_jd_minimality_index 1 -pep_nev 4 -n 24 -pep_ncv 18 -pep_target -50 -terse
152: requires: !single
153: filter: sed -e "s/[+-]0\.0*i//g"
155: test:
156: suffix: 3
157: args: -n 300 -pep_hermitian -pep_interval -10.1,-9.5 -pep_type stoar -st_type sinvert -st_pc_type cholesky -terse
158: filter: sed -e "s/52565/52566/" | sed -e "s/90758/90759/"
159: requires: !single
161: test:
162: suffix: 4
163: args: -n 300 -pep_hyperbolic -pep_interval -9.6,-.527 -pep_type stoar -st_type sinvert -st_pc_type cholesky -terse
164: requires: !single
165: timeoutfactor: 2
167: test:
168: suffix: 5
169: args: -n 300 -pep_hyperbolic -pep_interval -.506,-.3 -pep_type stoar -st_type sinvert -st_pc_type cholesky -pep_stoar_nev 11 -terse
170: requires: !single
172: test:
173: suffix: 6
174: args: -n 24 -pep_ncv 18 -pep_target -.5 -terse -pep_type jd -pep_jd_restart .6 -pep_jd_fix .001
175: requires: !single
177: TEST*/