Actual source code: iccbs.c
1: /*
2: Defines a Cholesky factorization preconditioner with BlockSolve95 interface.
4: Note that BlockSolve95 works with a scaled and permuted preconditioning matrix.
5: If the linear system matrix and preconditioning matrix are the same, we then
6: work directly with the permuted and scaled linear system:
7: - original system: Ax = b
8: - permuted and scaled system: Cz = f, where
9: C = P D^{-1/2} A D^{-1/2}
10: z = P D^{1/2} x
11: f = P D^{-1/2} b
12: D = diagonal of A
13: P = permutation matrix determined by coloring
14: In this case, we use pre-solve and post-solve phases to handle scaling and
15: permutation, and by default the scaled residual norm is monitored for the
16: ILU/ICC preconditioners. Use the option
17: -ksp_truemonitor
18: to print both the scaled and unscaled residual norms.
19: */
21: #include petsc.h
23: #include src/mat/impls/rowbs/mpi/mpirowbs.h
27: PetscErrorCode MatScaleSystem_MPIRowbs(Mat mat,Vec x,Vec rhs)
28: {
29: Mat_MPIRowbs *bsif = (Mat_MPIRowbs*)mat->data;
30: Vec v = bsif->xwork;
31: PetscScalar *xa,*rhsa,*va;
35: /* Permute and scale RHS and solution vectors */
36: if (x) {
37: VecGetArray(x,&xa);
38: VecGetArray(v,&va);
39: BSperm_dvec(xa,va,bsif->pA->perm);CHKERRBS(0);
40: VecRestoreArray(x,&xa);
41: VecRestoreArray(v,&va);
42: VecPointwiseDivide(v,bsif->diag,x);
43: }
45: if (rhs) {
46: VecGetArray(rhs,&rhsa);
47: VecGetArray(v,&va);
48: BSperm_dvec(rhsa,va,bsif->pA->perm);CHKERRBS(0);
49: VecRestoreArray(rhs,&rhsa);
50: VecRestoreArray(v,&va);
51: VecPointwiseMult(v,bsif->diag,rhs);
52: }
53: return(0);
54: }
58: PetscErrorCode MatUnScaleSystem_MPIRowbs(Mat mat,Vec x,Vec rhs)
59: {
60: Mat_MPIRowbs *bsif = (Mat_MPIRowbs*)mat->data;
61: Vec v = bsif->xwork;
62: PetscScalar *xa,*va,*rhsa;
66: /* Unpermute and unscale the solution and RHS vectors */
67: if (x) {
68: VecPointwiseMult(x,bsif->diag,v);
69: VecGetArray(v,&va);
70: VecGetArray(x,&xa);
71: BSiperm_dvec(va,xa,bsif->pA->perm);CHKERRBS(0);
72: VecRestoreArray(x,&xa);
73: VecRestoreArray(v,&va);
74: }
75: if (rhs) {
76: VecPointwiseDivide(rhs,bsif->diag,v);
77: VecGetArray(rhs,&rhsa);
78: VecGetArray(v,&va);
79: BSiperm_dvec(va,rhsa,bsif->pA->perm);CHKERRBS(0);
80: VecRestoreArray(rhs,&rhsa);
81: VecRestoreArray(v,&va);
82: }
83: return(0);
84: }
88: PetscErrorCode MatUseScaledForm_MPIRowbs(Mat mat,PetscTruth scale)
89: {
90: Mat_MPIRowbs *bsif = (Mat_MPIRowbs*)mat->data;
93: bsif->vecs_permscale = scale;
94: return(0);
95: }