Actual source code: test6.c

slepc-3.17.0 2022-03-31
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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 */
 37:   PetscInt       n=30,Istart,Iend,i,nev;
 38:   PetscReal      mu=1.0,tau=10.0,kappa=5.0;
 39:   PetscBool      terse=PETSC_FALSE;

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

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

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

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

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

 65:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 66:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 68:   /* C is a tridiagonal */
 69:   MatCreate(PETSC_COMM_WORLD,&C);
 70:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 71:   MatSetFromOptions(C);
 72:   MatSetUp(C);

 74:   MatGetOwnershipRange(C,&Istart,&Iend);
 75:   for (i=Istart;i<Iend;i++) {
 76:     if (i>0) MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
 77:     MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
 78:     if (i<n-1) MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
 79:   }

 81:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 82:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 84:   /* M is a diagonal matrix */
 85:   MatCreate(PETSC_COMM_WORLD,&M);
 86:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
 87:   MatSetFromOptions(M);
 88:   MatSetUp(M);
 89:   MatGetOwnershipRange(M,&Istart,&Iend);
 90:   for (i=Istart;i<Iend;i++) MatSetValue(M,i,i,mu,INSERT_VALUES);
 91:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
 92:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 95:                 Create the eigensolver and set various options
 96:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 98:   PEPCreate(PETSC_COMM_WORLD,&pep);
 99:   A[0] = K; A[1] = C; A[2] = M;
100:   PEPSetOperators(pep,3,A);
101:   PEPSetProblemType(pep,PEP_GENERAL);
102:   PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT);
103:   PEPSetFromOptions(pep);

105:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106:                       Solve the eigensystem
107:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

109:   PEPSolve(pep);
110:   PEPGetDimensions(pep,&nev,NULL,NULL);
111:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);

113:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114:                       Display solution of first solve
115:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
117:   if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
118:   else {
119:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
120:     PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
121:     PEPErrorView(pep,PEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
122:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
123:   }
124:   MatDestroy(&M);
125:   MatDestroy(&C);
126:   MatDestroy(&K);

128:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129:      Compute the eigensystem, (k^2*M+k*C+K)x=0 for bigger n
130:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

132:   n *= 2;
133:   /* K is a tridiagonal */
134:   MatCreate(PETSC_COMM_WORLD,&K);
135:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
136:   MatSetFromOptions(K);
137:   MatSetUp(K);

139:   MatGetOwnershipRange(K,&Istart,&Iend);
140:   for (i=Istart;i<Iend;i++) {
141:     if (i>0) MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
142:     MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
143:     if (i<n-1) MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
144:   }

146:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
147:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

149:   /* C is a tridiagonal */
150:   MatCreate(PETSC_COMM_WORLD,&C);
151:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
152:   MatSetFromOptions(C);
153:   MatSetUp(C);

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

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

165:   /* M is a diagonal matrix */
166:   MatCreate(PETSC_COMM_WORLD,&M);
167:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
168:   MatSetFromOptions(M);
169:   MatSetUp(M);
170:   MatGetOwnershipRange(M,&Istart,&Iend);
171:   for (i=Istart;i<Iend;i++) MatSetValue(M,i,i,mu,INSERT_VALUES);
172:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
173:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

175:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176:        Solve again, calling PEPReset() since matrix size has changed
177:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178:   /*PEPReset(pep);*/  /* not required, will be called in PEPSetOperators() */
179:   A[0] = K; A[1] = C; A[2] = M;
180:   PEPSetOperators(pep,3,A);
181:   PEPSolve(pep);

183:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184:                     Display solution and clean up
185:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
186:   if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
187:   else {
188:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
189:     PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
190:     PEPErrorView(pep,PEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
191:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
192:   }
193:   PEPDestroy(&pep);
194:   MatDestroy(&M);
195:   MatDestroy(&C);
196:   MatDestroy(&K);
197:   SlepcFinalize();
198:   return 0;
199: }

201: /*TEST

203:    test:
204:       suffix: 1
205:       args: -pep_type {{toar qarnoldi linear}} -pep_nev 4 -terse
206:       requires: !single

208:    test:
209:       suffix: 2
210:       args: -pep_type stoar -pep_hermitian -pep_nev 4 -terse
211:       requires: !single

213: TEST*/