Actual source code: test9.c

slepc-3.19.0 2023-03-31
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 four matrices and split preconditioner.\n\n";

 13: #include <slepcst.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat            A,B,C,D,Pa,Pb,Pc,Pd,Pmat,mat[4];
 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,"\nTest ST with four matrices, n=%" PetscInt_FMT "\n\n",n));
 30:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 31:      Compute the operator matrices
 32:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

 44:   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
 45:   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n));
 46:   PetscCall(MatSetFromOptions(C));
 47:   PetscCall(MatSetUp(C));

 49:   PetscCall(MatCreate(PETSC_COMM_WORLD,&D));
 50:   PetscCall(MatSetSizes(D,PETSC_DECIDE,PETSC_DECIDE,n,n));
 51:   PetscCall(MatSetFromOptions(D));
 52:   PetscCall(MatSetUp(D));

 54:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 55:   for (i=Istart;i<Iend;i++) {
 56:     PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
 57:     if (i>0) {
 58:       PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
 59:       PetscCall(MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES));
 60:     } else PetscCall(MatSetValue(B,i,i,-1.0,INSERT_VALUES));
 61:     if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
 62:     PetscCall(MatSetValue(C,i,n-i-1,1.0,INSERT_VALUES));
 63:     PetscCall(MatSetValue(D,i,i,i*.1,INSERT_VALUES));
 64:     if (i==0) PetscCall(MatSetValue(D,0,n-1,1.0,INSERT_VALUES));
 65:     if (i==n-1) PetscCall(MatSetValue(D,n-1,0,1.0,INSERT_VALUES));
 66:   }

 68:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 69:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
 70:   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
 71:   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
 72:   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
 73:   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
 74:   PetscCall(MatAssemblyBegin(D,MAT_FINAL_ASSEMBLY));
 75:   PetscCall(MatAssemblyEnd(D,MAT_FINAL_ASSEMBLY));
 76:   PetscCall(MatCreateVecs(A,&v,&w));
 77:   PetscCall(VecSet(v,1.0));

 79:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 80:      Compute the split preconditioner matrices (four diagonals)
 81:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 83:   PetscCall(MatCreate(PETSC_COMM_WORLD,&Pa));
 84:   PetscCall(MatSetSizes(Pa,PETSC_DECIDE,PETSC_DECIDE,n,n));
 85:   PetscCall(MatSetFromOptions(Pa));
 86:   PetscCall(MatSetUp(Pa));

 88:   PetscCall(MatCreate(PETSC_COMM_WORLD,&Pb));
 89:   PetscCall(MatSetSizes(Pb,PETSC_DECIDE,PETSC_DECIDE,n,n));
 90:   PetscCall(MatSetFromOptions(Pb));
 91:   PetscCall(MatSetUp(Pb));

 93:   PetscCall(MatCreate(PETSC_COMM_WORLD,&Pc));
 94:   PetscCall(MatSetSizes(Pc,PETSC_DECIDE,PETSC_DECIDE,n,n));
 95:   PetscCall(MatSetFromOptions(Pc));
 96:   PetscCall(MatSetUp(Pc));

 98:   PetscCall(MatCreate(PETSC_COMM_WORLD,&Pd));
 99:   PetscCall(MatSetSizes(Pd,PETSC_DECIDE,PETSC_DECIDE,n,n));
100:   PetscCall(MatSetFromOptions(Pd));
101:   PetscCall(MatSetUp(Pd));

103:   PetscCall(MatGetOwnershipRange(Pa,&Istart,&Iend));
104:   for (i=Istart;i<Iend;i++) {
105:     PetscCall(MatSetValue(Pa,i,i,2.0,INSERT_VALUES));
106:     if (i>0) PetscCall(MatSetValue(Pb,i,i,(PetscScalar)i,INSERT_VALUES));
107:     else PetscCall(MatSetValue(Pb,i,i,-1.0,INSERT_VALUES));
108:     PetscCall(MatSetValue(Pd,i,i,i*.1,INSERT_VALUES));
109:   }

111:   PetscCall(MatAssemblyBegin(Pa,MAT_FINAL_ASSEMBLY));
112:   PetscCall(MatAssemblyEnd(Pa,MAT_FINAL_ASSEMBLY));
113:   PetscCall(MatAssemblyBegin(Pb,MAT_FINAL_ASSEMBLY));
114:   PetscCall(MatAssemblyEnd(Pb,MAT_FINAL_ASSEMBLY));
115:   PetscCall(MatAssemblyBegin(Pc,MAT_FINAL_ASSEMBLY));
116:   PetscCall(MatAssemblyEnd(Pc,MAT_FINAL_ASSEMBLY));
117:   PetscCall(MatAssemblyBegin(Pd,MAT_FINAL_ASSEMBLY));
118:   PetscCall(MatAssemblyEnd(Pd,MAT_FINAL_ASSEMBLY));

120:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121:                 Create the spectral transformation object
122:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

124:   PetscCall(STCreate(PETSC_COMM_WORLD,&st));
125:   mat[0] = A;
126:   mat[1] = B;
127:   mat[2] = C;
128:   mat[3] = D;
129:   PetscCall(STSetMatrices(st,4,mat));
130:   mat[0] = Pa;
131:   mat[1] = Pb;
132:   mat[2] = Pc;
133:   mat[3] = Pd;
134:   PetscCall(STSetSplitPreconditioner(st,4,mat,SUBSET_NONZERO_PATTERN));
135:   PetscCall(STGetKSP(st,&ksp));
136:   PetscCall(KSPSetTolerances(ksp,100*PETSC_MACHINE_EPSILON,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
137:   PetscCall(STSetTransform(st,PETSC_TRUE));
138:   PetscCall(STSetFromOptions(st));
139:   PetscCall(STGetKSP(st,&ksp));
140:   PetscCall(KSPGetPC(ksp,&pc));

142:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143:                    Apply the operator
144:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

146:   /* sigma=0.0 */
147:   PetscCall(STSetUp(st));
148:   PetscCall(STGetType(st,&type));
149:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
150:   PetscCall(PCGetOperators(pc,NULL,&Pmat));
151:   PetscCall(MatView(Pmat,NULL));
152:   PetscCall(STMatSolve(st,v,w));
153:   PetscCall(VecView(w,NULL));

155:   /* sigma=0.1 */
156:   sigma = 0.1;
157:   PetscCall(STSetShift(st,sigma));
158:   PetscCall(STGetShift(st,&sigma));
159:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
160:   PetscCall(PCGetOperators(pc,NULL,&Pmat));
161:   PetscCall(MatView(Pmat,NULL));
162:   PetscCall(STMatSolve(st,v,w));
163:   PetscCall(VecView(w,NULL));

165:   PetscCall(STDestroy(&st));
166:   PetscCall(MatDestroy(&A));
167:   PetscCall(MatDestroy(&B));
168:   PetscCall(MatDestroy(&C));
169:   PetscCall(MatDestroy(&D));
170:   PetscCall(MatDestroy(&Pa));
171:   PetscCall(MatDestroy(&Pb));
172:   PetscCall(MatDestroy(&Pc));
173:   PetscCall(MatDestroy(&Pd));
174:   PetscCall(VecDestroy(&v));
175:   PetscCall(VecDestroy(&w));
176:   PetscCall(SlepcFinalize());
177:   return 0;
178: }

180: /*TEST

182:    test:
183:       suffix: 1
184:       args: -st_type {{shift sinvert}separate output} -st_pc_type jacobi
185:       requires: !single

187: TEST*/