Actual source code: sacusppoly.cu

petsc-3.4.2 2013-07-02
  2: /*  -------------------------------------------------------------------- */

  4: /*
  5:    Include files needed for the CUSP Smoothed Aggregation preconditioner with Chebyshev polynomial smoothing:
  6:      pcimpl.h - private include file intended for use by all preconditioners
  7: */

  9: #include <petsc-private/pcimpl.h>   /*I "petscpc.h" I*/
 10: #include <../src/mat/impls/aij/seq/aij.h>
 11: #include <cusp/monitor.h>
 12: #define USE_POLY_SMOOTHER 1
 13: #include <cusp/precond/smoothed_aggregation.h>
 14: #undef USE_POLY_SMOOTHER
 15: #include <../src/vec/vec/impls/dvecimpl.h>
 16: #include <../src/mat/impls/aij/seq/seqcusp/cuspmatimpl.h>

 18: #define cuspsaprecond cusp::precond::smoothed_aggregation<PetscInt,PetscScalar,cusp::device_memory>

 20: /*
 21:    Private context (data structure) for the SACUSPPoly preconditioner.
 22: */
 23: typedef struct {
 24:   cuspsaprecond * SACUSPPoly;
 25:   /*int cycles; */
 26: } PC_SACUSPPoly;


 29: /* -------------------------------------------------------------------------- */
 30: /*
 31:    PCSetUp_SACUSPPoly - Prepares for the use of the SACUSPPoly preconditioner
 32:                     by setting data structures and options.

 34:    Input Parameter:
 35: .  pc - the preconditioner context

 37:    Application Interface Routine: PCSetUp()

 39:    Notes:
 40:    The interface routine PCSetUp() is not usually called directly by
 41:    the user, but instead is called by PCApply() if necessary.
 42: */
 45: static PetscErrorCode PCSetUp_SACUSPPoly(PC pc)
 46: {
 47:   PC_SACUSPPoly  *sa = (PC_SACUSPPoly*)pc->data;
 48:   PetscBool      flg = PETSC_FALSE;
 50: #if !defined(PETSC_USE_COMPLEX)
 51:   // protect these in order to avoid compiler warnings. This preconditioner does
 52:   // not work for complex types.
 53:   Mat_SeqAIJCUSP *gpustruct;
 54:   CUSPMATRIX     *mat;
 55: #endif
 57:   PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJCUSP,&flg);
 58:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Currently only handles CUSP matrices");
 59:   if (pc->setupcalled != 0) {
 60:     try {
 61:       delete sa->SACUSPPoly;
 62:     } catch(char *ex) {
 63:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
 64:     }
 65:   }
 66:   try {
 67: #if defined(PETSC_USE_COMPLEX)
 68:     sa->SACUSPPoly = 0;CHKERRQ(1); /* TODO */
 69: #else
 70:     MatCUSPCopyToGPU(pc->pmat);
 71:     gpustruct = (Mat_SeqAIJCUSP*)(pc->pmat->spptr);
 72: #if defined(PETSC_HAVE_TXPETSCGPU)
 73:     gpustruct->mat->getCsrMatrix(&mat);CHKERRCUSP(ierr);
 74: #else
 75:     mat = (CUSPMATRIX*)gpustruct->mat;
 76: #endif

 78:     sa->SACUSPPoly = new cuspsaprecond(*mat);
 79: #endif
 80:   } catch(char *ex) {
 81:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
 82:   }
 83:   /*PetscOptionsInt("-pc_sacusp_cycles","Number of v-cycles to perform","PCSACUSPSetCycles",sa->cycles,
 84:     &sa->cycles,NULL);*/
 85:   return(0);
 86: }

 90: static PetscErrorCode PCApplyRichardson_SACUSPPoly(PC pc, Vec b, Vec y, Vec w,PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt its, PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
 91: {
 92: #if !defined(PETSC_USE_COMPLEX)
 93:   // protect these in order to avoid compiler warnings. This preconditioner does
 94:   // not work for complex types.
 95:   PC_SACUSPPoly *sac = (PC_SACUSPPoly*)pc->data;
 96: #endif
 98:   CUSPARRAY      *barray,*yarray;

101:   /* how to incorporate dtol, guesszero, w?*/
103:   VecCUSPGetArrayRead(b,&barray);
104:   VecCUSPGetArrayReadWrite(y,&yarray);
105:   cusp::default_monitor<PetscReal> monitor(*barray,its,rtol,abstol);
106: #if defined(PETSC_USE_COMPLEX)
107:   CHKERRQ(1); /* TODO */
108: #else
109:   sac->SACUSPPoly->solve(*barray,*yarray,monitor);
110: #endif
111:   *outits = monitor.iteration_count();
112:   if (monitor.converged()) *reason = PCRICHARDSON_CONVERGED_RTOL; /* how to discern between converging from RTOL or ATOL?*/
113:   else *reason = PCRICHARDSON_CONVERGED_ITS;

115:   PetscObjectStateIncrease((PetscObject)y);
116:   VecCUSPRestoreArrayRead(b,&barray);
117:   VecCUSPRestoreArrayReadWrite(y,&yarray);
118:   return(0);
119: }

121: /* -------------------------------------------------------------------------- */
122: /*
123:    PCApply_SACUSPPoly - Applies the SACUSPPoly preconditioner to a vector.

125:    Input Parameters:
126: .  pc - the preconditioner context
127: .  x - input vector

129:    Output Parameter:
130: .  y - output vector

132:    Application Interface Routine: PCApply()
133:  */
136: static PetscErrorCode PCApply_SACUSPPoly(PC pc,Vec x,Vec y)
137: {
138:   PC_SACUSPPoly  *sac = (PC_SACUSPPoly*)pc->data;
140:   PetscBool      flg1,flg2;
141:   CUSPARRAY      *xarray=NULL,*yarray=NULL;

144:   /*how to apply a certain fixed number of iterations?*/
145:   PetscObjectTypeCompare((PetscObject)x,VECSEQCUSP,&flg1);
146:   PetscObjectTypeCompare((PetscObject)y,VECSEQCUSP,&flg2);
147:   if (!(flg1 && flg2)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP, "Currently only handles CUSP vectors");
148:   if (!sac->SACUSPPoly) {
149:     PCSetUp_SACUSPPoly(pc);
150:   }
151:   VecSet(y,0.0);
152:   VecCUSPGetArrayRead(x,&xarray);
153:   VecCUSPGetArrayWrite(y,&yarray);
154:   try {
155: #if defined(PETSC_USE_COMPLEX)
156:     CHKERRQ(1); /* TODO */
157: #else
158:     cusp::multiply(*sac->SACUSPPoly,*xarray,*yarray);
159: #endif
160:   } catch(char *ex) {
161:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
162:   }
163:   VecCUSPRestoreArrayRead(x,&xarray);
164:   VecCUSPRestoreArrayWrite(y,&yarray);
165:   PetscObjectStateIncrease((PetscObject)y);
166:   return(0);
167: }
168: /* -------------------------------------------------------------------------- */
169: /*
170:    PCDestroy_SACUSPPoly - Destroys the private context for the SACUSPPoly preconditioner
171:    that was created with PCCreate_SACUSPPoly().

173:    Input Parameter:
174: .  pc - the preconditioner context

176:    Application Interface Routine: PCDestroy()
177: */
180: static PetscErrorCode PCDestroy_SACUSPPoly(PC pc)
181: {
182:   PC_SACUSPPoly  *sac = (PC_SACUSPPoly*)pc->data;

186:   if (sac->SACUSPPoly) {
187:     try {
188:       delete sac->SACUSPPoly;
189:     } catch(char *ex) {
190:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
191:     }
192:   }

194:   /*
195:       Free the private data structure that was hanging off the PC
196:   */
197:   PetscFree(sac);
198:   return(0);
199: }

203: static PetscErrorCode PCSetFromOptions_SACUSPPoly(PC pc)
204: {

208:   PetscOptionsHead("SACUSPPoly options");
209:   PetscOptionsTail();
210:   return(0);
211: }

213: /* -------------------------------------------------------------------------- */

217: PETSC_EXTERN PetscErrorCode PCCreate_SACUSPPoly(PC pc)
218: {
219:   PC_SACUSPPoly  *sac;

223:   /*
224:      Creates the private data structure for this preconditioner and
225:      attach it to the PC object.
226:   */
227:   PetscNewLog(pc,PC_SACUSPPoly,&sac);
228:   pc->data = (void*)sac;

230:   /*
231:      Initialize the pointer to zero
232:      Initialize number of v-cycles to default (1)
233:   */
234:   sac->SACUSPPoly = 0;
235:   /*sac->cycles=1;*/


238:   /*
239:       Set the pointers for the functions that are provided above.
240:       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
241:       are called, they will automatically call these functions.  Note we
242:       choose not to provide a couple of these functions since they are
243:       not needed.
244:   */
245:   pc->ops->apply               = PCApply_SACUSPPoly;
246:   pc->ops->applytranspose      = 0;
247:   pc->ops->setup               = PCSetUp_SACUSPPoly;
248:   pc->ops->destroy             = PCDestroy_SACUSPPoly;
249:   pc->ops->setfromoptions      = PCSetFromOptions_SACUSPPoly;
250:   pc->ops->view                = 0;
251:   pc->ops->applyrichardson     = PCApplyRichardson_SACUSPPoly;
252:   pc->ops->applysymmetricleft  = 0;
253:   pc->ops->applysymmetricright = 0;
254:   return(0);
255: }