Actual source code: test1.c
slepc-3.5.2 2014-10-10
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Test ST with shell matrices.\n\n";
24: #include <slepcst.h>
26: static PetscErrorCode MatGetDiagonal_Shell(Mat S,Vec diag);
27: static PetscErrorCode MatMultTranspose_Shell(Mat S,Vec x,Vec y);
28: static PetscErrorCode MatMult_Shell(Mat S,Vec x,Vec y);
29: static PetscErrorCode MatDuplicate_Shell(Mat S,MatDuplicateOption op,Mat *M);
33: static PetscErrorCode MyShellMatCreate(Mat *A,Mat *M)
34: {
36: MPI_Comm comm;
37: PetscInt n;
40: MatGetSize(*A,&n,NULL);
41: PetscObjectGetComm((PetscObject)*A,&comm);
42: MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,n,n,A,M);
43: MatSetFromOptions(*M);
44: MatShellSetOperation(*M,MATOP_MULT,(void(*)())MatMult_Shell);
45: MatShellSetOperation(*M,MATOP_MULT_TRANSPOSE,(void(*)())MatMultTranspose_Shell);
46: MatShellSetOperation(*M,MATOP_GET_DIAGONAL,(void(*)())MatGetDiagonal_Shell);
47: MatShellSetOperation(*M,MATOP_DUPLICATE,(void(*)())MatDuplicate_Shell);
48: return(0);
49: }
53: int main(int argc,char **argv)
54: {
55: Mat A,S,mat[1];
56: ST st;
57: Vec v,w;
58: STType type;
59: KSP ksp;
60: PC pc;
61: PetscScalar value[3],sigma;
62: PetscInt n=10,i,Istart,Iend,col[3];
63: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
66: SlepcInitialize(&argc,&argv,(char*)0,help);
67: PetscOptionsGetInt(NULL,"-n",&n,NULL);
68: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian with shell matrices, n=%D\n\n",n);
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Compute the operator matrix for the 1-D Laplacian
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74: MatCreate(PETSC_COMM_WORLD,&A);
75: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
76: MatSetFromOptions(A);
77: MatSetUp(A);
79: MatGetOwnershipRange(A,&Istart,&Iend);
80: if (Istart==0) FirstBlock=PETSC_TRUE;
81: if (Iend==n) LastBlock=PETSC_TRUE;
82: value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
83: for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
84: col[0]=i-1; col[1]=i; col[2]=i+1;
85: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
86: }
87: if (LastBlock) {
88: i=n-1; col[0]=n-2; col[1]=n-1;
89: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
90: }
91: if (FirstBlock) {
92: i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
93: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
94: }
96: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
97: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
99: /* create the shell version of A */
100: MyShellMatCreate(&A,&S);
102: /* work vectors */
103: MatGetVecs(A,&v,&w);
104: VecSet(v,1.0);
106: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107: Create the spectral transformation object
108: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110: STCreate(PETSC_COMM_WORLD,&st);
111: mat[0] = S;
112: STSetOperators(st,1,mat);
113: STSetFromOptions(st);
116: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117: Apply the transformed operator for several ST's
118: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120: /* shift, sigma=0.0 */
121: STSetUp(st);
122: STGetType(st,&type);
123: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
124: STApply(st,v,w);
125: VecView(w,NULL);
127: /* shift, sigma=0.1 */
128: sigma = 0.1;
129: STSetShift(st,sigma);
130: STGetShift(st,&sigma);
131: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
132: STApply(st,v,w);
133: VecView(w,NULL);
135: /* sinvert, sigma=0.1 */
136: STPostSolve(st); /* undo changes if inplace */
137: STSetType(st,STSINVERT);
138: STGetKSP(st,&ksp);
139: KSPSetType(ksp,KSPGMRES);
140: KSPGetPC(ksp,&pc);
141: PCSetType(pc,PCJACOBI);
142: STGetType(st,&type);
143: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
144: STGetShift(st,&sigma);
145: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
146: STApply(st,v,w);
147: VecView(w,NULL);
149: /* sinvert, sigma=-0.5 */
150: sigma = -0.5;
151: STSetShift(st,sigma);
152: STGetShift(st,&sigma);
153: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
154: STApply(st,v,w);
155: VecView(w,NULL);
157: STDestroy(&st);
158: MatDestroy(&A);
159: MatDestroy(&S);
160: VecDestroy(&v);
161: VecDestroy(&w);
162: SlepcFinalize();
163: return 0;
164: }
168: static PetscErrorCode MatMult_Shell(Mat S,Vec x,Vec y)
169: {
170: PetscErrorCode ierr;
171: Mat *A;
174: MatShellGetContext(S,&A);
175: MatMult(*A,x,y);
176: return(0);
177: }
181: static PetscErrorCode MatMultTranspose_Shell(Mat S,Vec x,Vec y)
182: {
183: PetscErrorCode ierr;
184: Mat *A;
187: MatShellGetContext(S,&A);
188: MatMultTranspose(*A,x,y);
189: return(0);
190: }
194: static PetscErrorCode MatGetDiagonal_Shell(Mat S,Vec diag)
195: {
196: PetscErrorCode ierr;
197: Mat *A;
200: MatShellGetContext(S,&A);
201: MatGetDiagonal(*A,diag);
202: return(0);
203: }
207: static PetscErrorCode MatDuplicate_Shell(Mat S,MatDuplicateOption op,Mat *M)
208: {
210: Mat *A;
213: MatShellGetContext(S,&A);
214: MyShellMatCreate(A,M);
215: return(0);
216: }