Actual source code: jacobic.c

  1: #ifdef PETSC_RCS_HEADER
  2: static char vcid[] = "$Id: jacobic.c,v 1.1 1999/05/20 21:56:24 knepley Exp $";
  3: #endif
  4: /*
  5:    Defines a Jacobi preconditioner for any constrained Mat implementation
  6: */
 7:  #include src/sles/pc/pcimpl.h
 8:  #include gvec.h

 10: /* We are cheating here since the interface is deficient */
 11: typedef struct {
 12:   Vec        diag;
 13:   Vec        diagsqrt;
 14:   PetscTruth userowmax;
 15: } PC_Jacobi;

 17: typedef struct {
 18:   PC jacobi;
 19: } PC_Jacobi_Constrained;

 21: #undef __FUNCT__  
 23: static int PCSetUp_Jacobi_Constrained(PC pc) {
 24:   PC_Jacobi_Constrained *jacConst     = (PC_Jacobi_Constrained *) pc->data;
 25:   PC_Jacobi             *jac          = (PC_Jacobi *) jacConst->jacobi->data;
 26:   PetscTruth             zeroDetected = PETSC_FALSE;
 27:   Grid                   grid;
 28:   Vec                    diag, diagsqrt;
 29:   PetscScalar           *x;
 30:   int                    i, n;
 31:   int                    ierr;

 34:   PCSetVector(jacConst->jacobi, pc->vec);
 35:   PCSetOperators(jacConst->jacobi, pc->mat, pc->pmat, pc->flag);
 36:   GMatGetGrid(pc->pmat, &grid);
 37:   if (grid == PETSC_NULL) SETERRQ(PETSC_ERR_SUP, "Constrained Jacobi is only for Grid matrices");

 39:   jacConst->jacobi->setupcalled = 2;
 40:   diag     = jac->diag;
 41:   diagsqrt = jac->diagsqrt;
 42:   if (jac->userowmax) SETERRQ(PETSC_ERR_SUP, "Not currently supported for constrained matrices");
 43:   if (diag) {
 44:     GMatGetDiagonalConstrained(pc->pmat, diag);
 45:     VecReciprocal(diag);
 46:     VecGetLocalSize(diag, &n);
 47:     VecGetArray(diag, &x);
 48:     for(i = 0; i < n; i++) {
 49:       if (x[i] == 0.0) {
 50:         x[i] = 1.0;
 51:         zeroDetected = PETSC_TRUE;
 52:       }
 53:     }
 54:     VecRestoreArray(diag, &x);
 55:   }
 56:   if (diagsqrt) {
 57:     GMatGetDiagonalConstrained(pc->pmat, diagsqrt);
 58:     VecGetLocalSize(diagsqrt, &n);
 59:     VecGetArray(diagsqrt, &x);
 60:     for(i = 0; i < n; i++) {
 61:       if (x[i] != 0.0) {
 62:         x[i] = 1.0/sqrt(PetscAbsScalar(x[i]));
 63:       } else {
 64:         x[i] = 1.0;
 65:         zeroDetected = PETSC_TRUE;
 66:       }
 67:     }
 68:     VecRestoreArray(diagsqrt, &x);
 69:   }

 71:   if (zeroDetected == PETSC_TRUE) {
 72:     PetscLogInfo(pc, "PCSetUp_Jacobi_Constrained:WARNING: Zero detected in diagonal while building Jacobi preconditionern");
 73:   }
 74:   return(0);
 75: }

 77: #undef __FUNCT__  
 79: static int PCDestroy_Jacobi_Constrained(PC pc) {
 80:   PC_Jacobi_Constrained *jac = (PC_Jacobi_Constrained *) pc->data;
 81:   int                    ierr;

 84:   PCDestroy(jac->jacobi);
 85:   PetscFree(jac);
 86:   return(0);
 87: }

 89: #undef __FUNCT__  
 91: static int PCApply_Jacobi_Constrained(PC pc, Vec x, Vec y) {
 92:   PC_Jacobi_Constrained *jac = (PC_Jacobi_Constrained *) pc->data;
 93:   PC_Jacobi             *j   = (PC_Jacobi *) jac->jacobi->data;
 94:   int                    ierr;

 96:   if (j->diag == PETSC_NULL) {
 97:     VecDuplicate(pc->vec, &j->diag);
 98:     PetscLogObjectParent(pc, j->diag);
 99:     PCSetUp_Jacobi_Constrained(pc);
100:   }
101:   return PCApply(jac->jacobi, x, y);
102: }

104: #undef __FUNCT__  
106: static int PCSetFromOptions_Jacobi_Constrained(PC pc) {
107:   PC_Jacobi_Constrained *jac = (PC_Jacobi_Constrained *) pc->data;
108:   return PCSetFromOptions(jac->jacobi);
109: }

111: #undef __FUNCT__  
113: static int PCApplySymmetricLeftOrRight_Jacobi_Constrained(PC pc, Vec x, Vec y) {
114:   PC_Jacobi_Constrained *jac = (PC_Jacobi_Constrained *) pc->data;
115:   PC_Jacobi             *j   = (PC_Jacobi *) jac->jacobi->data;
116:   int                    ierr;

118:   if (j->diagsqrt == PETSC_NULL) {
119:     VecDuplicate(pc->vec, &j->diagsqrt);
120:     PetscLogObjectParent(pc, j->diagsqrt);
121:     PCSetUp_Jacobi_Constrained(pc);
122:   }
123:   return PCApplySymmetricLeft(jac->jacobi, x, y);
124: }

126: EXTERN_C_BEGIN
127: #undef __FUNCT__  
129: int PCCreate_Jacobi_Constrained(PC pc) {
130:   PC_Jacobi_Constrained *jac;
131:   int                    ierr;

134:   PetscNew(PC_Jacobi_Constrained, &jac);
135:   PetscLogObjectMemory(pc,sizeof(PC_Jacobi_Constrained));
136:   pc->data = (void*) jac;

138:   PCCreate(pc->comm, &jac->jacobi);
139:   PCSetOptionsPrefix(jac->jacobi, "jacobic");
140:   PCSetType(jac->jacobi, PCJACOBI);
141:   pc->ops->apply               = PCApply_Jacobi_Constrained;
142:   pc->ops->applytranspose      = PCApply_Jacobi_Constrained;
143:   pc->ops->setup               = PCSetUp_Jacobi_Constrained;
144:   pc->ops->destroy             = PCDestroy_Jacobi_Constrained;
145:   pc->ops->setfromoptions      = PCSetFromOptions_Jacobi_Constrained;
146:   pc->ops->view                = PETSC_NULL;
147:   pc->ops->applyrichardson     = PETSC_NULL;
148:   pc->ops->applysymmetricleft  = PCApplySymmetricLeftOrRight_Jacobi_Constrained;
149:   pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi_Constrained;
150:   return(0);
151: }
152: EXTERN_C_END