Actual source code: spring.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */
 21: /*
 22:    This example implements one of the problems found at
 23:        NLEVP: A Collection of Nonlinear Eigenvalue Problems,
 24:        The University of Manchester.
 25:    The details of the collection can be found at:
 26:        [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
 27:            Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.

 29:    The spring problem is a QEP from the finite element model of a damped
 30:    mass-spring system. This implementation supports only scalar parameters,
 31:    that is all masses, dampers and springs have the same constants.
 32:    Furthermore, this implementation does not consider different constants
 33:    for dampers and springs connecting adjacent masses or masses to the ground.
 34: */

 36: static char help[] = "NLEVP problem: spring.\n\n"
 37:   "The command line options are:\n"
 38:   "  -n <n> ... dimension of the matrices.\n"
 39:   "  -mu <value> ... mass (default 1).\n"
 40:   "  -tau <value> ... damping constant of the dampers (default 10).\n"
 41:   "  -kappa <value> ... damping constant of the springs (default 5).\n\n";

 43: #include <slepcpep.h>

 47: int main(int argc,char **argv)
 48: {
 49:   Mat            M,C,K,A[3];      /* problem matrices */
 50:   PEP            pep;             /* polynomial eigenproblem solver context */
 51:   PetscInt       n=5,Istart,Iend,i;
 52:   PetscReal      mu=1,tau=10,kappa=5;

 55:   SlepcInitialize(&argc,&argv,(char*)0,help);

 57:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 58:   PetscOptionsGetReal(NULL,"-mu",&mu,NULL);
 59:   PetscOptionsGetReal(NULL,"-tau",&tau,NULL);
 60:   PetscOptionsGetReal(NULL,"-kappa",&kappa,NULL);
 61:   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);

 63:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 64:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 65:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 67:   /* K is a tridiagonal */
 68:   MatCreate(PETSC_COMM_WORLD,&K);
 69:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
 70:   MatSetFromOptions(K);
 71:   MatSetUp(K);
 72:   
 73:   MatGetOwnershipRange(K,&Istart,&Iend);
 74:   for (i=Istart;i<Iend;i++) {
 75:     if (i>0) {
 76:       MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
 77:     }
 78:     MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
 79:     if (i<n-1) {
 80:       MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
 81:     }
 82:   }

 84:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 85:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 87:   /* C is a tridiagonal */
 88:   MatCreate(PETSC_COMM_WORLD,&C);
 89:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 90:   MatSetFromOptions(C);
 91:   MatSetUp(C);
 92:   
 93:   MatGetOwnershipRange(C,&Istart,&Iend);
 94:   for (i=Istart;i<Iend;i++) {
 95:     if (i>0) {
 96:       MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
 97:     }
 98:     MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
 99:     if (i<n-1) {
100:       MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
101:     }
102:   }

104:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
105:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
106:   
107:   /* M is a diagonal matrix */
108:   MatCreate(PETSC_COMM_WORLD,&M);
109:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
110:   MatSetFromOptions(M);
111:   MatSetUp(M);
112:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
113:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
114:   MatShift(M,mu);

116:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
117:                 Create the eigensolver and solve the problem
118:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

120:   PEPCreate(PETSC_COMM_WORLD,&pep);
121:   A[0] = K; A[1] = C; A[2] = M;
122:   PEPSetOperators(pep,3,A);
123:   PEPSetFromOptions(pep);
124:   PEPSolve(pep);

126:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127:                     Display solution and clean up
128:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129:   
130:   PEPPrintSolution(pep,NULL);
131:   PEPDestroy(&pep);
132:   MatDestroy(&M);
133:   MatDestroy(&C);
134:   MatDestroy(&K);
135:   SlepcFinalize();
136:   return 0;
137: }