Actual source code: sleeper.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 sleeper problem is a proportionally damped QEP describing the
 30:    oscillations of a rail track resting on sleepers.
 31: */

 33: static char help[] = "NLEVP problem: sleeper.\n\n"
 34:   "The command line options are:\n"
 35:   "  -n <n>, where <n> = dimension of the matrices.\n\n";

 37: #include <slepcpep.h>

 41: int main(int argc,char **argv)
 42: {
 43:   Mat            M,C,K,A[3];      /* problem matrices */
 44:   PEP            pep;             /* polynomial eigenproblem solver context */
 45:   PetscInt       n=10,Istart,Iend,i;

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

 50:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 51:   PetscPrintf(PETSC_COMM_WORLD,"\nRailtrack resting on sleepers, n=%D\n\n",n);

 53:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 54:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 55:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 57:   /* K is a tridiagonal */
 58:   MatCreate(PETSC_COMM_WORLD,&K);
 59:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
 60:   MatSetFromOptions(K);
 61:   MatSetUp(K);
 62:   
 63:   MatGetOwnershipRange(K,&Istart,&Iend);
 64:   for (i=Istart;i<Iend;i++) {
 65:     if (i==0) { 
 66:       MatSetValue(K,i,n-1,-3.0,INSERT_VALUES);
 67:       MatSetValue(K,i,n-2,1.0,INSERT_VALUES);
 68:     }
 69:     if (i==1) { MatSetValue(K,i,n-1,1.0,INSERT_VALUES); }
 70:     if (i>0) { MatSetValue(K,i,i-1,-3.0,INSERT_VALUES); }
 71:     if (i>1) { MatSetValue(K,i,i-2,1.0,INSERT_VALUES); }
 72:     MatSetValue(K,i,i,5.0,INSERT_VALUES);
 73:     if (i==n-1) { 
 74:       MatSetValue(K,i,0,-3.0,INSERT_VALUES);
 75:       MatSetValue(K,i,1,1.0,INSERT_VALUES);
 76:     }
 77:     if (i==n-2) { MatSetValue(K,i,0,1.0,INSERT_VALUES); }
 78:     if (i<n-1) { MatSetValue(K,i,i+1,-3.0,INSERT_VALUES); }
 79:     if (i<n-2) { MatSetValue(K,i,i+2,1.0,INSERT_VALUES); }
 80:   }

 82:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 83:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 85:   /* C is a circulant matrix */
 86:   MatCreate(PETSC_COMM_WORLD,&C);
 87:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 88:   MatSetFromOptions(C);
 89:   MatSetUp(C);
 90:   
 91:   MatGetOwnershipRange(C,&Istart,&Iend);
 92:   for (i=Istart;i<Iend;i++) {
 93:     if (i==0) { 
 94:       MatSetValue(C,i,n-1,-4.0,INSERT_VALUES);
 95:       MatSetValue(C,i,n-2,1.0,INSERT_VALUES);
 96:     }
 97:     if (i==1) { MatSetValue(C,i,n-1,1.0,INSERT_VALUES); }
 98:     if (i>0) { MatSetValue(C,i,i-1,-4.0,INSERT_VALUES); }
 99:     if (i>1) { MatSetValue(C,i,i-2,1.0,INSERT_VALUES); }
100:     MatSetValue(C,i,i,7.0,INSERT_VALUES);
101:     if (i==n-1) { 
102:       MatSetValue(C,i,0,-4.0,INSERT_VALUES);
103:       MatSetValue(C,i,1,1.0,INSERT_VALUES);
104:     }
105:     if (i==n-2) { MatSetValue(C,i,0,1.0,INSERT_VALUES); }
106:     if (i<n-1) { MatSetValue(C,i,i+1,-4.0,INSERT_VALUES); }
107:     if (i<n-2) { MatSetValue(C,i,i+2,1.0,INSERT_VALUES); }
108:   }

110:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
111:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
112:   
113:   /* M is the identity matrix */
114:   MatCreate(PETSC_COMM_WORLD,&M);
115:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
116:   MatSetFromOptions(M);
117:   MatSetUp(M);
118:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
119:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
120:   MatShift(M,1.0);

122:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
123:                 Create the eigensolver and solve the problem
124:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

126:   PEPCreate(PETSC_COMM_WORLD,&pep);
127:   A[0] = K; A[1] = C; A[2] = M;
128:   PEPSetOperators(pep,3,A);
129:   PEPSetFromOptions(pep);
130:   PEPSolve(pep);

132:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133:                     Display solution and clean up
134:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135:   
136:   PEPPrintSolution(pep,NULL);
137:   PEPDestroy(&pep);
138:   MatDestroy(&M);
139:   MatDestroy(&C);
140:   MatDestroy(&K);
141:   SlepcFinalize();
142:   return 0;
143: }