Actual source code: test6.c

slepc-3.14.2 2021-02-01
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 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;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

113:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114:                              Clean up
115:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116:   STDestroy(&st);
117:   MatDestroy(&A);
118:   MatDestroy(&P);
119:   VecDestroy(&v);
120:   VecDestroy(&w);
121:   SlepcFinalize();
122:   return ierr;
123: }

125: /*TEST

127:    test:
128:       suffix: 1
129:       args: -st_matmode {{copy inplace}}
130:       requires: !single

132:    test:
133:       suffix: 2
134:       args: -st_matmode shell

136: TEST*/