Actual source code: test7.c

slepc-3.17.2 2022-08-09
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 ST with one matrix and split preconditioner.\n\n";

 13: #include <slepcst.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat            A,Pa,Pmat,mat[1];
 18:   ST             st;
 19:   KSP            ksp;
 20:   PC             pc;
 21:   Vec            v,w;
 22:   STType         type;
 23:   PetscBool      flg;
 24:   PetscScalar    sigma;
 25:   PetscInt       n=10,i,Istart,Iend;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);
 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian, n=%" PetscInt_FMT "\n\n",n);

 31:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 32:      Compute the operator matrix for the 1-D Laplacian
 33:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 35:   MatCreate(PETSC_COMM_WORLD,&A);
 36:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 37:   MatSetFromOptions(A);
 38:   MatSetUp(A);

 40:   MatGetOwnershipRange(A,&Istart,&Iend);
 41:   for (i=Istart;i<Iend;i++) {
 42:     if (i>0) MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
 43:     if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
 44:     MatSetValue(A,i,i,2.0,INSERT_VALUES);
 45:   }
 46:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 47:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 48:   MatCreateVecs(A,&v,&w);
 49:   VecSet(v,1.0);

 51:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 52:      Compute the split preconditioner matrix (one diagonal)
 53:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 55:   MatCreate(PETSC_COMM_WORLD,&Pa);
 56:   MatSetSizes(Pa,PETSC_DECIDE,PETSC_DECIDE,n,n);
 57:   MatSetFromOptions(Pa);
 58:   MatSetUp(Pa);

 60:   MatGetOwnershipRange(Pa,&Istart,&Iend);
 61:   for (i=Istart;i<Iend;i++) MatSetValue(Pa,i,i,2.0,INSERT_VALUES);
 62:   MatAssemblyBegin(Pa,MAT_FINAL_ASSEMBLY);
 63:   MatAssemblyEnd(Pa,MAT_FINAL_ASSEMBLY);

 65:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 66:                 Create the spectral transformation object
 67:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 69:   STCreate(PETSC_COMM_WORLD,&st);
 70:   mat[0] = A;
 71:   STSetMatrices(st,1,mat);
 72:   mat[0] = Pa;
 73:   STSetSplitPreconditioner(st,1,mat,SAME_NONZERO_PATTERN);
 74:   STSetTransform(st,PETSC_TRUE);
 75:   STSetFromOptions(st);
 76:   STCayleySetAntishift(st,-0.3);   /* only relevant for cayley */

 78:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 79:                Form the preconditioner matrix and print it
 80:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 82:   PetscObjectTypeCompareAny((PetscObject)st,&flg,STSINVERT,STCAYLEY,"");
 83:   if (flg) {
 84:     STGetKSP(st,&ksp);
 85:     KSPGetPC(ksp,&pc);
 86:     STGetOperator(st,NULL);
 87:     PCGetOperators(pc,NULL,&Pmat);
 88:     MatView(Pmat,NULL);
 89:   }

 91:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 92:                     Apply the operator
 93:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 95:   /* sigma=0.0 */
 96:   STSetUp(st);
 97:   STGetType(st,&type);
 98:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
 99:   STApply(st,v,w);
100:   VecView(w,NULL);

102:   /* sigma=0.1 */
103:   sigma = 0.1;
104:   STSetShift(st,sigma);
105:   STGetShift(st,&sigma);
106:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
107:   if (flg) {
108:     STGetOperator(st,NULL);
109:     PCGetOperators(pc,NULL,&Pmat);
110:     MatView(Pmat,NULL);
111:   }
112:   STApply(st,v,w);
113:   VecView(w,NULL);

115:   STDestroy(&st);
116:   MatDestroy(&A);
117:   MatDestroy(&Pa);
118:   VecDestroy(&v);
119:   VecDestroy(&w);
120:   SlepcFinalize();
121:   return 0;
122: }

124: /*TEST

126:    test:
127:       suffix: 1
128:       args: -st_type {{cayley shift sinvert}separate output}
129:       requires: !single

131: TEST*/