Actual source code: sacusppoly.cu
petsc-3.5.4 2015-05-23
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: #include <cusp/version.h>
13: #define USE_POLY_SMOOTHER 1
14: #if CUSP_VERSION >= 400
15: #include <cusp/precond/aggregation/smoothed_aggregation.h>
16: #define cuspsaprecond cusp::precond::aggregation::smoothed_aggregation<PetscInt,PetscScalar,cusp::device_memory>
17: #else
18: #include <cusp/precond/smoothed_aggregation.h>
19: #define cuspsaprecond cusp::precond::smoothed_aggregation<PetscInt,PetscScalar,cusp::device_memory>
20: #endif
21: #undef USE_POLY_SMOOTHER
22: #include <../src/vec/vec/impls/dvecimpl.h>
23: #include <../src/mat/impls/aij/seq/seqcusp/cuspmatimpl.h>
25: /*
26: Private context (data structure) for the SACUSPPoly preconditioner.
27: */
28: typedef struct {
29: cuspsaprecond * SACUSPPoly;
30: /*int cycles; */
31: } PC_SACUSPPoly;
34: /* -------------------------------------------------------------------------- */
35: /*
36: PCSetUp_SACUSPPoly - Prepares for the use of the SACUSPPoly preconditioner
37: by setting data structures and options.
39: Input Parameter:
40: . pc - the preconditioner context
42: Application Interface Routine: PCSetUp()
44: Notes:
45: The interface routine PCSetUp() is not usually called directly by
46: the user, but instead is called by PCApply() if necessary.
47: */
50: static PetscErrorCode PCSetUp_SACUSPPoly(PC pc)
51: {
52: PC_SACUSPPoly *sa = (PC_SACUSPPoly*)pc->data;
53: PetscBool flg = PETSC_FALSE;
55: #if !defined(PETSC_USE_COMPLEX)
56: // protect these in order to avoid compiler warnings. This preconditioner does
57: // not work for complex types.
58: Mat_SeqAIJCUSP *gpustruct;
59: #endif
61: PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJCUSP,&flg);
62: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Currently only handles CUSP matrices");
63: if (pc->setupcalled != 0) {
64: try {
65: delete sa->SACUSPPoly;
66: } catch(char *ex) {
67: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
68: }
69: }
70: try {
71: #if defined(PETSC_USE_COMPLEX)
72: sa->SACUSPPoly = 0;CHKERRQ(1); /* TODO */
73: #else
74: MatCUSPCopyToGPU(pc->pmat);
75: gpustruct = (Mat_SeqAIJCUSP*)(pc->pmat->spptr);
76:
77: if (gpustruct->format==MAT_CUSP_ELL) {
78: CUSPMATRIXELL *mat = (CUSPMATRIXELL*)gpustruct->mat;
79: sa->SACUSPPoly = new cuspsaprecond(*mat);
80: } else if (gpustruct->format==MAT_CUSP_DIA) {
81: CUSPMATRIXDIA *mat = (CUSPMATRIXDIA*)gpustruct->mat;
82: sa->SACUSPPoly = new cuspsaprecond(*mat);
83: } else {
84: CUSPMATRIX *mat = (CUSPMATRIX*)gpustruct->mat;
85: sa->SACUSPPoly = new cuspsaprecond(*mat);
86: }
87: #endif
88: } catch(char *ex) {
89: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
90: }
91: /*PetscOptionsInt("-pc_sacusp_cycles","Number of v-cycles to perform","PCSACUSPSetCycles",sa->cycles,
92: &sa->cycles,NULL);*/
93: return(0);
94: }
98: 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)
99: {
100: #if !defined(PETSC_USE_COMPLEX)
101: // protect these in order to avoid compiler warnings. This preconditioner does
102: // not work for complex types.
103: PC_SACUSPPoly *sac = (PC_SACUSPPoly*)pc->data;
104: #endif
106: CUSPARRAY *barray,*yarray;
109: /* how to incorporate dtol, guesszero, w?*/
111: VecCUSPGetArrayRead(b,&barray);
112: VecCUSPGetArrayReadWrite(y,&yarray);
113: cusp::default_monitor<PetscReal> monitor(*barray,its,rtol,abstol);
114: #if defined(PETSC_USE_COMPLEX)
115: CHKERRQ(1); /* TODO */
116: #else
117: sac->SACUSPPoly->solve(*barray,*yarray,monitor);
118: #endif
119: *outits = monitor.iteration_count();
120: if (monitor.converged()) *reason = PCRICHARDSON_CONVERGED_RTOL; /* how to discern between converging from RTOL or ATOL?*/
121: else *reason = PCRICHARDSON_CONVERGED_ITS;
123: PetscObjectStateIncrease((PetscObject)y);
124: VecCUSPRestoreArrayRead(b,&barray);
125: VecCUSPRestoreArrayReadWrite(y,&yarray);
126: return(0);
127: }
129: /* -------------------------------------------------------------------------- */
130: /*
131: PCApply_SACUSPPoly - Applies the SACUSPPoly preconditioner to a vector.
133: Input Parameters:
134: . pc - the preconditioner context
135: . x - input vector
137: Output Parameter:
138: . y - output vector
140: Application Interface Routine: PCApply()
141: */
144: static PetscErrorCode PCApply_SACUSPPoly(PC pc,Vec x,Vec y)
145: {
146: PC_SACUSPPoly *sac = (PC_SACUSPPoly*)pc->data;
148: PetscBool flg1,flg2;
149: CUSPARRAY *xarray=NULL,*yarray=NULL;
152: /*how to apply a certain fixed number of iterations?*/
153: PetscObjectTypeCompare((PetscObject)x,VECSEQCUSP,&flg1);
154: PetscObjectTypeCompare((PetscObject)y,VECSEQCUSP,&flg2);
155: if (!(flg1 && flg2)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP, "Currently only handles CUSP vectors");
156: if (!sac->SACUSPPoly) {
157: PCSetUp_SACUSPPoly(pc);
158: }
159: VecSet(y,0.0);
160: VecCUSPGetArrayRead(x,&xarray);
161: VecCUSPGetArrayWrite(y,&yarray);
162: try {
163: #if defined(PETSC_USE_COMPLEX)
164: CHKERRQ(1); /* TODO */
165: #else
166: cusp::multiply(*sac->SACUSPPoly,*xarray,*yarray);
167: #endif
168: } catch(char *ex) {
169: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
170: }
171: VecCUSPRestoreArrayRead(x,&xarray);
172: VecCUSPRestoreArrayWrite(y,&yarray);
173: PetscObjectStateIncrease((PetscObject)y);
174: return(0);
175: }
176: /* -------------------------------------------------------------------------- */
177: /*
178: PCDestroy_SACUSPPoly - Destroys the private context for the SACUSPPoly preconditioner
179: that was created with PCCreate_SACUSPPoly().
181: Input Parameter:
182: . pc - the preconditioner context
184: Application Interface Routine: PCDestroy()
185: */
188: static PetscErrorCode PCDestroy_SACUSPPoly(PC pc)
189: {
190: PC_SACUSPPoly *sac = (PC_SACUSPPoly*)pc->data;
194: if (sac->SACUSPPoly) {
195: try {
196: delete sac->SACUSPPoly;
197: } catch(char *ex) {
198: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
199: }
200: }
202: /*
203: Free the private data structure that was hanging off the PC
204: */
205: PetscFree(sac);
206: return(0);
207: }
211: static PetscErrorCode PCSetFromOptions_SACUSPPoly(PC pc)
212: {
216: PetscOptionsHead("SACUSPPoly options");
217: PetscOptionsTail();
218: return(0);
219: }
221: /* -------------------------------------------------------------------------- */
225: PETSC_EXTERN PetscErrorCode PCCreate_SACUSPPoly(PC pc)
226: {
227: PC_SACUSPPoly *sac;
231: /*
232: Creates the private data structure for this preconditioner and
233: attach it to the PC object.
234: */
235: PetscNewLog(pc,&sac);
236: pc->data = (void*)sac;
238: /*
239: Initialize the pointer to zero
240: Initialize number of v-cycles to default (1)
241: */
242: sac->SACUSPPoly = 0;
243: /*sac->cycles=1;*/
246: /*
247: Set the pointers for the functions that are provided above.
248: Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
249: are called, they will automatically call these functions. Note we
250: choose not to provide a couple of these functions since they are
251: not needed.
252: */
253: pc->ops->apply = PCApply_SACUSPPoly;
254: pc->ops->applytranspose = 0;
255: pc->ops->setup = PCSetUp_SACUSPPoly;
256: pc->ops->destroy = PCDestroy_SACUSPPoly;
257: pc->ops->setfromoptions = PCSetFromOptions_SACUSPPoly;
258: pc->ops->view = 0;
259: pc->ops->applyrichardson = PCApplyRichardson_SACUSPPoly;
260: pc->ops->applysymmetricleft = 0;
261: pc->ops->applysymmetricright = 0;
262: return(0);
263: }