Actual source code: test8.c

slepc-3.21.2 2024-09-25
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 two matrices and split preconditioner.\n\n";

 13: #include <slepcst.h>

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

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

 31:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 32:      Compute the operator matrices (1-D Laplacian and diagonal)
 33:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

 39:   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
 40:   PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
 41:   PetscCall(MatSetFromOptions(B));

 43:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 44:   for (i=Istart;i<Iend;i++) {
 45:     PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
 46:     if (i>0) {
 47:       PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
 48:       PetscCall(MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES));
 49:     } else PetscCall(MatSetValue(B,i,i,-1.0,INSERT_VALUES));
 50:     if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
 51:   }
 52:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 53:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
 54:   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
 55:   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
 56:   PetscCall(MatCreateVecs(A,&v,&w));
 57:   PetscCall(VecSet(v,1.0));

 59:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 60:      Compute the split preconditioner matrices (two diagonals)
 61:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 63:   PetscCall(MatCreate(PETSC_COMM_WORLD,&Pa));
 64:   PetscCall(MatSetSizes(Pa,PETSC_DECIDE,PETSC_DECIDE,n,n));
 65:   PetscCall(MatSetFromOptions(Pa));

 67:   PetscCall(MatCreate(PETSC_COMM_WORLD,&Pb));
 68:   PetscCall(MatSetSizes(Pb,PETSC_DECIDE,PETSC_DECIDE,n,n));
 69:   PetscCall(MatSetFromOptions(Pb));

 71:   PetscCall(MatGetOwnershipRange(Pa,&Istart,&Iend));
 72:   for (i=Istart;i<Iend;i++) {
 73:     PetscCall(MatSetValue(Pa,i,i,2.0,INSERT_VALUES));
 74:     if (i>0) PetscCall(MatSetValue(Pb,i,i,(PetscScalar)i,INSERT_VALUES));
 75:     else PetscCall(MatSetValue(Pb,i,i,-1.0,INSERT_VALUES));
 76:   }
 77:   PetscCall(MatAssemblyBegin(Pa,MAT_FINAL_ASSEMBLY));
 78:   PetscCall(MatAssemblyEnd(Pa,MAT_FINAL_ASSEMBLY));
 79:   PetscCall(MatAssemblyBegin(Pb,MAT_FINAL_ASSEMBLY));
 80:   PetscCall(MatAssemblyEnd(Pb,MAT_FINAL_ASSEMBLY));

 82:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 83:                 Create the spectral transformation object
 84:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 86:   PetscCall(STCreate(PETSC_COMM_WORLD,&st));
 87:   mat[0] = A;
 88:   mat[1] = B;
 89:   PetscCall(STSetMatrices(st,2,mat));
 90:   mat[0] = Pa;
 91:   mat[1] = Pb;
 92:   PetscCall(STSetSplitPreconditioner(st,2,mat,SAME_NONZERO_PATTERN));
 93:   PetscCall(STSetTransform(st,PETSC_TRUE));
 94:   PetscCall(STSetFromOptions(st));
 95:   PetscCall(STCayleySetAntishift(st,-0.2));   /* only relevant for cayley */

 97:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 98:                Form the preconditioner matrix and print it
 99:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

101:   PetscCall(STGetKSP(st,&ksp));
102:   PetscCall(KSPGetPC(ksp,&pc));
103:   PetscCall(STGetOperator(st,NULL));
104:   PetscCall(PCGetOperators(pc,NULL,&Pmat));
105:   PetscCall(MatView(Pmat,NULL));

107:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108:                    Apply the operator
109:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

111:   /* sigma=0.0 */
112:   PetscCall(STSetUp(st));
113:   PetscCall(STGetType(st,&type));
114:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
115:   PetscCall(STApply(st,v,w));
116:   PetscCall(VecView(w,NULL));

118:   /* sigma=0.1 */
119:   sigma = 0.1;
120:   PetscCall(STSetShift(st,sigma));
121:   PetscCall(STGetShift(st,&sigma));
122:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
123:   PetscCall(STGetOperator(st,NULL));
124:   PetscCall(PCGetOperators(pc,NULL,&Pmat));
125:   PetscCall(MatView(Pmat,NULL));
126:   PetscCall(STApply(st,v,w));
127:   PetscCall(VecView(w,NULL));

129:   PetscCall(STDestroy(&st));
130:   PetscCall(MatDestroy(&A));
131:   PetscCall(MatDestroy(&B));
132:   PetscCall(MatDestroy(&Pa));
133:   PetscCall(MatDestroy(&Pb));
134:   PetscCall(VecDestroy(&v));
135:   PetscCall(VecDestroy(&w));
136:   PetscCall(SlepcFinalize());
137:   return 0;
138: }

140: /*TEST

142:    test:
143:       suffix: 1
144:       args: -st_type {{cayley shift sinvert}separate output}
145:       requires: !single

147: TEST*/