Actual source code: ex16.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: static char help[] = "Quadratic eigenproblem for testing the QEP object.\n\n"
 23:   "The command line options are:\n"
 24:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 25:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";

 27: #include <slepcqep.h>

 31: int main(int argc,char **argv)
 32: {
 33:   Mat            M,C,K;           /* problem matrices */
 34:   QEP            qep;             /* quadratic eigenproblem solver context */
 35:   QEPType        type;
 36:   PetscInt       N,n=10,m,Istart,Iend,II,nev,i,j;
 37:   PetscBool      flag;

 40:   SlepcInitialize(&argc,&argv,(char*)0,help);

 42:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 43:   PetscOptionsGetInt(NULL,"-m",&m,&flag);
 44:   if (!flag) m=n;
 45:   N = n*m;
 46:   PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%D (%Dx%D grid)\n\n",N,n,m);

 48:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 49:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 50:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 52:   /* K is the 2-D Laplacian */
 53:   MatCreate(PETSC_COMM_WORLD,&K);
 54:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N);
 55:   MatSetFromOptions(K);
 56:   MatSetUp(K);

 58:   MatGetOwnershipRange(K,&Istart,&Iend);
 59:   for (II=Istart;II<Iend;II++) {
 60:     i = II/n; j = II-i*n;
 61:     if (i>0) { MatSetValue(K,II,II-n,-1.0,INSERT_VALUES); }
 62:     if (i<m-1) { MatSetValue(K,II,II+n,-1.0,INSERT_VALUES); }
 63:     if (j>0) { MatSetValue(K,II,II-1,-1.0,INSERT_VALUES); }
 64:     if (j<n-1) { MatSetValue(K,II,II+1,-1.0,INSERT_VALUES); }
 65:     MatSetValue(K,II,II,4.0,INSERT_VALUES);
 66:   }

 68:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 69:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 71:   /* C is the zero matrix */
 72:   MatCreate(PETSC_COMM_WORLD,&C);
 73:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
 74:   MatSetFromOptions(C);
 75:   MatSetUp(C);
 76:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 77:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 79:   /* M is the identity matrix */
 80:   MatCreate(PETSC_COMM_WORLD,&M);
 81:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
 82:   MatSetFromOptions(M);
 83:   MatSetUp(M);
 84:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
 85:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
 86:   MatShift(M,1.0);

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:                 Create the eigensolver and set various options
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 92:   /*
 93:      Create eigensolver context
 94:   */
 95:   QEPCreate(PETSC_COMM_WORLD,&qep);

 97:   /*
 98:      Set matrices and problem type
 99:   */
100:   QEPSetOperators(qep,M,C,K);
101:   QEPSetProblemType(qep,QEP_HERMITIAN);

103:   /*
104:      Set solver parameters at runtime
105:   */
106:   QEPSetFromOptions(qep);

108:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109:                       Solve the eigensystem
110:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

112:   QEPSolve(qep);

114:   /*
115:      Optional: Get some information from the solver and display it
116:   */
117:   QEPGetType(qep,&type);
118:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
119:   QEPGetDimensions(qep,&nev,NULL,NULL);
120:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

122:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123:                     Display solution and clean up
124:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

126:   QEPPrintSolution(qep,NULL);
127:   QEPDestroy(&qep);
128:   MatDestroy(&M);
129:   MatDestroy(&C);
130:   MatDestroy(&K);
131:   SlepcFinalize();
132:   return 0;
133: }