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