Actual source code: loaded_string.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 loaded_string problem is a rational eigenvalue problem for the
 30:    finite element model of a loaded vibrating string.
 31: */

 33: static char help[] = "NLEVP problem: loaded_string.\n\n"
 34:   "The command line options are:\n"
 35:   "  -n <n>, dimension of the matrices.\n"
 36:   "  -kappa <kappa>, stiffness of elastic spring.\n"
 37:   "  -mass <m>, mass of the attached load.\n\n";

 39: #include <slepcnep.h>

 41: #define NMAT 3

 45: int main(int argc,char **argv)
 46: {
 47:   Mat            A[NMAT];         /* problem matrices */
 48:   FN             f[NMAT];         /* functions to define the nonlinear operator */
 49:   NEP            nep;             /* nonlinear eigensolver context */
 50:   PetscInt       n=20,Istart,Iend,i,nconv;
 51:   PetscReal      kappa=1.0,m=1.0,re,im,norm;
 52:   PetscScalar    kr,ki,sigma,numer[2],denom[2];

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

 57:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 58:   PetscOptionsGetReal(NULL,"-kappa",&kappa,NULL);
 59:   PetscOptionsGetReal(NULL,"-mass",&m,NULL);
 60:   sigma = kappa/m;
 61:   PetscPrintf(PETSC_COMM_WORLD,"Loaded vibrating string, n=%D kappa=%g m=%g\n\n",n,(double)kappa,(double)m);

 63:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 64:                        Build the problem matrices 
 65:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 67:   /* initialize matrices */
 68:   for (i=0;i<NMAT;i++) {
 69:     MatCreate(PETSC_COMM_WORLD,&A[i]);
 70:     MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n);
 71:     MatSetFromOptions(A[i]);
 72:     MatSetUp(A[i]);
 73:   }
 74:   MatGetOwnershipRange(A[0],&Istart,&Iend);

 76:   /* A0 */
 77:   for (i=Istart;i<Iend;i++) {
 78:     MatSetValue(A[0],i,i,(i==n-1)?1.0*n:2.0*n,INSERT_VALUES);
 79:     if (i>0) { MatSetValue(A[0],i,i-1,-1.0*n,INSERT_VALUES); }
 80:     if (i<n-1) { MatSetValue(A[0],i,i+1,-1.0*n,INSERT_VALUES); }
 81:   }

 83:   /* A1 */

 85:   for (i=Istart;i<Iend;i++) {
 86:     MatSetValue(A[1],i,i,(i==n-1)?2.0/(6.0*n):4.0/(6.0*n),INSERT_VALUES);
 87:     if (i>0) { MatSetValue(A[1],i,i-1,1.0/(6.0*n),INSERT_VALUES); }
 88:     if (i<n-1) { MatSetValue(A[1],i,i+1,1.0/(6.0*n),INSERT_VALUES); }
 89:   }

 91:   /* A2 */
 92:   if (Istart<=n-1 && n-1<Iend) {
 93:     MatSetValue(A[2],n-1,n-1,kappa,INSERT_VALUES); 
 94:   }

 96:   /* assemble matrices */
 97:   for (i=0;i<NMAT;i++) {
 98:     MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY);
 99:   }
100:   for (i=0;i<NMAT;i++) {
101:     MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY);
102:   }

104:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
105:                        Create the problem functions
106:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

108:   /* f1=1 */
109:   FNCreate(PETSC_COMM_WORLD,&f[0]);
110:   FNSetType(f[0],FNRATIONAL);
111:   numer[0] = 1.0;
112:   FNSetParameters(f[0],1,numer,0,NULL);

114:   /* f2=-lambda */
115:   FNCreate(PETSC_COMM_WORLD,&f[1]);
116:   FNSetType(f[1],FNRATIONAL);
117:   numer[0] = -1.0; numer[1] = 0.0;
118:   FNSetParameters(f[1],2,numer,0,NULL);

120:   /* f3=lambda/(lambda-sigma) */
121:   FNCreate(PETSC_COMM_WORLD,&f[2]);
122:   FNSetType(f[2],FNRATIONAL);
123:   numer[0] = 1.0; numer[1] = 0.0;
124:   denom[0] = 1.0; denom[1] = -sigma;
125:   FNSetParameters(f[2],2,numer,2,denom);

127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
128:                 Create the eigensolver and solve the problem
129:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

131:   NEPCreate(PETSC_COMM_WORLD,&nep);
132:   NEPSetSplitOperator(nep,3,A,f,SUBSET_NONZERO_PATTERN);
133:   NEPSetFromOptions(nep);
134:   NEPSolve(nep);

136:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137:                     Display solution and clean up
138:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139:   
140:   /*
141:      Get number of converged approximate eigenpairs
142:   */
143:   NEPGetConverged(nep,&nconv);
144:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %D\n\n",nconv);

146:   if (nconv>0) {
147:     /*
148:        Display eigenvalues and relative errors
149:     */
150:     PetscPrintf(PETSC_COMM_WORLD,
151:          "           k              ||T(k)x||\n"
152:          "   ----------------- ------------------\n");
153:     for (i=0;i<nconv;i++) {
154:       NEPGetEigenpair(nep,i,&kr,&ki,NULL,NULL);
155:       NEPComputeRelativeError(nep,i,&norm);
156: #if defined(PETSC_USE_COMPLEX)
157:       re = PetscRealPart(kr);
158:       im = PetscImaginaryPart(kr);
159: #else
160:       re = kr;
161:       im = ki;
162: #endif
163:       if (im!=0.0) {
164:         PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g\n",(double)re,(double)im,(double)norm);
165:       } else {
166:         PetscPrintf(PETSC_COMM_WORLD,"   %12f         %12g\n",(double)re,(double)norm);
167:       }
168:     }
169:     PetscPrintf(PETSC_COMM_WORLD,"\n");
170:   }

172:   NEPDestroy(&nep);
173:   for (i=0;i<NMAT;i++) {
174:     MatDestroy(&A[i]);
175:     FNDestroy(&f[i]);
176:   }
177:   SlepcFinalize();
178:   return 0;
179: }