Actual source code: test5.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 various ST interface functions.\n\n";

 13: #include <slepceps.h>

 15: int main (int argc,char **argv)
 16: {
 17:   ST             st;
 18:   KSP            ksp;
 19:   PC             pc;
 20:   Mat            A,mat[1],Op;
 21:   Vec            v,w;
 22:   PetscInt       N,n=4,i,j,II,Istart,Iend;
 23:   PetscScalar    d;

 25:   PetscFunctionBeginUser;
 26:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 27:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 28:   N = n*n;

 30:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 31:      Create non-symmetric matrix
 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(MatGetOwnershipRange(A,&Istart,&Iend));
 40:   for (II=Istart;II<Iend;II++) {
 41:     i = II/n; j = II-i*n;
 42:     d = 0.0;
 43:     if (i>0) { PetscCall(MatSetValue(A,II,II-n,1.0,INSERT_VALUES)); d=d+1.0; }
 44:     if (i<n-1) { PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES)); d=d+1.0; }
 45:     if (j>0) { PetscCall(MatSetValue(A,II,II-1,1.0,INSERT_VALUES)); d=d+1.0; }
 46:     if (j<n-1) { PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES)); d=d+1.0; }
 47:     PetscCall(MatSetValue(A,II,II,d,INSERT_VALUES));
 48:   }

 50:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 51:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
 52:   PetscCall(MatCreateVecs(A,&v,&w));
 53:   PetscCall(VecSetValue(v,0,-.5,INSERT_VALUES));
 54:   PetscCall(VecSetValue(v,1,1.5,INSERT_VALUES));
 55:   PetscCall(VecSetValue(v,2,2,INSERT_VALUES));
 56:   PetscCall(VecAssemblyBegin(v));
 57:   PetscCall(VecAssemblyEnd(v));
 58:   PetscCall(VecView(v,NULL));

 60:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 61:      Create the spectral transformation object
 62:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 63:   PetscCall(STCreate(PETSC_COMM_WORLD,&st));
 64:   mat[0] = A;
 65:   PetscCall(STSetMatrices(st,1,mat));
 66:   PetscCall(STSetType(st,STCAYLEY));
 67:   PetscCall(STSetShift(st,2.0));
 68:   PetscCall(STCayleySetAntishift(st,1.0));
 69:   PetscCall(STSetTransform(st,PETSC_TRUE));

 71:   PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp));
 72:   PetscCall(KSPSetType(ksp,KSPPREONLY));
 73:   PetscCall(KSPGetPC(ksp,&pc));
 74:   PetscCall(PCSetType(pc,PCLU));
 75:   PetscCall(KSPSetTolerances(ksp,100*PETSC_MACHINE_EPSILON,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
 76:   PetscCall(KSPSetFromOptions(ksp));
 77:   PetscCall(STSetKSP(st,ksp));
 78:   PetscCall(KSPDestroy(&ksp));

 80:   PetscCall(STSetFromOptions(st));

 82:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 83:      Apply the operator
 84:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 85:   PetscCall(STApply(st,v,w));
 86:   PetscCall(VecView(w,NULL));

 88:   PetscCall(STApplyTranspose(st,v,w));
 89:   PetscCall(VecView(w,NULL));

 91:   PetscCall(STMatMult(st,1,v,w));
 92:   PetscCall(VecView(w,NULL));

 94:   PetscCall(STMatMultTranspose(st,1,v,w));
 95:   PetscCall(VecView(w,NULL));

 97:   PetscCall(STMatSolve(st,v,w));
 98:   PetscCall(VecView(w,NULL));

100:   PetscCall(STMatSolveTranspose(st,v,w));
101:   PetscCall(VecView(w,NULL));

103:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104:      Get the operator matrix
105:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106:   PetscCall(STGetOperator(st,&Op));
107:   PetscCall(MatMult(Op,v,w));
108:   PetscCall(VecView(w,NULL));
109:   PetscCall(MatMultTranspose(Op,v,w));
110:   PetscCall(VecView(w,NULL));
111:   PetscCall(STRestoreOperator(st,&Op));

113:   PetscCall(STDestroy(&st));
114:   PetscCall(MatDestroy(&A));
115:   PetscCall(VecDestroy(&v));
116:   PetscCall(VecDestroy(&w));
117:   PetscCall(SlepcFinalize());
118:   return 0;
119: }

121: /*TEST

123:    testset:
124:       output_file: output/test5_1.out
125:       requires: !single
126:       test:
127:          args: -st_matmode {{copy inplace}}
128:       test:
129:          args: -st_matmode shell -ksp_type bcgs -pc_type jacobi

131: TEST*/