Actual source code: sacusppoly.cu
petsc-dev 2014-02-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: #endif
56: PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJCUSP,&flg);
57: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Currently only handles CUSP matrices");
58: if (pc->setupcalled != 0) {
59: try {
60: delete sa->SACUSPPoly;
61: } catch(char *ex) {
62: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
63: }
64: }
65: try {
66: #if defined(PETSC_USE_COMPLEX)
67: sa->SACUSPPoly = 0;CHKERRQ(1); /* TODO */
68: #else
69: MatCUSPCopyToGPU(pc->pmat);
70: gpustruct = (Mat_SeqAIJCUSP*)(pc->pmat->spptr);
71:
72: if (gpustruct->format==MAT_CUSP_ELL) {
73: CUSPMATRIXELL *mat = (CUSPMATRIXELL*)gpustruct->mat;
74: sa->SACUSPPoly = new cuspsaprecond(*mat);
75: } else if (gpustruct->format==MAT_CUSP_DIA) {
76: CUSPMATRIXDIA *mat = (CUSPMATRIXDIA*)gpustruct->mat;
77: sa->SACUSPPoly = new cuspsaprecond(*mat);
78: } else {
79: CUSPMATRIX *mat = (CUSPMATRIX*)gpustruct->mat;
80: sa->SACUSPPoly = new cuspsaprecond(*mat);
81: }
82: #endif
83: } catch(char *ex) {
84: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
85: }
86: /*PetscOptionsInt("-pc_sacusp_cycles","Number of v-cycles to perform","PCSACUSPSetCycles",sa->cycles,
87: &sa->cycles,NULL);*/
88: return(0);
89: }
93: 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)
94: {
95: #if !defined(PETSC_USE_COMPLEX)
96: // protect these in order to avoid compiler warnings. This preconditioner does
97: // not work for complex types.
98: PC_SACUSPPoly *sac = (PC_SACUSPPoly*)pc->data;
99: #endif
101: CUSPARRAY *barray,*yarray;
104: /* how to incorporate dtol, guesszero, w?*/
106: VecCUSPGetArrayRead(b,&barray);
107: VecCUSPGetArrayReadWrite(y,&yarray);
108: cusp::default_monitor<PetscReal> monitor(*barray,its,rtol,abstol);
109: #if defined(PETSC_USE_COMPLEX)
110: CHKERRQ(1); /* TODO */
111: #else
112: sac->SACUSPPoly->solve(*barray,*yarray,monitor);
113: #endif
114: *outits = monitor.iteration_count();
115: if (monitor.converged()) *reason = PCRICHARDSON_CONVERGED_RTOL; /* how to discern between converging from RTOL or ATOL?*/
116: else *reason = PCRICHARDSON_CONVERGED_ITS;
118: PetscObjectStateIncrease((PetscObject)y);
119: VecCUSPRestoreArrayRead(b,&barray);
120: VecCUSPRestoreArrayReadWrite(y,&yarray);
121: return(0);
122: }
124: /* -------------------------------------------------------------------------- */
125: /*
126: PCApply_SACUSPPoly - Applies the SACUSPPoly preconditioner to a vector.
128: Input Parameters:
129: . pc - the preconditioner context
130: . x - input vector
132: Output Parameter:
133: . y - output vector
135: Application Interface Routine: PCApply()
136: */
139: static PetscErrorCode PCApply_SACUSPPoly(PC pc,Vec x,Vec y)
140: {
141: PC_SACUSPPoly *sac = (PC_SACUSPPoly*)pc->data;
143: PetscBool flg1,flg2;
144: CUSPARRAY *xarray=NULL,*yarray=NULL;
147: /*how to apply a certain fixed number of iterations?*/
148: PetscObjectTypeCompare((PetscObject)x,VECSEQCUSP,&flg1);
149: PetscObjectTypeCompare((PetscObject)y,VECSEQCUSP,&flg2);
150: if (!(flg1 && flg2)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP, "Currently only handles CUSP vectors");
151: if (!sac->SACUSPPoly) {
152: PCSetUp_SACUSPPoly(pc);
153: }
154: VecSet(y,0.0);
155: VecCUSPGetArrayRead(x,&xarray);
156: VecCUSPGetArrayWrite(y,&yarray);
157: try {
158: #if defined(PETSC_USE_COMPLEX)
159: CHKERRQ(1); /* TODO */
160: #else
161: cusp::multiply(*sac->SACUSPPoly,*xarray,*yarray);
162: #endif
163: } catch(char *ex) {
164: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
165: }
166: VecCUSPRestoreArrayRead(x,&xarray);
167: VecCUSPRestoreArrayWrite(y,&yarray);
168: PetscObjectStateIncrease((PetscObject)y);
169: return(0);
170: }
171: /* -------------------------------------------------------------------------- */
172: /*
173: PCDestroy_SACUSPPoly - Destroys the private context for the SACUSPPoly preconditioner
174: that was created with PCCreate_SACUSPPoly().
176: Input Parameter:
177: . pc - the preconditioner context
179: Application Interface Routine: PCDestroy()
180: */
183: static PetscErrorCode PCDestroy_SACUSPPoly(PC pc)
184: {
185: PC_SACUSPPoly *sac = (PC_SACUSPPoly*)pc->data;
189: if (sac->SACUSPPoly) {
190: try {
191: delete sac->SACUSPPoly;
192: } catch(char *ex) {
193: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
194: }
195: }
197: /*
198: Free the private data structure that was hanging off the PC
199: */
200: PetscFree(sac);
201: return(0);
202: }
206: static PetscErrorCode PCSetFromOptions_SACUSPPoly(PC pc)
207: {
211: PetscOptionsHead("SACUSPPoly options");
212: PetscOptionsTail();
213: return(0);
214: }
216: /* -------------------------------------------------------------------------- */
220: PETSC_EXTERN PetscErrorCode PCCreate_SACUSPPoly(PC pc)
221: {
222: PC_SACUSPPoly *sac;
226: /*
227: Creates the private data structure for this preconditioner and
228: attach it to the PC object.
229: */
230: PetscNewLog(pc,&sac);
231: pc->data = (void*)sac;
233: /*
234: Initialize the pointer to zero
235: Initialize number of v-cycles to default (1)
236: */
237: sac->SACUSPPoly = 0;
238: /*sac->cycles=1;*/
241: /*
242: Set the pointers for the functions that are provided above.
243: Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
244: are called, they will automatically call these functions. Note we
245: choose not to provide a couple of these functions since they are
246: not needed.
247: */
248: pc->ops->apply = PCApply_SACUSPPoly;
249: pc->ops->applytranspose = 0;
250: pc->ops->setup = PCSetUp_SACUSPPoly;
251: pc->ops->destroy = PCDestroy_SACUSPPoly;
252: pc->ops->setfromoptions = PCSetFromOptions_SACUSPPoly;
253: pc->ops->view = 0;
254: pc->ops->applyrichardson = PCApplyRichardson_SACUSPPoly;
255: pc->ops->applysymmetricleft = 0;
256: pc->ops->applysymmetricright = 0;
257: return(0);
258: }