Actual source code: test8.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 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;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);
 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   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:   MatCreate(PETSC_COMM_WORLD,&A);
 36:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 37:   MatSetFromOptions(A);
 38:   MatSetUp(A);

 40:   MatCreate(PETSC_COMM_WORLD,&B);
 41:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 42:   MatSetFromOptions(B);
 43:   MatSetUp(B);

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

 61:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 62:      Compute the split preconditioner matrices (two diagonals)
 63:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 65:   MatCreate(PETSC_COMM_WORLD,&Pa);
 66:   MatSetSizes(Pa,PETSC_DECIDE,PETSC_DECIDE,n,n);
 67:   MatSetFromOptions(Pa);
 68:   MatSetUp(Pa);

 70:   MatCreate(PETSC_COMM_WORLD,&Pb);
 71:   MatSetSizes(Pb,PETSC_DECIDE,PETSC_DECIDE,n,n);
 72:   MatSetFromOptions(Pb);
 73:   MatSetUp(Pb);

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

 86:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 87:                 Create the spectral transformation object
 88:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

101:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102:                Form the preconditioner matrix and print it
103:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

105:   STGetKSP(st,&ksp);
106:   KSPGetPC(ksp,&pc);
107:   STGetOperator(st,NULL);
108:   PCGetOperators(pc,NULL,&Pmat);
109:   MatView(Pmat,NULL);

111:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112:                    Apply the operator
113:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

115:   /* sigma=0.0 */
116:   STSetUp(st);
117:   STGetType(st,&type);
118:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
119:   STApply(st,v,w);
120:   VecView(w,NULL);

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

133:   STDestroy(&st);
134:   MatDestroy(&A);
135:   MatDestroy(&B);
136:   MatDestroy(&Pa);
137:   MatDestroy(&Pb);
138:   VecDestroy(&v);
139:   VecDestroy(&w);
140:   SlepcFinalize();
141:   return 0;
142: }

144: /*TEST

146:    test:
147:       suffix: 1
148:       args: -st_type {{cayley shift sinvert}separate output}
149:       requires: !single

151: TEST*/