Actual source code: test12.c

slepc-3.18.3 2023-03-24
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: */

 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);

 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=%" PetscInt_FMT " 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) MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
 50:     MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
 51:     if (i<n-1) MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
 52:   }

 54:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 55:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 57:   /* C is a tridiagonal */
 58:   MatCreate(PETSC_COMM_WORLD,&C);
 59:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 60:   MatSetFromOptions(C);
 61:   MatSetUp(C);

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

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

 73:   /* M is a diagonal matrix */
 74:   MatCreate(PETSC_COMM_WORLD,&M);
 75:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
 76:   MatSetFromOptions(M);
 77:   MatSetUp(M);
 78:   MatGetOwnershipRange(M,&Istart,&Iend);
 79:   for (i=Istart;i<Iend;i++) MatSetValue(M,i,i,mu,INSERT_VALUES);
 80:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
 81:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

 83:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 84:                     Create a region for filtering
 85:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 87:   RGCreate(PETSC_COMM_WORLD,&rg);
 88:   RGSetType(rg,RGINTERVAL);
 89: #if defined(PETSC_USE_COMPLEX)
 90:   RGIntervalSetEndpoints(rg,-47.0,-35.0,-0.001,0.001);
 91: #else
 92:   RGIntervalSetEndpoints(rg,-47.0,-35.0,-0.0,0.0);
 93: #endif

 95:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 96:                 Create the eigensolver and solve the problem
 97:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 99:   PEPCreate(PETSC_COMM_WORLD,&pep);
100:   PEPSetRG(pep,rg);
101:   A[0] = K; A[1] = C; A[2] = M;
102:   PEPSetOperators(pep,3,A);
103:   PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT);
104:   PEPSetFromOptions(pep);
105:   PEPSolve(pep);

107:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108:                     Display solution and clean up
109:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

111:   PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
112:   PEPDestroy(&pep);
113:   MatDestroy(&M);
114:   MatDestroy(&C);
115:   MatDestroy(&K);
116:   RGDestroy(&rg);
117:   SlepcFinalize();
118:   return 0;
119: }

121: /*TEST

123:    test:
124:       args: -pep_nev 8 -pep_type {{toar linear qarnoldi}}
125:       suffix: 1
126:       requires: !single
127:       output_file: output/test12_1.out

129: TEST*/