Actual source code: spring.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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,tau=10,kappa=5;
 40:   PetscBool      terse;

 43:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 45:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 46:   PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL);
 47:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 48:   PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL);
 49:   PetscPrintf(PETSC_COMM_WORLD,"\nDamped mass-spring system, n=%D 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:   MatCreate(PETSC_COMM_WORLD,&K);
 57:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
 58:   MatSetFromOptions(K);
 59:   MatSetUp(K);

 61:   MatGetOwnershipRange(K,&Istart,&Iend);
 62:   for (i=Istart;i<Iend;i++) {
 63:     if (i>0) {
 64:       MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
 65:     }
 66:     MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
 67:     if (i<n-1) {
 68:       MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
 69:     }
 70:   }

 72:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 73:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 75:   /* C is a tridiagonal */
 76:   MatCreate(PETSC_COMM_WORLD,&C);
 77:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 78:   MatSetFromOptions(C);
 79:   MatSetUp(C);

 81:   MatGetOwnershipRange(C,&Istart,&Iend);
 82:   for (i=Istart;i<Iend;i++) {
 83:     if (i>0) {
 84:       MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
 85:     }
 86:     MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
 87:     if (i<n-1) {
 88:       MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
 89:     }
 90:   }

 92:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 93:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 95:   /* M is a diagonal matrix */
 96:   MatCreate(PETSC_COMM_WORLD,&M);
 97:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
 98:   MatSetFromOptions(M);
 99:   MatSetUp(M);
100:   MatGetOwnershipRange(M,&Istart,&Iend);
101:   for (i=Istart;i<Iend;i++) {
102:     MatSetValue(M,i,i,mu,INSERT_VALUES);
103:   }
104:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
105:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

107:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108:                 Create the eigensolver and solve the problem
109:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

111:   PEPCreate(PETSC_COMM_WORLD,&pep);
112:   A[0] = K; A[1] = C; A[2] = M;
113:   PEPSetOperators(pep,3,A);
114:   PEPSetFromOptions(pep);
115:   PEPSolve(pep);

117:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118:                     Display solution and clean up
119:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

121:   /* show detailed info unless -terse option is given by user */
122:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
123:   if (terse) {
124:     PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
125:   } else {
126:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
127:     PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
128:     PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
129:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
130:   }
131:   PEPDestroy(&pep);
132:   MatDestroy(&M);
133:   MatDestroy(&C);
134:   MatDestroy(&K);
135:   SlepcFinalize();
136:   return ierr;
137: }

139: /*TEST

141:    testset:
142:       args: -pep_nev 4 -n 24 -pep_ncv 18 -pep_target -.5 -st_type sinvert -pep_scale diagonal -terse
143:       requires: !single
144:       output_file: output/spring_1.out
145:       test:
146:          suffix: 1
147:          args: -pep_type {{toar linear}} -pep_conv_norm
148:       test:
149:          suffix: 1_stoar
150:          args: -pep_type stoar -pep_hermitian -pep_conv_rel
151:       test:
152:          suffix: 1_qarnoldi
153:          args: -pep_type qarnoldi -pep_conv_rel

155:    test:
156:       suffix: 2
157:       args: -pep_type jd -pep_jd_minimality_index 1 -pep_nev 4 -n 24 -pep_ncv 18 -pep_target -50 -terse
158:       requires: !single

160:    test:
161:       suffix: 3
162:       args: -n 300 -pep_hermitian -pep_interval -10.1,-9.5 -pep_type stoar -st_type sinvert -st_pc_type cholesky -terse
163:       requires: !single

165:    test:
166:       suffix: 4
167:       args: -n 300 -pep_hyperbolic -pep_interval -9.6,-.527 -pep_type stoar -st_type sinvert -st_pc_type cholesky -terse
168:       requires: !single

170:    test:
171:       suffix: 5
172:       args: -n 300 -pep_hyperbolic -pep_interval -.506,-.3 -pep_type stoar -st_type sinvert -st_pc_type cholesky -pep_stoar_nev 10 -terse
173:       requires: !single !complex

175:    test:
176:       suffix: 6
177:       args: -n 24 -pep_ncv 18 -pep_target -.5 -terse -pep_type jd -pep_jd_restart .6 -pep_jd_fix .001
178:       requires: !single
179: TEST*/