Actual source code: test1.c

slepc-3.20.2 2024-03-15
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 shell matrices.\n\n";

 13: #include <slepcst.h>

 15: static PetscErrorCode MatGetDiagonal_Shell(Mat S,Vec diag);
 16: static PetscErrorCode MatMultTranspose_Shell(Mat S,Vec x,Vec y);
 17: static PetscErrorCode MatMult_Shell(Mat S,Vec x,Vec y);
 18: static PetscErrorCode MatDuplicate_Shell(Mat S,MatDuplicateOption op,Mat *M);

 20: static PetscErrorCode MyShellMatCreate(Mat *A,Mat *M)
 21: {
 22:   MPI_Comm       comm;
 23:   PetscInt       n;

 25:   PetscFunctionBeginUser;
 26:   PetscCall(MatGetSize(*A,&n,NULL));
 27:   PetscCall(PetscObjectGetComm((PetscObject)*A,&comm));
 28:   PetscCall(MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,n,n,A,M));
 29:   PetscCall(MatShellSetOperation(*M,MATOP_MULT,(void(*)(void))MatMult_Shell));
 30:   PetscCall(MatShellSetOperation(*M,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell));
 31:   PetscCall(MatShellSetOperation(*M,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Shell));
 32:   PetscCall(MatShellSetOperation(*M,MATOP_DUPLICATE,(void(*)(void))MatDuplicate_Shell));
 33:   PetscFunctionReturn(PETSC_SUCCESS);
 34: }

 36: int main(int argc,char **argv)
 37: {
 38:   Mat            A,S,mat[1];
 39:   ST             st;
 40:   Vec            v,w;
 41:   STType         type;
 42:   KSP            ksp;
 43:   PC             pc;
 44:   PetscScalar    sigma;
 45:   PetscInt       n=10,i,Istart,Iend;

 47:   PetscFunctionBeginUser;
 48:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 49:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 50:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian with shell matrices, n=%" PetscInt_FMT "\n\n",n));

 52:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 53:      Compute the operator matrix for the 1-D Laplacian
 54:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 56:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 57:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
 58:   PetscCall(MatSetFromOptions(A));
 59:   PetscCall(MatSetUp(A));

 61:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 62:   for (i=Istart;i<Iend;i++) {
 63:     if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
 64:     if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
 65:     PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
 66:   }
 67:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 68:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));

 70:   /* create the shell version of A */
 71:   PetscCall(MyShellMatCreate(&A,&S));

 73:   /* work vectors */
 74:   PetscCall(MatCreateVecs(A,&v,&w));
 75:   PetscCall(VecSet(v,1.0));

 77:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 78:                 Create the spectral transformation object
 79:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 81:   PetscCall(STCreate(PETSC_COMM_WORLD,&st));
 82:   mat[0] = S;
 83:   PetscCall(STSetMatrices(st,1,mat));
 84:   PetscCall(STSetTransform(st,PETSC_TRUE));
 85:   PetscCall(STSetFromOptions(st));

 87:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 88:               Apply the transformed operator for several ST's
 89:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 91:   /* shift, sigma=0.0 */
 92:   PetscCall(STSetUp(st));
 93:   PetscCall(STGetType(st,&type));
 94:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
 95:   PetscCall(STApply(st,v,w));
 96:   PetscCall(VecView(w,NULL));
 97:   PetscCall(STApplyTranspose(st,v,w));
 98:   PetscCall(VecView(w,NULL));

100:   /* shift, sigma=0.1 */
101:   sigma = 0.1;
102:   PetscCall(STSetShift(st,sigma));
103:   PetscCall(STGetShift(st,&sigma));
104:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
105:   PetscCall(STApply(st,v,w));
106:   PetscCall(VecView(w,NULL));

108:   /* sinvert, sigma=0.1 */
109:   PetscCall(STPostSolve(st));   /* undo changes if inplace */
110:   PetscCall(STSetType(st,STSINVERT));
111:   PetscCall(STGetKSP(st,&ksp));
112:   PetscCall(KSPSetType(ksp,KSPGMRES));
113:   PetscCall(KSPGetPC(ksp,&pc));
114:   PetscCall(PCSetType(pc,PCJACOBI));
115:   PetscCall(STGetType(st,&type));
116:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
117:   PetscCall(STGetShift(st,&sigma));
118:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
119:   PetscCall(STApply(st,v,w));
120:   PetscCall(VecView(w,NULL));

122:   /* sinvert, sigma=-0.5 */
123:   sigma = -0.5;
124:   PetscCall(STSetShift(st,sigma));
125:   PetscCall(STGetShift(st,&sigma));
126:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
127:   PetscCall(STApply(st,v,w));
128:   PetscCall(VecView(w,NULL));

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

139: static PetscErrorCode MatMult_Shell(Mat S,Vec x,Vec y)
140: {
141:   Mat               *A;

143:   PetscFunctionBeginUser;
144:   PetscCall(MatShellGetContext(S,&A));
145:   PetscCall(MatMult(*A,x,y));
146:   PetscFunctionReturn(PETSC_SUCCESS);
147: }

149: static PetscErrorCode MatMultTranspose_Shell(Mat S,Vec x,Vec y)
150: {
151:   Mat               *A;

153:   PetscFunctionBeginUser;
154:   PetscCall(MatShellGetContext(S,&A));
155:   PetscCall(MatMultTranspose(*A,x,y));
156:   PetscFunctionReturn(PETSC_SUCCESS);
157: }

159: static PetscErrorCode MatGetDiagonal_Shell(Mat S,Vec diag)
160: {
161:   Mat               *A;

163:   PetscFunctionBeginUser;
164:   PetscCall(MatShellGetContext(S,&A));
165:   PetscCall(MatGetDiagonal(*A,diag));
166:   PetscFunctionReturn(PETSC_SUCCESS);
167: }

169: static PetscErrorCode MatDuplicate_Shell(Mat S,MatDuplicateOption op,Mat *M)
170: {
171:   Mat            *A;

173:   PetscFunctionBeginUser;
174:   PetscCall(MatShellGetContext(S,&A));
175:   PetscCall(MyShellMatCreate(A,M));
176:   PetscFunctionReturn(PETSC_SUCCESS);
177: }

179: /*TEST

181:    test:
182:       suffix: 1
183:       args: -st_matmode {{inplace shell}}
184:       output_file: output/test1_1.out
185:       requires: !single

187: TEST*/