Actual source code: test7.c

slepc-3.17.1 2022-04-11
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:   SlepcInitialize(&argc,&argv,(char*)0,help);

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

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

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

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

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

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

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

 67:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 68:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

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

 80:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 81:                 Create the eigensolver and solve the problem
 82:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

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

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

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

131:   PEPSetFromOptions(pep);

133:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134:              Compute all eigenvalues in interval and display info
135:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

157:   PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);

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

170: /*TEST

172:    test:
173:       requires: !single
174:       args: -showinertia 0

176: TEST*/