Actual source code: spring.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  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:   SlepcInitialize(&argc,&argv,(char*)0,help);

 44:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 45:   PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL);
 46:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 47:   PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL);
 48:   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);

 50:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 51:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 52:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 54:   /* K is a tridiagonal */
 55:   MatCreate(PETSC_COMM_WORLD,&K);
 56:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
 57:   MatSetFromOptions(K);
 58:   MatSetUp(K);

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

 67:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 68:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 70:   /* C is a tridiagonal */
 71:   MatCreate(PETSC_COMM_WORLD,&C);
 72:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 73:   MatSetFromOptions(C);
 74:   MatSetUp(C);

 76:   MatGetOwnershipRange(C,&Istart,&Iend);
 77:   for (i=Istart;i<Iend;i++) {
 78:     if (i>0) MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
 79:     MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
 80:     if (i<n-1) MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
 81:   }

 83:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 84:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 86:   /* M is a diagonal matrix */
 87:   MatCreate(PETSC_COMM_WORLD,&M);
 88:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
 89:   MatSetFromOptions(M);
 90:   MatSetUp(M);
 91:   MatGetOwnershipRange(M,&Istart,&Iend);
 92:   for (i=Istart;i<Iend;i++) MatSetValue(M,i,i,mu,INSERT_VALUES);
 93:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
 94:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

 96:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97:                 Create the eigensolver and solve the problem
 98:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

100:   PEPCreate(PETSC_COMM_WORLD,&pep);
101:   A[0] = K; A[1] = C; A[2] = M;
102:   PEPSetOperators(pep,3,A);
103:   PEPSetFromOptions(pep);
104:   PEPSolve(pep);

106:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107:                     Display solution and clean up
108:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

110:   /* show detailed info unless -terse option is given by user */
111:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
112:   if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
113:   else {
114:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
115:     PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
116:     PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
117:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
118:   }
119:   PEPDestroy(&pep);
120:   MatDestroy(&M);
121:   MatDestroy(&C);
122:   MatDestroy(&K);
123:   SlepcFinalize();
124:   return 0;
125: }

127: /*TEST

129:    testset:
130:       args: -pep_nev 4 -n 24 -pep_ncv 18 -pep_target -.5 -st_type sinvert -pep_scale diagonal -terse
131:       output_file: output/spring_1.out
132:       filter: sed -e "s/[+-]0\.0*i//g"
133:       test:
134:          suffix: 1
135:          args: -pep_type {{toar linear}} -pep_conv_norm
136:       test:
137:          suffix: 1_stoar
138:          args: -pep_type stoar -pep_hermitian -pep_conv_rel
139:       test:
140:          suffix: 1_qarnoldi
141:          args: -pep_type qarnoldi -pep_conv_rel
142:       test:
143:          suffix: 1_cuda
144:          args: -mat_type aijcusparse
145:          requires: cuda

147:    test:
148:       suffix: 2
149:       args: -pep_type jd -pep_jd_minimality_index 1 -pep_nev 4 -n 24 -pep_ncv 18 -pep_target -50 -terse
150:       requires: !single
151:       filter: sed -e "s/[+-]0\.0*i//g"

153:    test:
154:       suffix: 3
155:       args: -n 300 -pep_hermitian -pep_interval -10.1,-9.5 -pep_type stoar -st_type sinvert -st_pc_type cholesky -terse
156:       filter: sed -e "s/52565/52566/" | sed -e "s/90758/90759/"

158:    test:
159:       suffix: 4
160:       args: -n 300 -pep_hyperbolic -pep_interval -9.6,-.527 -pep_type stoar -st_type sinvert -st_pc_type cholesky -terse
161:       requires: !single
162:       timeoutfactor: 2

164:    test:
165:       suffix: 5
166:       args: -n 300 -pep_hyperbolic -pep_interval -.506,-.3 -pep_type stoar -st_type sinvert -st_pc_type cholesky -pep_stoar_nev 11 -terse
167:       requires: !single

169:    test:
170:       suffix: 6
171:       args: -n 24 -pep_ncv 18 -pep_target -.5 -terse -pep_type jd -pep_jd_restart .6 -pep_jd_fix .001
172:       requires: !single

174: TEST*/