Actual source code: test5.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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;
26: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
27: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
28: N = n*n;
30: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31: Create non-symmetric matrix
32: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
34: MatCreate(PETSC_COMM_WORLD,&A);
35: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
36: MatSetFromOptions(A);
37: MatSetUp(A);
39: 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) { MatSetValue(A,II,II-n,1.0,INSERT_VALUES); d=d+1.0; }
44: if (i<n-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); d=d+1.0; }
45: if (j>0) { MatSetValue(A,II,II-1,1.0,INSERT_VALUES); d=d+1.0; }
46: if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); d=d+1.0; }
47: MatSetValue(A,II,II,d,INSERT_VALUES);
48: }
50: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
51: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
52: MatCreateVecs(A,&v,&w);
53: VecSetValue(v,0,-.5,INSERT_VALUES);
54: VecSetValue(v,1,1.5,INSERT_VALUES);
55: VecSetValue(v,2,2,INSERT_VALUES);
56: VecAssemblyBegin(v);
57: VecAssemblyEnd(v);
58: VecView(v,NULL);
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Create the spectral transformation object
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: STCreate(PETSC_COMM_WORLD,&st);
64: mat[0] = A;
65: STSetMatrices(st,1,mat);
66: STSetType(st,STCAYLEY);
67: STSetShift(st,2.0);
68: STCayleySetAntishift(st,1.0);
69: STSetTransform(st,PETSC_TRUE);
71: KSPCreate(PETSC_COMM_WORLD,&ksp);
72: KSPSetType(ksp,KSPPREONLY);
73: KSPGetPC(ksp,&pc);
74: PCSetType(pc,PCLU);
75: KSPSetTolerances(ksp,100*PETSC_MACHINE_EPSILON,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
76: KSPSetFromOptions(ksp);
77: STSetKSP(st,ksp);
78: KSPDestroy(&ksp);
80: STSetFromOptions(st);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Apply the operator
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85: STApply(st,v,w);
86: VecView(w,NULL);
88: STApplyTranspose(st,v,w);
89: VecView(w,NULL);
91: STMatMult(st,1,v,w);
92: VecView(w,NULL);
94: STMatMultTranspose(st,1,v,w);
95: VecView(w,NULL);
97: STMatSolve(st,v,w);
98: VecView(w,NULL);
100: STMatSolveTranspose(st,v,w);
101: VecView(w,NULL);
103: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: Get the operator matrix
105: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106: STGetOperator(st,&Op);
107: MatMult(Op,v,w);
108: VecView(w,NULL);
109: MatMultTranspose(Op,v,w);
110: VecView(w,NULL);
111: MatDestroy(&Op);
113: STDestroy(&st);
114: MatDestroy(&A);
115: VecDestroy(&v);
116: VecDestroy(&w);
117: SlepcFinalize();
118: return ierr;
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*/