Actual source code: test6.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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:    Example based on spring problem in NLEVP collection [1]. See the parameters
 12:    meaning at Example 2 in [2].

 14:    [1] T. Betcke, N. J. Higham, V. Mehrmann, C. Schroder, and F. Tisseur,
 15:        NLEVP: A Collection of Nonlinear Eigenvalue Problems, MIMS EPrint
 16:        2010.98, November 2010.
 17:    [2] F. Tisseur, Backward error and condition of polynomial eigenvalue
 18:        problems, Linear Algebra and its Applications, 309 (2000), pp. 339--361,
 19:        April 2000.
 20: */

 22: static char help[] = "Tests multiple calls to PEPSolve with different matrix of different size.\n\n"
 23:   "This is based on the spring problem from NLEVP collection.\n\n"
 24:   "The command line options are:\n"
 25:   "  -n <n> ... number of grid subdivisions.\n"
 26:   "  -mu <value> ... mass (default 1).\n"
 27:   "  -tau <value> ... damping constant of the dampers (default 10).\n"
 28:   "  -kappa <value> ... damping constant of the springs (default 5).\n"
 29:   "  -initv ... set an initial vector.\n\n";

 31: #include <slepcpep.h>

 33: int main(int argc,char **argv)
 34: {
 35:   Mat            M,C,K,A[3];      /* problem matrices */
 36:   PEP            pep;             /* polynomial eigenproblem solver context */
 38:   PetscInt       n=30,Istart,Iend,i,nev;
 39:   PetscScalar    mu=1.0,tau=10.0,kappa=5.0;
 40:   PetscBool      terse=PETSC_FALSE;

 42:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 44:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 45:   PetscOptionsGetScalar(NULL,NULL,"-mu",&mu,NULL);
 46:   PetscOptionsGetScalar(NULL,NULL,"-tau",&tau,NULL);
 47:   PetscOptionsGetScalar(NULL,NULL,"-kappa",&kappa,NULL);

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

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

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

 70:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 71:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 73:   /* C is a tridiagonal */
 74:   MatCreate(PETSC_COMM_WORLD,&C);
 75:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 76:   MatSetFromOptions(C);
 77:   MatSetUp(C);

 79:   MatGetOwnershipRange(C,&Istart,&Iend);
 80:   for (i=Istart;i<Iend;i++) {
 81:     if (i>0) {
 82:       MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
 83:     }
 84:     MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
 85:     if (i<n-1) {
 86:       MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
 87:     }
 88:   }

 90:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 91:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 93:   /* M is a diagonal matrix */
 94:   MatCreate(PETSC_COMM_WORLD,&M);
 95:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
 96:   MatSetFromOptions(M);
 97:   MatSetUp(M);
 98:   MatGetOwnershipRange(M,&Istart,&Iend);
 99:   for (i=Istart;i<Iend;i++) {
100:     MatSetValue(M,i,i,mu,INSERT_VALUES);
101:   }
102:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
103:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

105:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106:                 Create the eigensolver and set various options
107:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

109:   PEPCreate(PETSC_COMM_WORLD,&pep);
110:   A[0] = K; A[1] = C; A[2] = M;
111:   PEPSetOperators(pep,3,A);
112:   PEPSetProblemType(pep,PEP_GENERAL);
113:   PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT);
114:   PEPSetFromOptions(pep);

116:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117:                       Solve the eigensystem
118:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

120:   PEPSolve(pep);
121:   PEPGetDimensions(pep,&nev,NULL,NULL);
122:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

124:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125:                       Display solution of first solve            
126:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
128:   if (terse) {
129:     PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
130:   } else {
131:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
132:     PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
133:     PEPErrorView(pep,PEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
134:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
135:   }
136:   MatDestroy(&M);
137:   MatDestroy(&C);
138:   MatDestroy(&K);

140:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141:      Compute the eigensystem, (k^2*M+k*C+K)x=0 for bigger n
142:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143:   
144:   n *= 2;
145:   /* K is a tridiagonal */
146:   MatCreate(PETSC_COMM_WORLD,&K);
147:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
148:   MatSetFromOptions(K);
149:   MatSetUp(K);

151:   MatGetOwnershipRange(K,&Istart,&Iend);
152:   for (i=Istart;i<Iend;i++) {
153:     if (i>0) {
154:       MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
155:     }
156:     MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
157:     if (i<n-1) {
158:       MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
159:     }
160:   }

162:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
163:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

165:   /* C is a tridiagonal */
166:   MatCreate(PETSC_COMM_WORLD,&C);
167:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
168:   MatSetFromOptions(C);
169:   MatSetUp(C);

171:   MatGetOwnershipRange(C,&Istart,&Iend);
172:   for (i=Istart;i<Iend;i++) {
173:     if (i>0) {
174:       MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
175:     }
176:     MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
177:     if (i<n-1) {
178:       MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
179:     }
180:   }

182:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
183:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

185:   /* M is a diagonal matrix */
186:   MatCreate(PETSC_COMM_WORLD,&M);
187:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
188:   MatSetFromOptions(M);
189:   MatSetUp(M);
190:   MatGetOwnershipRange(M,&Istart,&Iend);
191:   for (i=Istart;i<Iend;i++) {
192:     MatSetValue(M,i,i,mu,INSERT_VALUES);
193:   }
194:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
195:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

197:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198:        Solve again, calling PEPReset() since matrix size has changed
199:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200:   PEPReset(pep);  /* if this is omitted, it will be called in PEPSetOperators() */
201:   A[0] = K; A[1] = C; A[2] = M;
202:   PEPSetOperators(pep,3,A);
203:   PEPSolve(pep);

205:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206:                     Display solution and clean up
207:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
208:   if (terse) {
209:     PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
210:   } else {
211:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
212:     PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
213:     PEPErrorView(pep,PEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
214:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
215:   }
216:   PEPDestroy(&pep);
217:   MatDestroy(&M);
218:   MatDestroy(&C);
219:   MatDestroy(&K);
220:   SlepcFinalize();
221:   return ierr;
222: }

224: /*TEST

226:    test:
227:       suffix: 1
228:       args: -pep_type {{toar qarnoldi linear}} -pep_nev 4 -terse
229:       requires: !single

231:    test:
232:       suffix: 2
233:       args: -pep_type stoar -pep_hermitian -pep_nev 4 -terse
234:       requires: !single

236: TEST*/