Actual source code: test6.c
slepc-3.17.1 2022-04-11
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 STPRECOND operations.\n\n";
13: #include <slepcst.h>
15: int main(int argc,char **argv)
16: {
17: Mat A,P,mat[1];
18: ST st;
19: KSP ksp;
20: Vec v,w;
21: PetscScalar sigma;
22: PetscInt n=10,i,Istart,Iend;
23: STMatMode matmode;
25: SlepcInitialize(&argc,&argv,(char*)0,help);
26: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
27: PetscPrintf(PETSC_COMM_WORLD,"\nPreconditioner for 1-D Laplacian, n=%" PetscInt_FMT "\n\n",n);
29: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30: Compute the operator matrix for the 1-D Laplacian
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 (i=Istart;i<Iend;i++) {
40: if (i>0) MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
41: if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
42: MatSetValue(A,i,i,2.0,INSERT_VALUES);
43: }
44: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
45: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
46: MatCreateVecs(A,&v,&w);
47: VecSet(v,1.0);
49: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50: Create the spectral transformation object
51: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53: STCreate(PETSC_COMM_WORLD,&st);
54: mat[0] = A;
55: STSetMatrices(st,1,mat);
56: STSetType(st,STPRECOND);
57: STGetKSP(st,&ksp);
58: KSPSetType(ksp,KSPPREONLY);
59: STSetFromOptions(st);
61: /* set up */
62: /* - the transform flag is necessary so that A-sigma*I is built explicitly */
63: STSetTransform(st,PETSC_TRUE);
64: /* - the ksphasmat flag is necessary when using STApply(), otherwise can only use PCApply() */
65: STPrecondSetKSPHasMat(st,PETSC_TRUE);
66: /* no need to call STSetUp() explicitly */
67: STSetUp(st);
69: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
70: Apply the preconditioner
71: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73: /* default shift */
74: PetscPrintf(PETSC_COMM_WORLD,"With default shift\n");
75: STApply(st,v,w);
76: VecView(w,NULL);
78: /* change shift */
79: sigma = 0.1;
80: STSetShift(st,sigma);
81: STGetShift(st,&sigma);
82: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
83: STApply(st,v,w);
84: VecView(w,NULL);
85: STPostSolve(st); /* undo changes if inplace */
87: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: Test a user-provided preconditioner matrix
89: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91: MatCreate(PETSC_COMM_WORLD,&P);
92: MatSetSizes(P,PETSC_DECIDE,PETSC_DECIDE,n,n);
93: MatSetFromOptions(P);
94: MatSetUp(P);
96: MatGetOwnershipRange(P,&Istart,&Iend);
97: for (i=Istart;i<Iend;i++) MatSetValue(P,i,i,2.0,INSERT_VALUES);
98: if (Istart==0) {
99: MatSetValue(P,1,0,-1.0,INSERT_VALUES);
100: MatSetValue(P,0,1,-1.0,INSERT_VALUES);
101: }
102: MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);
103: MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY);
105: /* apply new preconditioner */
106: STSetPreconditionerMat(st,P);
107: PetscPrintf(PETSC_COMM_WORLD,"With user-provided matrix\n");
108: STApply(st,v,w);
109: VecView(w,NULL);
111: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: Test a user-provided preconditioner in split form
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115: STGetMatMode(st,&matmode);
116: if (matmode==ST_MATMODE_COPY) {
117: STSetPreconditionerMat(st,NULL);
118: mat[0] = P;
119: STSetSplitPreconditioner(st,1,mat,SAME_NONZERO_PATTERN);
121: /* apply new preconditioner */
122: PetscPrintf(PETSC_COMM_WORLD,"With split preconditioner\n");
123: STApply(st,v,w);
124: VecView(w,NULL);
125: }
127: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: Clean up
129: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: STDestroy(&st);
131: MatDestroy(&A);
132: MatDestroy(&P);
133: VecDestroy(&v);
134: VecDestroy(&w);
135: SlepcFinalize();
136: return 0;
137: }
139: /*TEST
141: test:
142: suffix: 1
143: args: -st_matmode {{copy inplace shell}separate output}
144: requires: !single
146: TEST*/