Actual source code: test6.c

slepc-3.18.1 2022-11-02
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 STPRECOND operations.\n\n";

 13: #include <slepcst.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat            A,P,mat[1];
 18:   ST             st;
 19:   KSP            ksp;
 20:   Vec            v,w;
 21:   PetscScalar    sigma;
 22:   PetscInt       n=10,i,Istart,Iend;
 23:   STMatMode      matmode;

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

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

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

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

 50:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 51:                 Create the spectral transformation object
 52:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 54:   STCreate(PETSC_COMM_WORLD,&st);
 55:   mat[0] = A;
 56:   STSetMatrices(st,1,mat);
 57:   STSetType(st,STPRECOND);
 58:   STGetKSP(st,&ksp);
 59:   KSPSetType(ksp,KSPPREONLY);
 60:   STSetFromOptions(st);

 62:   /* set up */
 63:   /* - the transform flag is necessary so that A-sigma*I is built explicitly */
 64:   STSetTransform(st,PETSC_TRUE);
 65:   /* - the ksphasmat flag is necessary when using STApply(), otherwise can only use PCApply() */
 66:   STPrecondSetKSPHasMat(st,PETSC_TRUE);
 67:   /* no need to call STSetUp() explicitly */
 68:   STSetUp(st);

 70:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 71:                        Apply the preconditioner
 72:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 74:   /* default shift */
 75:   PetscPrintf(PETSC_COMM_WORLD,"With default shift\n");
 76:   STApply(st,v,w);
 77:   VecView(w,NULL);

 79:   /* change shift */
 80:   sigma = 0.1;
 81:   STSetShift(st,sigma);
 82:   STGetShift(st,&sigma);
 83:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
 84:   STApply(st,v,w);
 85:   VecView(w,NULL);
 86:   STPostSolve(st);   /* undo changes if inplace */

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:                  Test a user-provided preconditioner matrix
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 92:   MatCreate(PETSC_COMM_WORLD,&P);
 93:   MatSetSizes(P,PETSC_DECIDE,PETSC_DECIDE,n,n);
 94:   MatSetFromOptions(P);
 95:   MatSetUp(P);

 97:   MatGetOwnershipRange(P,&Istart,&Iend);
 98:   for (i=Istart;i<Iend;i++) MatSetValue(P,i,i,2.0,INSERT_VALUES);
 99:   if (Istart==0) {
100:     MatSetValue(P,1,0,-1.0,INSERT_VALUES);
101:     MatSetValue(P,0,1,-1.0,INSERT_VALUES);
102:   }
103:   MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);
104:   MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY);

106:   /* apply new preconditioner */
107:   STSetPreconditionerMat(st,P);
108:   PetscPrintf(PETSC_COMM_WORLD,"With user-provided matrix\n");
109:   STApply(st,v,w);
110:   VecView(w,NULL);

112:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113:             Test a user-provided preconditioner in split form
114:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

116:   STGetMatMode(st,&matmode);
117:   if (matmode==ST_MATMODE_COPY) {
118:     STSetPreconditionerMat(st,NULL);
119:     mat[0] = P;
120:     STSetSplitPreconditioner(st,1,mat,SAME_NONZERO_PATTERN);

122:     /* apply new preconditioner */
123:     PetscPrintf(PETSC_COMM_WORLD,"With split preconditioner\n");
124:     STApply(st,v,w);
125:     VecView(w,NULL);
126:   }

128:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129:                              Clean up
130:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131:   STDestroy(&st);
132:   MatDestroy(&A);
133:   MatDestroy(&P);
134:   VecDestroy(&v);
135:   VecDestroy(&w);
136:   SlepcFinalize();
137:   return 0;
138: }

140: /*TEST

142:    test:
143:       suffix: 1
144:       args: -st_matmode {{copy inplace shell}separate output}
145:       requires: !single

147: TEST*/