Actual source code: iccbs.c

  1: /*$Id: iccbs.c,v 1.45 2001/08/07 03:02:51 balay Exp $*/
  2: /*
  3:    Defines a Cholesky factorization preconditioner with BlockSolve95 interface.

  5:    Note that BlockSolve95 works with a scaled and permuted preconditioning matrix.
  6:    If the linear system matrix and preconditioning matrix are the same, we then
  7:    work directly with the permuted and scaled linear system:
  8:       - original system:  Ax = b
  9:       - permuted and scaled system:   Cz = f, where
 10:              C = P D^{-1/2} A D^{-1/2}
 11:              z = P D^{1/2} x
 12:              f = P D^{-1/2} b
 13:              D = diagonal of A
 14:              P = permutation matrix determined by coloring
 15:    In this case, we use pre-solve and post-solve phases to handle scaling and
 16:    permutation, and by default the scaled residual norm is monitored for the
 17:    ILU/ICC preconditioners.  Use the option
 18:      -ksp_truemonitor
 19:    to print both the scaled and unscaled residual norms.
 20: */

 22:  #include petsc.h

 24:  #include src/mat/impls/rowbs/mpi/mpirowbs.h

 28: int MatScaleSystem_MPIRowbs(Mat mat,Vec x,Vec rhs)
 29: {
 30:   Mat_MPIRowbs *bsif  = (Mat_MPIRowbs*)mat->data;
 31:   Vec          v = bsif->xwork;
 32:   PetscScalar  *xa,*rhsa,*va;
 33:   int          ierr;

 36:   /* Permute and scale RHS and solution vectors */
 37:   if (x) {
 38:     VecGetArray(x,&xa);
 39:     VecGetArray(v,&va);
 40:     BSperm_dvec(xa,va,bsif->pA->perm);CHKERRBS(0);
 41:     VecRestoreArray(x,&xa);
 42:     VecRestoreArray(v,&va);
 43:     VecPointwiseDivide(v,bsif->diag,x);
 44:   }

 46:   if (rhs) {
 47:     VecGetArray(rhs,&rhsa);
 48:     VecGetArray(v,&va);
 49:     BSperm_dvec(rhsa,va,bsif->pA->perm);CHKERRBS(0);
 50:     VecRestoreArray(rhs,&rhsa);
 51:     VecRestoreArray(v,&va);
 52:     VecPointwiseMult(v,bsif->diag,rhs);
 53:   }
 54:   return(0);
 55: }

 59: int MatUnScaleSystem_MPIRowbs(Mat mat,Vec x,Vec rhs)
 60: {
 61:   Mat_MPIRowbs *bsif  = (Mat_MPIRowbs*)mat->data;
 62:   Vec          v = bsif->xwork;
 63:   PetscScalar  *xa,*va,*rhsa;
 64:   int          ierr;

 67:   /* Unpermute and unscale the solution and RHS vectors */
 68:   if (x) {
 69:     VecPointwiseMult(x,bsif->diag,v);
 70:     VecGetArray(v,&va);
 71:     VecGetArray(x,&xa);
 72:     BSiperm_dvec(va,xa,bsif->pA->perm);CHKERRBS(0);
 73:     VecRestoreArray(x,&xa);
 74:     VecRestoreArray(v,&va);
 75:   }
 76:   if (rhs) {
 77:     VecPointwiseDivide(rhs,bsif->diag,v);
 78:     VecGetArray(rhs,&rhsa);
 79:     VecGetArray(v,&va);
 80:     BSiperm_dvec(va,rhsa,bsif->pA->perm);CHKERRBS(0);
 81:     VecRestoreArray(rhs,&rhsa);
 82:     VecRestoreArray(v,&va);
 83:   }
 84:   return(0);
 85: }

 89: int MatUseScaledForm_MPIRowbs(Mat mat,PetscTruth scale)
 90: {
 91:   Mat_MPIRowbs *bsif  = (Mat_MPIRowbs*)mat->data;

 94:   bsif->vecs_permscale = scale;
 95:   return(0);
 96: }