Actual source code: test5.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 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: SlepcInitialize(&argc,&argv,(char*)0,help);
26: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
27: N = n*n;
29: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30: Create non-symmetric matrix
31: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
33: MatCreate(PETSC_COMM_WORLD,&A);
34: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
35: MatSetFromOptions(A);
36: MatSetUp(A);
38: MatGetOwnershipRange(A,&Istart,&Iend);
39: for (II=Istart;II<Iend;II++) {
40: i = II/n; j = II-i*n;
41: d = 0.0;
42: if (i>0) { MatSetValue(A,II,II-n,1.0,INSERT_VALUES); d=d+1.0; }
43: if (i<n-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); d=d+1.0; }
44: if (j>0) { MatSetValue(A,II,II-1,1.0,INSERT_VALUES); d=d+1.0; }
45: if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); d=d+1.0; }
46: MatSetValue(A,II,II,d,INSERT_VALUES);
47: }
49: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
50: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
51: MatCreateVecs(A,&v,&w);
52: VecSetValue(v,0,-.5,INSERT_VALUES);
53: VecSetValue(v,1,1.5,INSERT_VALUES);
54: VecSetValue(v,2,2,INSERT_VALUES);
55: VecAssemblyBegin(v);
56: VecAssemblyEnd(v);
57: VecView(v,NULL);
59: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60: Create the spectral transformation object
61: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
62: STCreate(PETSC_COMM_WORLD,&st);
63: mat[0] = A;
64: STSetMatrices(st,1,mat);
65: STSetType(st,STCAYLEY);
66: STSetShift(st,2.0);
67: STCayleySetAntishift(st,1.0);
68: STSetTransform(st,PETSC_TRUE);
70: KSPCreate(PETSC_COMM_WORLD,&ksp);
71: KSPSetType(ksp,KSPPREONLY);
72: KSPGetPC(ksp,&pc);
73: PCSetType(pc,PCLU);
74: KSPSetTolerances(ksp,100*PETSC_MACHINE_EPSILON,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
75: KSPSetFromOptions(ksp);
76: STSetKSP(st,ksp);
77: KSPDestroy(&ksp);
79: STSetFromOptions(st);
81: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82: Apply the operator
83: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84: STApply(st,v,w);
85: VecView(w,NULL);
87: STApplyTranspose(st,v,w);
88: VecView(w,NULL);
90: STMatMult(st,1,v,w);
91: VecView(w,NULL);
93: STMatMultTranspose(st,1,v,w);
94: VecView(w,NULL);
96: STMatSolve(st,v,w);
97: VecView(w,NULL);
99: STMatSolveTranspose(st,v,w);
100: VecView(w,NULL);
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Get the operator matrix
104: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105: STGetOperator(st,&Op);
106: MatMult(Op,v,w);
107: VecView(w,NULL);
108: MatMultTranspose(Op,v,w);
109: VecView(w,NULL);
110: STRestoreOperator(st,&Op);
112: STDestroy(&st);
113: MatDestroy(&A);
114: VecDestroy(&v);
115: VecDestroy(&w);
116: SlepcFinalize();
117: return 0;
118: }
120: /*TEST
122: testset:
123: output_file: output/test5_1.out
124: requires: !single
125: test:
126: args: -st_matmode {{copy inplace}}
127: test:
128: args: -st_matmode shell -ksp_type bcgs -pc_type jacobi
130: TEST*/