Actual source code: test2.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:       
  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:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 21:    Example based on spring problem in NLEVP collection [1]. See the parameters
 22:    meaning at Example 2 in [2].

 24:    [1] T. Betcke, N. J. Higham, V. Mehrmann, C. Schroder, and F. Tisseur,
 25:        NLEVP: A Collection of Nonlinear Eigenvalue Problems, MIMS EPrint
 26:        2010.98, November 2010.
 27:    [2] F. Tisseur, Backward error and condition of polynomial eigenvalue
 28:        problems, Linear Algebra and its Applications, 309 (2000), pp. 339--361,
 29:        April 2000.
 30: */

 32: static char help[] = "Test the solution of a QEP from a finite element model of "
 33:   "damped mass-spring system (problem from NLEVP collection).\n\n"
 34:   "The command line options are:\n"
 35:   "  -n <n>, where <n> = number of grid subdivisions.\n"
 36:   "  -tau <tau>, where <tau> = tau parameter.\n"
 37:   "  -kappa <kappa>, where <kappa> = kappa paramter.\n\n";

 39:  #include slepcqep.h

 43: int main( int argc, char **argv )
 44: {
 45:   Mat                  M, C, K;         /* problem matrices */
 46:   QEP                  qep;             /* quadratic eigenproblem solver context */
 47:   const QEPType  type;
 49:   PetscInt             n=30,Istart,Iend,i,maxit,nev;
 50:   PetscScalar          mu=1.0,tau=10.0,kappa=5.0;

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

 54:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 55:   PetscOptionsGetScalar(PETSC_NULL,"-mu",&mu,PETSC_NULL);
 56:   PetscOptionsGetScalar(PETSC_NULL,"-tau",&tau,PETSC_NULL);
 57:   PetscOptionsGetScalar(PETSC_NULL,"-kappa",&kappa,PETSC_NULL);

 59:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 60:      Compute the matrices that define the eigensystem, (k^2*K+k*X+M)x=0
 61:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

 79:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 80:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 82:   /* C is a tridiagonal */
 83:   MatCreate(PETSC_COMM_WORLD,&C);
 84:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 85:   MatSetFromOptions(C);
 86: 
 87:   MatGetOwnershipRange(C,&Istart,&Iend);
 88:   for (i=Istart; i<Iend; i++) {
 89:     if (i>0) {
 90:       MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
 91:     }
 92:     MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
 93:     if (i<n-1) {
 94:       MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
 95:     }
 96:   }

 98:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 99:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
100: 
101:   /* M is a diagonal matrix */
102:   MatCreate(PETSC_COMM_WORLD,&M);
103:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
104:   MatSetFromOptions(M);
105:   MatGetOwnershipRange(M,&Istart,&Iend);
106:   for (i=Istart; i<Iend; i++) {
107:     MatSetValue(M,i,i,mu,INSERT_VALUES);
108:   }
109:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
110:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
111: 
112:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
113:                 Create the eigensolver and set various options
114:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

116:   /* 
117:      Create eigensolver context
118:   */
119:   QEPCreate(PETSC_COMM_WORLD,&qep);

121:   /* 
122:      Set matrices, the problem type and other settings
123:   */
124:   QEPSetOperators(qep,M,C,K);
125:   QEPSetProblemType(qep,QEP_GENERAL);
126:   QEPSetTolerances(qep,PETSC_SMALL,PETSC_DEFAULT);
127:   QEPSetFromOptions(qep);

129:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
130:                       Solve the eigensystem
131:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

133:   QEPSolve(qep);

135:   /*
136:      Optional: Get some information from the solver and display it
137:   */
138:   QEPGetType(qep,&type);
139:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
140:   QEPGetDimensions(qep,&nev,PETSC_NULL,PETSC_NULL);
141:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
142:   QEPGetTolerances(qep,PETSC_NULL,&maxit);
143:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: maxit=%D\n",maxit);

145:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
146:                     Display solution and clean up
147:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

149:   QEPPrintSolution(qep,PETSC_NULL);
150: 
151:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
152:      Free work space
153:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

155:   QEPDestroy(&qep);
156:   MatDestroy(&M);
157:   MatDestroy(&C);
158:   MatDestroy(&K);
159:   SlepcFinalize();
160:   return 0;
161: }