Actual source code: test7.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: */

 11: static char help[] = "Test interface functions of spectrum-slicing STOAR.\n\n"
 12:   "This is based on ex38.c. The command line options are:\n"
 13:   "  -n <n> ... dimension of the matrices.\n\n";

 15: #include <slepcpep.h>

 17: int main(int argc,char **argv)
 18: {
 19:   Mat            M,C,K,A[3]; /* problem matrices */
 20:   PEP            pep;        /* polynomial eigenproblem solver context */
 21:   ST             st;         /* spectral transformation context */
 22:   KSP            ksp;
 23:   PC             pc;
 24:   PetscBool      showinertia=PETSC_TRUE,lock,detect,checket;
 25:   PetscInt       n=100,Istart,Iend,i,*inertias,ns,nev,ncv,mpd;
 26:   PetscReal      mu=1.0,tau=10.0,kappa=5.0,int0,int1,*shifts;

 28:   PetscFunctionBeginUser;
 29:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

 31:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 32:   PetscCall(PetscOptionsGetBool(NULL,NULL,"-showinertia",&showinertia,NULL));
 33:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSpectrum slicing on PEP, n=%" PetscInt_FMT "\n\n",n));

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

 39:   /* K is a tridiagonal */
 40:   PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
 41:   PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
 42:   PetscCall(MatSetFromOptions(K));

 44:   PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
 45:   for (i=Istart;i<Iend;i++) {
 46:     if (i>0) PetscCall(MatSetValue(K,i,i-1,-kappa,INSERT_VALUES));
 47:     PetscCall(MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES));
 48:     if (i<n-1) PetscCall(MatSetValue(K,i,i+1,-kappa,INSERT_VALUES));
 49:   }

 51:   PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
 52:   PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));

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

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

 66:   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
 67:   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));

 69:   /* M is a diagonal matrix */
 70:   PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
 71:   PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
 72:   PetscCall(MatSetFromOptions(M));
 73:   PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
 74:   for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(M,i,i,mu,INSERT_VALUES));
 75:   PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
 76:   PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));

 78:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 79:                 Create the eigensolver and solve the problem
 80:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 82:   PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
 83:   A[0] = K; A[1] = C; A[2] = M;
 84:   PetscCall(PEPSetOperators(pep,3,A));
 85:   PetscCall(PEPSetProblemType(pep,PEP_HYPERBOLIC));
 86:   PetscCall(PEPSetType(pep,PEPSTOAR));

 88:   /*
 89:      Set interval and other settings for spectrum slicing
 90:   */
 91:   int0 = -11.3;
 92:   int1 = -9.5;
 93:   PetscCall(PEPSetInterval(pep,int0,int1));
 94:   PetscCall(PEPSetWhichEigenpairs(pep,PEP_ALL));
 95:   PetscCall(PEPGetST(pep,&st));
 96:   PetscCall(STSetType(st,STSINVERT));
 97:   PetscCall(STGetKSP(st,&ksp));
 98:   PetscCall(KSPSetType(ksp,KSPPREONLY));
 99:   PetscCall(KSPGetPC(ksp,&pc));
100:   PetscCall(PCSetType(pc,PCCHOLESKY));

102:   /*
103:      Test interface functions of STOAR solver
104:   */
105:   PetscCall(PEPSTOARGetDetectZeros(pep,&detect));
106:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Detect zeros before changing = %d",(int)detect));
107:   PetscCall(PEPSTOARSetDetectZeros(pep,PETSC_TRUE));
108:   PetscCall(PEPSTOARGetDetectZeros(pep,&detect));
109:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)detect));

111:   PetscCall(PEPSTOARGetLocking(pep,&lock));
112:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Locking flag before changing = %d",(int)lock));
113:   PetscCall(PEPSTOARSetLocking(pep,PETSC_TRUE));
114:   PetscCall(PEPSTOARGetLocking(pep,&lock));
115:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)lock));

117:   PetscCall(PEPSTOARGetCheckEigenvalueType(pep,&checket));
118:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Check eigenvalue type flag before changing = %d",(int)checket));
119:   PetscCall(PEPSTOARSetCheckEigenvalueType(pep,PETSC_FALSE));
120:   PetscCall(PEPSTOARGetCheckEigenvalueType(pep,&checket));
121:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)checket));

123:   PetscCall(PEPSTOARGetDimensions(pep,&nev,&ncv,&mpd));
124:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Sub-solve dimensions before changing = [%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]",nev,ncv,mpd));
125:   PetscCall(PEPSTOARSetDimensions(pep,30,60,60));
126:   PetscCall(PEPSTOARGetDimensions(pep,&nev,&ncv,&mpd));
127:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to [%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]\n",nev,ncv,mpd));

129:   PetscCall(PEPSetFromOptions(pep));

131:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132:              Compute all eigenvalues in interval and display info
133:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

135:   PetscCall(PEPSetUp(pep));
136:   PetscCall(PEPSTOARGetInertias(pep,&ns,&shifts,&inertias));
137:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Inertias (after setup):\n"));
138:   for (i=0;i<ns;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," .. %g (%" PetscInt_FMT ")\n",(double)shifts[i],inertias[i]));
139:   PetscCall(PetscFree(shifts));
140:   PetscCall(PetscFree(inertias));

142:   PetscCall(PEPSolve(pep));
143:   PetscCall(PEPGetDimensions(pep,&nev,NULL,NULL));
144:   PetscCall(PEPGetInterval(pep,&int0,&int1));
145:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Found %" PetscInt_FMT " eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1));

147:   if (showinertia) {
148:     PetscCall(PEPSTOARGetInertias(pep,&ns,&shifts,&inertias));
149:     PetscCall(PetscPrintf(PETSC_COMM_WORLD," Used %" PetscInt_FMT " shifts (inertia):\n",ns));
150:     for (i=0;i<ns;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," .. %g (%" PetscInt_FMT ")\n",(double)shifts[i],inertias[i]));
151:     PetscCall(PetscFree(shifts));
152:     PetscCall(PetscFree(inertias));
153:   }

155:   PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));

157:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158:                     Clean up
159:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
160:   PetscCall(PEPDestroy(&pep));
161:   PetscCall(MatDestroy(&M));
162:   PetscCall(MatDestroy(&C));
163:   PetscCall(MatDestroy(&K));
164:   PetscCall(SlepcFinalize());
165:   return 0;
166: }

168: /*TEST

170:    test:
171:       requires: !single
172:       args: -showinertia 0

174: TEST*/