Actual source code: sacusppoly.cu
petsc-3.8.3 2017-12-09
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: */
8: #define PETSC_SKIP_SPINLOCK
9: #include <petsc/private/pcimpl.h>
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>
24: #include <../src/vec/vec/impls/seq/seqcusp/cuspvecimpl.h>
26: /*
27: Private context (data structure) for the SACUSPPoly preconditioner.
28: */
29: typedef struct {
30: cuspsaprecond * SACUSPPoly;
31: /*int cycles; */
32: } PC_SACUSPPoly;
35: /* -------------------------------------------------------------------------- */
36: /*
37: PCSetUp_SACUSPPoly - Prepares for the use of the SACUSPPoly preconditioner
38: by setting data structures and options.
40: Input Parameter:
41: . pc - the preconditioner context
43: Application Interface Routine: PCSetUp()
45: Notes:
46: The interface routine PCSetUp() is not usually called directly by
47: the user, but instead is called by PCApply() if necessary.
48: */
49: static PetscErrorCode PCSetUp_SACUSPPoly(PC pc)
50: {
51: PC_SACUSPPoly *sa = (PC_SACUSPPoly*)pc->data;
52: PetscBool flg = PETSC_FALSE;
54: #if !defined(PETSC_USE_COMPLEX)
55: // protect these in order to avoid compiler warnings. This preconditioner does
56: // not work for complex types.
57: Mat_SeqAIJCUSP *gpustruct;
58: #endif
60: PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJCUSP,&flg);
61: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Currently only handles CUSP matrices");
62: if (pc->setupcalled != 0) {
63: try {
64: delete sa->SACUSPPoly;
65: } catch(char *ex) {
66: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
67: }
68: }
69: try {
70: #if defined(PETSC_USE_COMPLEX)
71: sa->SACUSPPoly = 0;CHKERRQ(1); /* TODO */
72: #else
73: MatCUSPCopyToGPU(pc->pmat);
74: gpustruct = (Mat_SeqAIJCUSP*)(pc->pmat->spptr);
75:
76: if (gpustruct->format==MAT_CUSP_ELL) {
77: CUSPMATRIXELL *mat = (CUSPMATRIXELL*)gpustruct->mat;
78: sa->SACUSPPoly = new cuspsaprecond(*mat);
79: } else if (gpustruct->format==MAT_CUSP_DIA) {
80: CUSPMATRIXDIA *mat = (CUSPMATRIXDIA*)gpustruct->mat;
81: sa->SACUSPPoly = new cuspsaprecond(*mat);
82: } else {
83: CUSPMATRIX *mat = (CUSPMATRIX*)gpustruct->mat;
84: sa->SACUSPPoly = new cuspsaprecond(*mat);
85: }
86: #endif
87: } catch(char *ex) {
88: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
89: }
90: /*PetscOptionsInt("-pc_sacusp_cycles","Number of v-cycles to perform","PCSACUSPSetCycles",sa->cycles,
91: &sa->cycles,NULL);*/
92: return(0);
93: }
95: 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)
96: {
97: #if !defined(PETSC_USE_COMPLEX)
98: // protect these in order to avoid compiler warnings. This preconditioner does
99: // not work for complex types.
100: PC_SACUSPPoly *sac = (PC_SACUSPPoly*)pc->data;
101: #endif
103: CUSPARRAY *barray,*yarray;
106: /* how to incorporate dtol, guesszero, w?*/
108: VecCUSPGetArrayRead(b,&barray);
109: VecCUSPGetArrayReadWrite(y,&yarray);
110: #if defined(CUSP_VERSION) && CUSP_VERSION >= 500
111: cusp::monitor<PetscReal> monitor(*barray,its,rtol,abstol);
112: #else
113: cusp::default_monitor<PetscReal> monitor(*barray,its,rtol,abstol);
114: #endif
115: #if defined(PETSC_USE_COMPLEX)
116: CHKERRQ(1); /* TODO */
117: #else
118: sac->SACUSPPoly->solve(*barray,*yarray,monitor);
119: #endif
120: *outits = monitor.iteration_count();
121: if (monitor.converged()) *reason = PCRICHARDSON_CONVERGED_RTOL; /* how to discern between converging from RTOL or ATOL?*/
122: else *reason = PCRICHARDSON_CONVERGED_ITS;
124: PetscObjectStateIncrease((PetscObject)y);
125: VecCUSPRestoreArrayRead(b,&barray);
126: VecCUSPRestoreArrayReadWrite(y,&yarray);
127: return(0);
128: }
130: /* -------------------------------------------------------------------------- */
131: /*
132: PCApply_SACUSPPoly - Applies the SACUSPPoly preconditioner to a vector.
134: Input Parameters:
135: . pc - the preconditioner context
136: . x - input vector
138: Output Parameter:
139: . y - output vector
141: Application Interface Routine: PCApply()
142: */
143: static PetscErrorCode PCApply_SACUSPPoly(PC pc,Vec x,Vec y)
144: {
145: PC_SACUSPPoly *sac = (PC_SACUSPPoly*)pc->data;
147: PetscBool flg1,flg2;
148: CUSPARRAY *xarray=NULL,*yarray=NULL;
151: /*how to apply a certain fixed number of iterations?*/
152: PetscObjectTypeCompare((PetscObject)x,VECSEQCUSP,&flg1);
153: PetscObjectTypeCompare((PetscObject)y,VECSEQCUSP,&flg2);
154: if (!(flg1 && flg2)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP, "Currently only handles CUSP vectors");
155: if (!sac->SACUSPPoly) {
156: PCSetUp_SACUSPPoly(pc);
157: }
158: VecSet(y,0.0);
159: VecCUSPGetArrayRead(x,&xarray);
160: VecCUSPGetArrayWrite(y,&yarray);
161: try {
162: #if defined(PETSC_USE_COMPLEX)
163: CHKERRQ(1); /* TODO */
164: #else
165: cusp::multiply(*sac->SACUSPPoly,*xarray,*yarray);
166: #endif
167: } catch(char *ex) {
168: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
169: }
170: VecCUSPRestoreArrayRead(x,&xarray);
171: VecCUSPRestoreArrayWrite(y,&yarray);
172: PetscObjectStateIncrease((PetscObject)y);
173: return(0);
174: }
175: /* -------------------------------------------------------------------------- */
176: /*
177: PCDestroy_SACUSPPoly - Destroys the private context for the SACUSPPoly preconditioner
178: that was created with PCCreate_SACUSPPoly().
180: Input Parameter:
181: . pc - the preconditioner context
183: Application Interface Routine: PCDestroy()
184: */
185: static PetscErrorCode PCDestroy_SACUSPPoly(PC pc)
186: {
187: PC_SACUSPPoly *sac = (PC_SACUSPPoly*)pc->data;
191: if (sac->SACUSPPoly) {
192: try {
193: delete sac->SACUSPPoly;
194: } catch(char *ex) {
195: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
196: }
197: }
199: /*
200: Free the private data structure that was hanging off the PC
201: */
202: PetscFree(sac);
203: return(0);
204: }
206: static PetscErrorCode PCSetFromOptions_SACUSPPoly(PetscOptionItems *PetscOptionsObject,PC pc)
207: {
211: PetscOptionsHead(PetscOptionsObject,"SACUSPPoly options");
212: PetscOptionsTail();
213: return(0);
214: }
216: /* -------------------------------------------------------------------------- */
218: PETSC_EXTERN PetscErrorCode PCCreate_SACUSPPoly(PC pc)
219: {
220: PC_SACUSPPoly *sac;
224: /*
225: Creates the private data structure for this preconditioner and
226: attach it to the PC object.
227: */
228: PetscNewLog(pc,&sac);
229: pc->data = (void*)sac;
231: /*
232: Initialize the pointer to zero
233: Initialize number of v-cycles to default (1)
234: */
235: sac->SACUSPPoly = 0;
236: /*sac->cycles=1;*/
239: /*
240: Set the pointers for the functions that are provided above.
241: Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
242: are called, they will automatically call these functions. Note we
243: choose not to provide a couple of these functions since they are
244: not needed.
245: */
246: pc->ops->apply = PCApply_SACUSPPoly;
247: pc->ops->applytranspose = 0;
248: pc->ops->setup = PCSetUp_SACUSPPoly;
249: pc->ops->destroy = PCDestroy_SACUSPPoly;
250: pc->ops->setfromoptions = PCSetFromOptions_SACUSPPoly;
251: pc->ops->view = 0;
252: pc->ops->applyrichardson = PCApplyRichardson_SACUSPPoly;
253: pc->ops->applysymmetricleft = 0;
254: pc->ops->applysymmetricright = 0;
255: return(0);
256: }