Actual source code: test12.c

slepc-3.16.3 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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: */

 11: static char help[] = "Illustrates region filtering in PEP (based on spring).\n"
 12:   "The command line options are:\n"
 13:   "  -n <n> ... number of grid subdivisions.\n"
 14:   "  -mu <value> ... mass (default 1).\n"
 15:   "  -tau <value> ... damping constant of the dampers (default 10).\n"
 16:   "  -kappa <value> ... damping constant of the springs (default 5).\n\n";

 18: #include <slepcpep.h>

 20: int main(int argc,char **argv)
 21: {
 22:   Mat            M,C,K,A[3];      /* problem matrices */
 23:   PEP            pep;             /* polynomial eigenproblem solver context */
 24:   RG             rg;
 25:   PetscInt       n=30,Istart,Iend,i;
 26:   PetscReal      mu=1.0,tau=10.0,kappa=5.0;

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

 31:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 32:   PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL);
 33:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 34:   PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL);
 35:   PetscPrintf(PETSC_COMM_WORLD,"\nDamped mass-spring system, n=%D mu=%g tau=%g kappa=%g\n\n",n,(double)mu,(double)tau,(double)kappa);

 37:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 38:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 39:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 41:   /* K is a tridiagonal */
 42:   MatCreate(PETSC_COMM_WORLD,&K);
 43:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
 44:   MatSetFromOptions(K);
 45:   MatSetUp(K);

 47:   MatGetOwnershipRange(K,&Istart,&Iend);
 48:   for (i=Istart;i<Iend;i++) {
 49:     if (i>0) {
 50:       MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
 51:     }
 52:     MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
 53:     if (i<n-1) {
 54:       MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
 55:     }
 56:   }

 58:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 59:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 61:   /* C is a tridiagonal */
 62:   MatCreate(PETSC_COMM_WORLD,&C);
 63:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 64:   MatSetFromOptions(C);
 65:   MatSetUp(C);

 67:   MatGetOwnershipRange(C,&Istart,&Iend);
 68:   for (i=Istart;i<Iend;i++) {
 69:     if (i>0) {
 70:       MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
 71:     }
 72:     MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
 73:     if (i<n-1) {
 74:       MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
 75:     }
 76:   }

 78:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 79:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

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

 93:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 94:                     Create a region for filtering
 95:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 97:   RGCreate(PETSC_COMM_WORLD,&rg);
 98:   RGSetType(rg,RGINTERVAL);
 99: #if defined(PETSC_USE_COMPLEX)
100:   RGIntervalSetEndpoints(rg,-47.0,-35.0,-0.001,0.001);
101: #else
102:   RGIntervalSetEndpoints(rg,-47.0,-35.0,-0.0,0.0);
103: #endif

105:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106:                 Create the eigensolver and solve the problem
107:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

117:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118:                     Display solution and clean up
119:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

121:   PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
122:   PEPDestroy(&pep);
123:   MatDestroy(&M);
124:   MatDestroy(&C);
125:   MatDestroy(&K);
126:   RGDestroy(&rg);
127:   SlepcFinalize();
128:   return ierr;
129: }

131: /*TEST

133:    test:
134:       args: -pep_nev 8 -pep_type {{toar linear qarnoldi}}
135:       suffix: 1
136:       requires: !single
137:       output_file: output/test12_1.out

139: TEST*/