Actual source code: sacusp.cu
petsc-3.4.2 2013-07-02
2: /* -------------------------------------------------------------------- */
4: /*
5: Include files needed for the CUSP Smoothed Aggregation preconditioner:
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/precond/smoothed_aggregation.h>
13: #include <../src/vec/vec/impls/dvecimpl.h>
14: #include <../src/mat/impls/aij/seq/seqcusp/cuspmatimpl.h>
16: #define cuspsaprecond cusp::precond::smoothed_aggregation<PetscInt,PetscScalar,cusp::device_memory>
18: /*
19: Private context (data structure) for the SACUSP preconditioner.
20: */
21: typedef struct {
22: cuspsaprecond * SACUSP;
23: /*int cycles; */
24: } PC_SACUSP;
28: static PetscErrorCode PCSACUSPSetCycles(PC pc, int n)
29: {
30: PC_SACUSP *sac = (PC_SACUSP*)pc->data;
33: sac->cycles = n;
34: return(0);
36: }*/
38: /* -------------------------------------------------------------------------- */
39: /*
40: PCSetUp_SACUSP - Prepares for the use of the SACUSP preconditioner
41: by setting data structures and options.
43: Input Parameter:
44: . pc - the preconditioner context
46: Application Interface Routine: PCSetUp()
48: Notes:
49: The interface routine PCSetUp() is not usually called directly by
50: the user, but instead is called by PCApply() if necessary.
51: */
54: static PetscErrorCode PCSetUp_SACUSP(PC pc)
55: {
56: PC_SACUSP *sa = (PC_SACUSP*)pc->data;
57: PetscBool flg = PETSC_FALSE;
59: #if !defined(PETSC_USE_COMPLEX)
60: // protect these in order to avoid compiler warnings. This preconditioner does
61: // not work for complex types.
62: Mat_SeqAIJCUSP *gpustruct;
63: CUSPMATRIX *mat;
64: #endif
67: PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJCUSP,&flg);
68: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Currently only handles CUSP matrices");
69: if (pc->setupcalled != 0) {
70: try {
71: delete sa->SACUSP;
72: } catch(char *ex) {
73: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
74: }
75: }
76: try {
77: #if defined(PETSC_USE_COMPLEX)
78: sa->SACUSP = 0;CHKERRQ(1); /* TODO */
79: #else
80: MatCUSPCopyToGPU(pc->pmat);
81: gpustruct = (Mat_SeqAIJCUSP*)(pc->pmat->spptr);
82: #if defined(PETSC_HAVE_TXPETSCGPU)
83: gpustruct->mat->getCsrMatrix(&mat);CHKERRCUSP(ierr);
84: #else
85: mat = (CUSPMATRIX*)gpustruct->mat;
86: #endif
87: sa->SACUSP = new cuspsaprecond(*mat);
88: #endif
90: } catch(char *ex) {
91: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
92: }
93: /*PetscOptionsInt("-pc_sacusp_cycles","Number of v-cycles to perform","PCSACUSPSetCycles",sa->cycles,
94: &sa->cycles,NULL);*/
95: return(0);
96: }
100: static PetscErrorCode PCApplyRichardson_SACUSP(PC pc, Vec b, Vec y, Vec w,PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt its, PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
101: {
102: #if !defined(PETSC_USE_COMPLEX)
103: // protect these in order to avoid compiler warnings. This preconditioner does
104: // not work for complex types.
105: PC_SACUSP *sac = (PC_SACUSP*)pc->data;
106: #endif
108: CUSPARRAY *barray,*yarray;
111: /* how to incorporate dtol, guesszero, w?*/
113: VecCUSPGetArrayRead(b,&barray);
114: VecCUSPGetArrayReadWrite(y,&yarray);
115: cusp::default_monitor<PetscReal> monitor(*barray,its,rtol,abstol);
116: #if defined(PETSC_USE_COMPLEX)
117: CHKERRQ(1);
118: /* TODO */
119: #else
120: sac->SACUSP->solve(*barray,*yarray,monitor);
121: *outits = monitor.iteration_count();
122: if (monitor.converged()) *reason = PCRICHARDSON_CONVERGED_RTOL; /* how to discern between converging from RTOL or ATOL?*/
123: else *reason = PCRICHARDSON_CONVERGED_ITS;
124: #endif
125: PetscObjectStateIncrease((PetscObject)y);
126: VecCUSPRestoreArrayRead(b,&barray);
127: VecCUSPRestoreArrayReadWrite(y,&yarray);
128: return(0);
129: }
131: /* -------------------------------------------------------------------------- */
132: /*
133: PCApply_SACUSP - Applies the SACUSP preconditioner to a vector.
135: Input Parameters:
136: . pc - the preconditioner context
137: . x - input vector
139: Output Parameter:
140: . y - output vector
142: Application Interface Routine: PCApply()
143: */
146: static PetscErrorCode PCApply_SACUSP(PC pc,Vec x,Vec y)
147: {
148: PC_SACUSP *sac = (PC_SACUSP*)pc->data;
150: PetscBool flg1,flg2;
151: CUSPARRAY *xarray=NULL,*yarray=NULL;
154: /*how to apply a certain fixed number of iterations?*/
155: PetscObjectTypeCompare((PetscObject)x,VECSEQCUSP,&flg1);
156: PetscObjectTypeCompare((PetscObject)y,VECSEQCUSP,&flg2);
157: if (!(flg1 && flg2)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP, "Currently only handles CUSP vectors");
158: if (!sac->SACUSP) {
159: PCSetUp_SACUSP(pc);
160: }
161: VecSet(y,0.0);
162: VecCUSPGetArrayRead(x,&xarray);
163: VecCUSPGetArrayWrite(y,&yarray);
164: try {
165: #if defined(PETSC_USE_COMPLEX)
167: #else
168: cusp::multiply(*sac->SACUSP,*xarray,*yarray);
169: #endif
170: } catch(char * ex) {
171: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
172: }
173: VecCUSPRestoreArrayRead(x,&xarray);
174: VecCUSPRestoreArrayWrite(y,&yarray);
175: PetscObjectStateIncrease((PetscObject)y);
176: return(0);
177: }
178: /* -------------------------------------------------------------------------- */
179: /*
180: PCDestroy_SACUSP - Destroys the private context for the SACUSP preconditioner
181: that was created with PCCreate_SACUSP().
183: Input Parameter:
184: . pc - the preconditioner context
186: Application Interface Routine: PCDestroy()
187: */
190: static PetscErrorCode PCDestroy_SACUSP(PC pc)
191: {
192: PC_SACUSP *sac = (PC_SACUSP*)pc->data;
196: if (sac->SACUSP) {
197: try {
198: delete sac->SACUSP;
199: } catch(char * ex) {
200: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
201: }
202: }
204: /*
205: Free the private data structure that was hanging off the PC
206: */
207: PetscFree(pc->data);
208: return(0);
209: }
213: static PetscErrorCode PCSetFromOptions_SACUSP(PC pc)
214: {
218: PetscOptionsHead("SACUSP options");
219: PetscOptionsTail();
220: return(0);
221: }
223: /* -------------------------------------------------------------------------- */
226: /*MC
227: PCSACUSP - A smoothed agglomeration algorithm that runs on the Nvidia GPU.
230: http://research.nvidia.com/sites/default/files/publications/nvr-2011-002.pdf
232: Level: advanced
234: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC
236: M*/
240: PETSC_EXTERN PetscErrorCode PCCreate_SACUSP(PC pc)
241: {
242: PC_SACUSP *sac;
246: /*
247: Creates the private data structure for this preconditioner and
248: attach it to the PC object.
249: */
250: PetscNewLog(pc,PC_SACUSP,&sac);
251: pc->data = (void*)sac;
253: /*
254: Initialize the pointer to zero
255: Initialize number of v-cycles to default (1)
256: */
257: sac->SACUSP = 0;
258: /*sac->cycles=1;*/
261: /*
262: Set the pointers for the functions that are provided above.
263: Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
264: are called, they will automatically call these functions. Note we
265: choose not to provide a couple of these functions since they are
266: not needed.
267: */
268: pc->ops->apply = PCApply_SACUSP;
269: pc->ops->applytranspose = 0;
270: pc->ops->setup = PCSetUp_SACUSP;
271: pc->ops->destroy = PCDestroy_SACUSP;
272: pc->ops->setfromoptions = PCSetFromOptions_SACUSP;
273: pc->ops->view = 0;
274: pc->ops->applyrichardson = PCApplyRichardson_SACUSP;
275: pc->ops->applysymmetricleft = 0;
276: pc->ops->applysymmetricright = 0;
277: return(0);
278: }