Actual source code: test7.c

slepc-3.13.2 2020-05-12
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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,tau=10,kappa=5,int0,int1,*shifts;

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

 31:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 32:   PetscOptionsGetBool(NULL,NULL,"-showinertia",&showinertia,NULL);
 33:   PetscPrintf(PETSC_COMM_WORLD,"\nSpectrum slicing on PEP, n=%D\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:   MatCreate(PETSC_COMM_WORLD,&K);
 41:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
 42:   MatSetFromOptions(K);
 43:   MatSetUp(K);

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

 56:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 57:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

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

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

 76:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 77:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

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

 91:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 92:                 Create the eigensolver and solve the problem
 93:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 95:   PEPCreate(PETSC_COMM_WORLD,&pep);
 96:   A[0] = K; A[1] = C; A[2] = M;
 97:   PEPSetOperators(pep,3,A);
 98:   PEPSetProblemType(pep,PEP_HYPERBOLIC);
 99:   PEPSetType(pep,PEPSTOAR);

101:   /*
102:      Set interval and other settings for spectrum slicing
103:   */
104:   int0 = -11.3;
105:   int1 = -9.5;
106:   PEPSetInterval(pep,int0,int1);
107:   PEPSetWhichEigenpairs(pep,PEP_ALL);
108:   PEPGetST(pep,&st);
109:   STSetType(st,STSINVERT);
110:   STGetKSP(st,&ksp);
111:   KSPSetType(ksp,KSPPREONLY);
112:   KSPGetPC(ksp,&pc);
113:   PCSetType(pc,PCCHOLESKY);

115:   /*
116:      Test interface functions of STOAR solver
117:   */
118:   PEPSTOARGetDetectZeros(pep,&detect);
119:   PetscPrintf(PETSC_COMM_WORLD," Detect zeros before changing = %d",(int)detect);
120:   PEPSTOARSetDetectZeros(pep,PETSC_TRUE);
121:   PEPSTOARGetDetectZeros(pep,&detect);
122:   PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)detect);

124:   PEPSTOARGetLocking(pep,&lock);
125:   PetscPrintf(PETSC_COMM_WORLD," Locking flag before changing = %d",(int)lock);
126:   PEPSTOARSetLocking(pep,PETSC_TRUE);
127:   PEPSTOARGetLocking(pep,&lock);
128:   PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)lock);

130:   PEPSTOARGetCheckEigenvalueType(pep,&checket);
131:   PetscPrintf(PETSC_COMM_WORLD," Check eigenvalue type flag before changing = %d",(int)checket);
132:   PEPSTOARSetCheckEigenvalueType(pep,PETSC_FALSE);
133:   PEPSTOARGetCheckEigenvalueType(pep,&checket);
134:   PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)checket);

136:   PEPSTOARGetDimensions(pep,&nev,&ncv,&mpd);
137:   PetscPrintf(PETSC_COMM_WORLD," Sub-solve dimensions before changing = [%D,%D,%D]",nev,ncv,mpd);
138:   PEPSTOARSetDimensions(pep,30,60,60);
139:   PEPSTOARGetDimensions(pep,&nev,&ncv,&mpd);
140:   PetscPrintf(PETSC_COMM_WORLD," ... changed to [%D,%D,%D]\n",nev,ncv,mpd);

142:   PEPSetFromOptions(pep);

144:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145:              Compute all eigenvalues in interval and display info
146:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

148:   PEPSetUp(pep);
149:   PEPSTOARGetInertias(pep,&ns,&shifts,&inertias);
150:   PetscPrintf(PETSC_COMM_WORLD," Inertias (after setup):\n");
151:   for (i=0;i<ns;i++) {
152:     PetscPrintf(PETSC_COMM_WORLD," .. %g (%D)\n",(double)shifts[i],inertias[i]);
153:   }
154:   PetscFree(shifts);
155:   PetscFree(inertias);

157:   PEPSolve(pep);
158:   PEPGetDimensions(pep,&nev,NULL,NULL);
159:   PEPGetInterval(pep,&int0,&int1);
160:   PetscPrintf(PETSC_COMM_WORLD," Found %D eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1);

162:   if (showinertia) {
163:     PEPSTOARGetInertias(pep,&ns,&shifts,&inertias);
164:     PetscPrintf(PETSC_COMM_WORLD," Used %D shifts (inertia):\n",ns);
165:     for (i=0;i<ns;i++) {
166:       PetscPrintf(PETSC_COMM_WORLD," .. %g (%D)\n",(double)shifts[i],inertias[i]);
167:     }
168:     PetscFree(shifts);
169:     PetscFree(inertias);
170:   }

172:   PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);

174:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175:                     Clean up
176:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
177:   PEPDestroy(&pep);
178:   MatDestroy(&M);
179:   MatDestroy(&C);
180:   MatDestroy(&K);
181:   SlepcFinalize();
182:   return ierr;
183: }

185: /*TEST

187:    test:
188:       requires: !single
189:       args: -showinertia 0

191: TEST*/