xref: /petsc/src/ksp/pc/impls/saviennacl/saviennacl.cxx (revision 7efe37a1cedd385a2f501b843d47cdf14dfb49ea)
1 
2 /*
3    Include files needed for the ViennaCL Smoothed Aggregation preconditioner:
4      pcimpl.h - private include file intended for use by all preconditioners
5 */
6 #define PETSC_SKIP_SPINLOCK
7 #define PETSC_SKIP_IMMINTRIN_H_CUDAWORKAROUND 1
8 #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
9 #include <../src/mat/impls/aij/seq/aij.h>
10 #include <../src/vec/vec/impls/dvecimpl.h>
11 #include <../src/mat/impls/aij/seq/seqviennacl/viennaclmatimpl.h>
12 #include <../src/vec/vec/impls/seq/seqviennacl/viennaclvecimpl.h>
13 #include <viennacl/linalg/amg.hpp>
14 
15 /*
16    Private context (data structure) for the SAVIENNACL preconditioner.
17 */
18 typedef struct {
19   viennacl::linalg::amg_precond<viennacl::compressed_matrix<PetscScalar>> *SAVIENNACL;
20 } PC_SAVIENNACL;
21 
22 /*
23    PCSetUp_SAVIENNACL - Prepares for the use of the SAVIENNACL preconditioner
24                         by setting data structures and options.
25 
26    Input Parameter:
27 .  pc - the preconditioner context
28 
29    Application Interface Routine: PCSetUp()
30 
31    Note:
32    The interface routine PCSetUp() is not usually called directly by
33    the user, but instead is called by PCApply() if necessary.
34 */
35 static PetscErrorCode PCSetUp_SAVIENNACL(PC pc)
36 {
37   PC_SAVIENNACL      *sa  = (PC_SAVIENNACL *)pc->data;
38   PetscBool           flg = PETSC_FALSE;
39   Mat_SeqAIJViennaCL *gpustruct;
40 
41   PetscFunctionBegin;
42   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATSEQAIJVIENNACL, &flg));
43   PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Currently only handles ViennaCL matrices");
44   if (pc->setupcalled != 0) {
45     try {
46       delete sa->SAVIENNACL;
47     } catch (char *ex) {
48       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "ViennaCL error: %s", ex);
49     }
50   }
51   try {
52 #if defined(PETSC_USE_COMPLEX)
53     gpustruct = NULL;
54     SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "No support for complex arithmetic in SAVIENNACL preconditioner");
55 #else
56     PetscCall(MatViennaCLCopyToGPU(pc->pmat));
57     gpustruct = (Mat_SeqAIJViennaCL *)(pc->pmat->spptr);
58 
59     viennacl::linalg::amg_tag amg_tag_sa_pmis;
60     amg_tag_sa_pmis.set_coarsening_method(viennacl::linalg::AMG_COARSENING_METHOD_MIS2_AGGREGATION);
61     amg_tag_sa_pmis.set_interpolation_method(viennacl::linalg::AMG_INTERPOLATION_METHOD_SMOOTHED_AGGREGATION);
62     ViennaCLAIJMatrix *mat = (ViennaCLAIJMatrix *)gpustruct->mat;
63     sa->SAVIENNACL         = new viennacl::linalg::amg_precond<viennacl::compressed_matrix<PetscScalar>>(*mat, amg_tag_sa_pmis);
64     sa->SAVIENNACL->setup();
65 #endif
66   } catch (char *ex) {
67     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "ViennaCL error: %s", ex);
68   }
69   PetscFunctionReturn(0);
70 }
71 
72 /*
73    PCApply_SAVIENNACL - Applies the SAVIENNACL preconditioner to a vector.
74 
75    Input Parameters:
76 .  pc - the preconditioner context
77 .  x - input vector
78 
79    Output Parameter:
80 .  y - output vector
81 
82    Application Interface Routine: PCApply()
83  */
84 static PetscErrorCode PCApply_SAVIENNACL(PC pc, Vec x, Vec y)
85 {
86   PC_SAVIENNACL                       *sac = (PC_SAVIENNACL *)pc->data;
87   PetscBool                            flg1, flg2;
88   viennacl::vector<PetscScalar> const *xarray = NULL;
89   viennacl::vector<PetscScalar>       *yarray = NULL;
90 
91   PetscFunctionBegin;
92   /*how to apply a certain fixed number of iterations?*/
93   PetscCall(PetscObjectTypeCompare((PetscObject)x, VECSEQVIENNACL, &flg1));
94   PetscCall(PetscObjectTypeCompare((PetscObject)y, VECSEQVIENNACL, &flg2));
95   PetscCheck((flg1 && flg2), PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Currently only handles ViennaCL vectors");
96   if (!sac->SAVIENNACL) PetscCall(PCSetUp_SAVIENNACL(pc));
97   PetscCall(VecViennaCLGetArrayRead(x, &xarray));
98   PetscCall(VecViennaCLGetArrayWrite(y, &yarray));
99   try {
100 #if !defined(PETSC_USE_COMPLEX)
101     *yarray = *xarray;
102     sac->SAVIENNACL->apply(*yarray);
103 #endif
104   } catch (char *ex) {
105     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "ViennaCL error: %s", ex);
106   }
107   PetscCall(VecViennaCLRestoreArrayRead(x, &xarray));
108   PetscCall(VecViennaCLRestoreArrayWrite(y, &yarray));
109   PetscCall(PetscObjectStateIncrease((PetscObject)y));
110   PetscFunctionReturn(0);
111 }
112 
113 /*
114    PCDestroy_SAVIENNACL - Destroys the private context for the SAVIENNACL preconditioner
115    that was created with PCCreate_SAVIENNACL().
116 
117    Input Parameter:
118 .  pc - the preconditioner context
119 
120    Application Interface Routine: PCDestroy()
121 */
122 static PetscErrorCode PCDestroy_SAVIENNACL(PC pc)
123 {
124   PC_SAVIENNACL *sac = (PC_SAVIENNACL *)pc->data;
125 
126   PetscFunctionBegin;
127   if (sac->SAVIENNACL) {
128     try {
129       delete sac->SAVIENNACL;
130     } catch (char *ex) {
131       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "ViennaCL error: %s", ex);
132     }
133   }
134 
135   /*
136       Free the private data structure that was hanging off the PC
137   */
138   PetscCall(PetscFree(pc->data));
139   PetscFunctionReturn(0);
140 }
141 
142 static PetscErrorCode PCSetFromOptions_SAVIENNACL(PC pc, PetscOptionItems *PetscOptionsObject)
143 {
144   PetscFunctionBegin;
145   PetscOptionsHeadBegin(PetscOptionsObject, "SAVIENNACL options");
146   PetscOptionsHeadEnd();
147   PetscFunctionReturn(0);
148 }
149 
150 /*MC
151      PCSAViennaCL  - A smoothed agglomeration algorithm that can be used via the CUDA, OpenCL, and OpenMP backends of ViennaCL
152 
153    Level: advanced
154 
155    Developer Note:
156    This `PCType` does not appear to be registered
157 
158 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`
159 M*/
160 
161 PETSC_EXTERN PetscErrorCode PCCreate_SAVIENNACL(PC pc)
162 {
163   PC_SAVIENNACL *sac;
164 
165   PetscFunctionBegin;
166   /*
167      Creates the private data structure for this preconditioner and
168      attach it to the PC object.
169   */
170   PetscCall(PetscNew(&sac));
171   pc->data = (void *)sac;
172 
173   /*
174      Initialize the pointer to zero
175      Initialize number of v-cycles to default (1)
176   */
177   sac->SAVIENNACL = 0;
178 
179   /*
180       Set the pointers for the functions that are provided above.
181       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
182       are called, they will automatically call these functions.  Note we
183       choose not to provide a couple of these functions since they are
184       not needed.
185   */
186   pc->ops->apply               = PCApply_SAVIENNACL;
187   pc->ops->applytranspose      = 0;
188   pc->ops->setup               = PCSetUp_SAVIENNACL;
189   pc->ops->destroy             = PCDestroy_SAVIENNACL;
190   pc->ops->setfromoptions      = PCSetFromOptions_SAVIENNACL;
191   pc->ops->view                = 0;
192   pc->ops->applyrichardson     = 0;
193   pc->ops->applysymmetricleft  = 0;
194   pc->ops->applysymmetricright = 0;
195   PetscFunctionReturn(0);
196 }
197