Actual source code: test13.c

slepc-3.16.0 2021-09-30
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[] = "Solve a quadratic problem with CISS.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 14:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";

 16: #include <slepcpep.h>

 18: int main(int argc,char **argv)
 19: {
 20:   Mat               M,C,K,A[3];
 21:   PEP               pep;
 22:   RG                rg;
 23:   KSP               *ksp;
 24:   PC                pc;
 25:   PEPCISSExtraction ext;
 26:   PetscInt          N,n=10,m,Istart,Iend,II,i,j,nsolve;
 27:   PetscBool         flg;
 28:   PetscErrorCode    ierr;

 30:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 31:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 32:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flg);
 33:   if (!flg) m=n;
 34:   N = n*m;
 35:   PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%D (%Dx%D grid)\n\n",N,n,m);

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

 41:   /* K is the 2-D Laplacian */
 42:   MatCreate(PETSC_COMM_WORLD,&K);
 43:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N);
 44:   MatSetFromOptions(K);
 45:   MatSetUp(K);
 46:   MatGetOwnershipRange(K,&Istart,&Iend);
 47:   for (II=Istart;II<Iend;II++) {
 48:     i = II/n; j = II-i*n;
 49:     if (i>0) { MatSetValue(K,II,II-n,-1.0,INSERT_VALUES); }
 50:     if (i<m-1) { MatSetValue(K,II,II+n,-1.0,INSERT_VALUES); }
 51:     if (j>0) { MatSetValue(K,II,II-1,-1.0,INSERT_VALUES); }
 52:     if (j<n-1) { MatSetValue(K,II,II+1,-1.0,INSERT_VALUES); }
 53:     MatSetValue(K,II,II,4.0,INSERT_VALUES);
 54:   }
 55:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 56:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 58:   /* C is the 1-D Laplacian on horizontal lines */
 59:   MatCreate(PETSC_COMM_WORLD,&C);
 60:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
 61:   MatSetFromOptions(C);
 62:   MatSetUp(C);
 63:   MatGetOwnershipRange(C,&Istart,&Iend);
 64:   for (II=Istart;II<Iend;II++) {
 65:     i = II/n; j = II-i*n;
 66:     if (j>0) { MatSetValue(C,II,II-1,-1.0,INSERT_VALUES); }
 67:     if (j<n-1) { MatSetValue(C,II,II+1,-1.0,INSERT_VALUES); }
 68:     MatSetValue(C,II,II,2.0,INSERT_VALUES);
 69:   }
 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 (II=Istart;II<Iend;II++) {
 80:     MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES);
 81:   }
 82:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
 83:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

 85:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 86:              Create the eigensolver and solve the eigensystem
 87:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 89:   PEPCreate(PETSC_COMM_WORLD,&pep);
 90:   A[0] = K; A[1] = C; A[2] = M;
 91:   PEPSetOperators(pep,3,A);
 92:   PEPSetProblemType(pep,PEP_GENERAL);

 94:   /* customize polynomial eigensolver; set runtime options */
 95:   PEPSetType(pep,PEPCISS);
 96:   PEPGetRG(pep,&rg);
 97:   RGSetType(rg,RGELLIPSE);
 98:   RGEllipseSetParameters(rg,PetscCMPLX(-0.1,0.3),0.1,0.25);
 99:   PEPCISSSetSizes(pep,24,PETSC_DEFAULT,PETSC_DEFAULT,1,PETSC_DEFAULT,PETSC_TRUE);
100:   PEPCISSGetKSPs(pep,&nsolve,&ksp);
101:   for (i=0;i<nsolve;i++) {
102:     KSPSetTolerances(ksp[i],1e-12,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
103:     KSPSetType(ksp[i],KSPPREONLY);
104:     KSPGetPC(ksp[i],&pc);
105:     PCSetType(pc,PCLU);
106:   }
107:   PEPSetFromOptions(pep);

109:   /* solve */
110:   PetscObjectTypeCompare((PetscObject)pep,PEPCISS,&flg);
111:   if (flg) {
112:     PEPCISSGetExtraction(pep,&ext);
113:     PetscPrintf(PETSC_COMM_WORLD," Running CISS with %D KSP solvers (%s extraction)\n",nsolve,PEPCISSExtractions[ext]);
114:   }
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:   SlepcFinalize();
127:   return ierr;
128: }

130: /*TEST

132:    build:
133:       requires: complex

135:    test:
136:       suffix: 1
137:       requires: complex

139: TEST*/