Actual source code: test3.c
slepc-3.20.2 2024-03-15
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 two matrices.\n\n";
13: #include <slepcst.h>
15: int main(int argc,char **argv)
16: {
17: Mat A,B,M,mat[2];
18: ST st;
19: Vec v,w;
20: STType type;
21: PetscScalar sigma,tau;
22: PetscInt n=10,i,Istart,Iend;
24: PetscFunctionBeginUser;
25: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
26: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
27: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian plus diagonal, n=%" PetscInt_FMT "\n\n",n));
29: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30: Compute the operator matrix for the 1-D Laplacian
31: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
33: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
34: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
35: PetscCall(MatSetFromOptions(A));
36: PetscCall(MatSetUp(A));
38: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
39: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
40: PetscCall(MatSetFromOptions(B));
41: PetscCall(MatSetUp(B));
43: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
44: for (i=Istart;i<Iend;i++) {
45: PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
46: if (i>0) {
47: PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
48: PetscCall(MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES));
49: } else PetscCall(MatSetValue(B,i,i,-1.0,INSERT_VALUES));
50: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
51: }
52: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
53: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
54: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
55: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
56: PetscCall(MatCreateVecs(A,&v,&w));
57: PetscCall(VecSet(v,1.0));
59: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60: Create the spectral transformation object
61: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: PetscCall(STCreate(PETSC_COMM_WORLD,&st));
64: mat[0] = A;
65: mat[1] = B;
66: PetscCall(STSetMatrices(st,2,mat));
67: PetscCall(STSetTransform(st,PETSC_TRUE));
68: PetscCall(STSetFromOptions(st));
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Apply the transformed operator for several ST's
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74: /* shift, sigma=0.0 */
75: PetscCall(STSetUp(st));
76: PetscCall(STGetType(st,&type));
77: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
78: PetscCall(STApply(st,v,w));
79: PetscCall(VecView(w,NULL));
81: /* shift, sigma=0.1 */
82: sigma = 0.1;
83: PetscCall(STSetShift(st,sigma));
84: PetscCall(STGetShift(st,&sigma));
85: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
86: PetscCall(STApply(st,v,w));
87: PetscCall(VecView(w,NULL));
89: /* sinvert, sigma=0.1 */
90: PetscCall(STPostSolve(st)); /* undo changes if inplace */
91: PetscCall(STSetType(st,STSINVERT));
92: PetscCall(STGetType(st,&type));
93: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
94: PetscCall(STGetShift(st,&sigma));
95: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
96: PetscCall(STApply(st,v,w));
97: PetscCall(VecView(w,NULL));
99: /* sinvert, sigma=-0.5 */
100: sigma = -0.5;
101: PetscCall(STSetShift(st,sigma));
102: PetscCall(STGetShift(st,&sigma));
103: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
104: PetscCall(STApply(st,v,w));
105: PetscCall(VecView(w,NULL));
107: /* cayley, sigma=-0.5, tau=-0.5 (equal to sigma by default) */
108: PetscCall(STPostSolve(st)); /* undo changes if inplace */
109: PetscCall(STSetType(st,STCAYLEY));
110: PetscCall(STSetUp(st));
111: PetscCall(STGetType(st,&type));
112: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
113: PetscCall(STGetShift(st,&sigma));
114: PetscCall(STCayleyGetAntishift(st,&tau));
115: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau)));
116: PetscCall(STApply(st,v,w));
117: PetscCall(VecView(w,NULL));
119: /* cayley, sigma=1.1, tau=1.1 (still equal to sigma) */
120: sigma = 1.1;
121: PetscCall(STSetShift(st,sigma));
122: PetscCall(STGetShift(st,&sigma));
123: PetscCall(STCayleyGetAntishift(st,&tau));
124: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau)));
125: PetscCall(STApply(st,v,w));
126: PetscCall(VecView(w,NULL));
128: /* cayley, sigma=1.1, tau=-1.0 */
129: tau = -1.0;
130: PetscCall(STCayleySetAntishift(st,tau));
131: PetscCall(STGetShift(st,&sigma));
132: PetscCall(STCayleyGetAntishift(st,&tau));
133: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau)));
134: PetscCall(STApply(st,v,w));
135: PetscCall(VecView(w,NULL));
137: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138: Check inner product matrix in Cayley
139: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140: PetscCall(STGetBilinearForm(st,&M));
141: PetscCall(MatMult(M,v,w));
142: PetscCall(VecView(w,NULL));
144: PetscCall(STDestroy(&st));
145: PetscCall(MatDestroy(&A));
146: PetscCall(MatDestroy(&B));
147: PetscCall(MatDestroy(&M));
148: PetscCall(VecDestroy(&v));
149: PetscCall(VecDestroy(&w));
150: PetscCall(SlepcFinalize());
151: return 0;
152: }
154: /*TEST
156: test:
157: suffix: 1
158: args: -st_matmode {{copy inplace shell}}
159: output_file: output/test3_1.out
160: requires: !single
162: TEST*/