Actual source code: test6.c

slepc-3.21.1 2024-04-26
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:   PetscFunctionBeginUser;
 42:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

 44:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 45:   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL));
 46:   PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
 47:   PetscCall(PetscOptionsGetReal(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:   PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
 55:   PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
 56:   PetscCall(MatSetFromOptions(K));

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

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

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

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

 80:   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
 81:   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));

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

 92:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 93:                 Create the eigensolver and set various options
 94:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

103:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104:                       Solve the eigensystem
105:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

130:   n *= 2;
131:   /* K is a tridiagonal */
132:   PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
133:   PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
134:   PetscCall(MatSetFromOptions(K));

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

143:   PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
144:   PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));

146:   /* C is a tridiagonal */
147:   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
148:   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n));
149:   PetscCall(MatSetFromOptions(C));

151:   PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
152:   for (i=Istart;i<Iend;i++) {
153:     if (i>0) PetscCall(MatSetValue(C,i,i-1,-tau,INSERT_VALUES));
154:     PetscCall(MatSetValue(C,i,i,tau*3.0,INSERT_VALUES));
155:     if (i<n-1) PetscCall(MatSetValue(C,i,i+1,-tau,INSERT_VALUES));
156:   }

158:   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
159:   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));

161:   /* M is a diagonal matrix */
162:   PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
163:   PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
164:   PetscCall(MatSetFromOptions(M));
165:   PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
166:   for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(M,i,i,mu,INSERT_VALUES));
167:   PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
168:   PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));

170:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171:        Solve again, calling PEPReset() since matrix size has changed
172:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173:   /* PetscCall(PEPReset(pep)); */  /* not required, will be called in PEPSetOperators() */
174:   A[0] = K; A[1] = C; A[2] = M;
175:   PetscCall(PEPSetOperators(pep,3,A));
176:   PetscCall(PEPSolve(pep));

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

196: /*TEST

198:    test:
199:       suffix: 1
200:       args: -pep_type {{toar qarnoldi linear}} -pep_nev 4 -terse
201:       requires: double

203:    test:
204:       suffix: 2
205:       args: -pep_type stoar -pep_hermitian -pep_nev 4 -terse
206:       requires: !single

208: TEST*/