Actual source code: test1.c
slepc-3.17.2 2022-08-09
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;
26: MatGetSize(*A,&n,NULL);
27: PetscObjectGetComm((PetscObject)*A,&comm);
28: MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,n,n,A,M);
29: MatShellSetOperation(*M,MATOP_MULT,(void(*)(void))MatMult_Shell);
30: MatShellSetOperation(*M,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell);
31: MatShellSetOperation(*M,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Shell);
32: MatShellSetOperation(*M,MATOP_DUPLICATE,(void(*)(void))MatDuplicate_Shell);
33: PetscFunctionReturn(0);
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: SlepcInitialize(&argc,&argv,(char*)0,help);
48: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
49: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian with shell matrices, n=%" PetscInt_FMT "\n\n",n);
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Compute the operator matrix for the 1-D Laplacian
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: MatCreate(PETSC_COMM_WORLD,&A);
56: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
57: MatSetFromOptions(A);
58: MatSetUp(A);
60: MatGetOwnershipRange(A,&Istart,&Iend);
61: for (i=Istart;i<Iend;i++) {
62: if (i>0) MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
63: if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
64: MatSetValue(A,i,i,2.0,INSERT_VALUES);
65: }
66: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
67: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
69: /* create the shell version of A */
70: MyShellMatCreate(&A,&S);
72: /* work vectors */
73: MatCreateVecs(A,&v,&w);
74: VecSet(v,1.0);
76: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77: Create the spectral transformation object
78: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
80: STCreate(PETSC_COMM_WORLD,&st);
81: mat[0] = S;
82: STSetMatrices(st,1,mat);
83: STSetTransform(st,PETSC_TRUE);
84: STSetFromOptions(st);
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Apply the transformed operator for several ST's
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: /* shift, sigma=0.0 */
91: STSetUp(st);
92: STGetType(st,&type);
93: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
94: STApply(st,v,w);
95: VecView(w,NULL);
96: STApplyTranspose(st,v,w);
97: VecView(w,NULL);
99: /* shift, sigma=0.1 */
100: sigma = 0.1;
101: STSetShift(st,sigma);
102: STGetShift(st,&sigma);
103: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
104: STApply(st,v,w);
105: VecView(w,NULL);
107: /* sinvert, sigma=0.1 */
108: STPostSolve(st); /* undo changes if inplace */
109: STSetType(st,STSINVERT);
110: STGetKSP(st,&ksp);
111: KSPSetType(ksp,KSPGMRES);
112: KSPGetPC(ksp,&pc);
113: PCSetType(pc,PCJACOBI);
114: STGetType(st,&type);
115: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
116: STGetShift(st,&sigma);
117: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
118: STApply(st,v,w);
119: VecView(w,NULL);
121: /* sinvert, sigma=-0.5 */
122: sigma = -0.5;
123: STSetShift(st,sigma);
124: STGetShift(st,&sigma);
125: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
126: STApply(st,v,w);
127: VecView(w,NULL);
129: STDestroy(&st);
130: MatDestroy(&A);
131: MatDestroy(&S);
132: VecDestroy(&v);
133: VecDestroy(&w);
134: SlepcFinalize();
135: return 0;
136: }
138: static PetscErrorCode MatMult_Shell(Mat S,Vec x,Vec y)
139: {
140: Mat *A;
143: MatShellGetContext(S,&A);
144: MatMult(*A,x,y);
145: PetscFunctionReturn(0);
146: }
148: static PetscErrorCode MatMultTranspose_Shell(Mat S,Vec x,Vec y)
149: {
150: Mat *A;
153: MatShellGetContext(S,&A);
154: MatMultTranspose(*A,x,y);
155: PetscFunctionReturn(0);
156: }
158: static PetscErrorCode MatGetDiagonal_Shell(Mat S,Vec diag)
159: {
160: Mat *A;
163: MatShellGetContext(S,&A);
164: MatGetDiagonal(*A,diag);
165: PetscFunctionReturn(0);
166: }
168: static PetscErrorCode MatDuplicate_Shell(Mat S,MatDuplicateOption op,Mat *M)
169: {
170: Mat *A;
173: MatShellGetContext(S,&A);
174: MyShellMatCreate(A,M);
175: PetscFunctionReturn(0);
176: }
178: /*TEST
180: test:
181: suffix: 1
182: args: -st_matmode {{inplace shell}}
183: output_file: output/test1_1.out
184: requires: !single
186: TEST*/