170baa948SKarl Rupp /* 270baa948SKarl Rupp Include files needed for the ViennaCL row-scaling preconditioner: 370baa948SKarl Rupp pcimpl.h - private include file intended for use by all preconditioners 470baa948SKarl Rupp */ 599acd6aaSStefano Zampini #define PETSC_SKIP_IMMINTRIN_H_CUDAWORKAROUND 1 699acd6aaSStefano Zampini 770baa948SKarl Rupp #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 870baa948SKarl Rupp #include <../src/mat/impls/aij/seq/aij.h> 970baa948SKarl Rupp #include <../src/vec/vec/impls/dvecimpl.h> 1070baa948SKarl Rupp #include <../src/mat/impls/aij/seq/seqviennacl/viennaclmatimpl.h> 1170baa948SKarl Rupp #include <../src/vec/vec/impls/seq/seqviennacl/viennaclvecimpl.h> 1270baa948SKarl Rupp #include <viennacl/linalg/row_scaling.hpp> 1370baa948SKarl Rupp 1470baa948SKarl Rupp /* 1570baa948SKarl Rupp Private context (data structure) for the ROWSCALINGVIENNACL preconditioner. 1670baa948SKarl Rupp */ 1770baa948SKarl Rupp typedef struct { 1870baa948SKarl Rupp viennacl::linalg::row_scaling<viennacl::compressed_matrix<PetscScalar>> *ROWSCALINGVIENNACL; 1970baa948SKarl Rupp } PC_ROWSCALINGVIENNACL; 2070baa948SKarl Rupp 2170baa948SKarl Rupp /* 2270baa948SKarl Rupp PCSetUp_ROWSCALINGVIENNACL - Prepares for the use of the ROWSCALINGVIENNACL preconditioner 2370baa948SKarl Rupp by setting data structures and options. 2470baa948SKarl Rupp 2570baa948SKarl Rupp Input Parameter: 2670baa948SKarl Rupp . pc - the preconditioner context 2770baa948SKarl Rupp 2870baa948SKarl Rupp Application Interface Routine: PCSetUp() 2970baa948SKarl Rupp 30f1580f4eSBarry Smith Note: 3170baa948SKarl Rupp The interface routine PCSetUp() is not usually called directly by 3270baa948SKarl Rupp the user, but instead is called by PCApply() if necessary. 3370baa948SKarl Rupp */ 34d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_ROWSCALINGVIENNACL(PC pc) 35d71ae5a4SJacob Faibussowitsch { 3670baa948SKarl Rupp PC_ROWSCALINGVIENNACL *rowscaling = (PC_ROWSCALINGVIENNACL *)pc->data; 3770baa948SKarl Rupp PetscBool flg = PETSC_FALSE; 3870baa948SKarl Rupp Mat_SeqAIJViennaCL *gpustruct; 3970baa948SKarl Rupp 4070baa948SKarl Rupp PetscFunctionBegin; 419566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATSEQAIJVIENNACL, &flg)); 4228b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Currently only handles ViennaCL matrices"); 43*371d2eb7SMartin Diehl if (pc->setupcalled) { 4470baa948SKarl Rupp try { 4570baa948SKarl Rupp delete rowscaling->ROWSCALINGVIENNACL; 46d71ae5a4SJacob Faibussowitsch } catch (char *ex) { 47d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "ViennaCL error: %s", ex); 48d71ae5a4SJacob Faibussowitsch } 4970baa948SKarl Rupp } 5070baa948SKarl Rupp try { 5170baa948SKarl Rupp #if defined(PETSC_USE_COMPLEX) 5270baa948SKarl Rupp gpustruct = NULL; 5370baa948SKarl Rupp SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "No support for complex arithmetic in ROWSCALINGVIENNACL preconditioner"); 5470baa948SKarl Rupp #else 559566063dSJacob Faibussowitsch PetscCall(MatViennaCLCopyToGPU(pc->pmat)); 56f4f49eeaSPierre Jolivet gpustruct = (Mat_SeqAIJViennaCL *)pc->pmat->spptr; 5770baa948SKarl Rupp 5870baa948SKarl Rupp viennacl::linalg::row_scaling_tag pc_tag(1); 5970baa948SKarl Rupp ViennaCLAIJMatrix *mat = (ViennaCLAIJMatrix *)gpustruct->mat; 6070baa948SKarl Rupp rowscaling->ROWSCALINGVIENNACL = new viennacl::linalg::row_scaling<viennacl::compressed_matrix<PetscScalar>>(*mat, pc_tag); 6170baa948SKarl Rupp #endif 62d71ae5a4SJacob Faibussowitsch } catch (char *ex) { 63d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "ViennaCL error: %s", ex); 64d71ae5a4SJacob Faibussowitsch } 653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6670baa948SKarl Rupp } 6770baa948SKarl Rupp 6870baa948SKarl Rupp /* 6970baa948SKarl Rupp PCApply_ROWSCALINGVIENNACL - Applies the ROWSCALINGVIENNACL preconditioner to a vector. 7070baa948SKarl Rupp 7170baa948SKarl Rupp Input Parameters: 7270baa948SKarl Rupp . pc - the preconditioner context 7370baa948SKarl Rupp . x - input vector 7470baa948SKarl Rupp 7570baa948SKarl Rupp Output Parameter: 7670baa948SKarl Rupp . y - output vector 7770baa948SKarl Rupp 7870baa948SKarl Rupp Application Interface Routine: PCApply() 7970baa948SKarl Rupp */ 80d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_ROWSCALINGVIENNACL(PC pc, Vec x, Vec y) 81d71ae5a4SJacob Faibussowitsch { 8270baa948SKarl Rupp PC_ROWSCALINGVIENNACL *ilu = (PC_ROWSCALINGVIENNACL *)pc->data; 8370baa948SKarl Rupp PetscBool flg1, flg2; 8470baa948SKarl Rupp viennacl::vector<PetscScalar> const *xarray = NULL; 8570baa948SKarl Rupp viennacl::vector<PetscScalar> *yarray = NULL; 8670baa948SKarl Rupp 8770baa948SKarl Rupp PetscFunctionBegin; 8870baa948SKarl Rupp /*how to apply a certain fixed number of iterations?*/ 899566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)x, VECSEQVIENNACL, &flg1)); 909566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)y, VECSEQVIENNACL, &flg2)); 9157508eceSPierre Jolivet PetscCheck(flg1 && flg2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Currently only handles ViennaCL vectors"); 9248a46eb9SPierre Jolivet if (!ilu->ROWSCALINGVIENNACL) PetscCall(PCSetUp_ROWSCALINGVIENNACL(pc)); 939566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 949566063dSJacob Faibussowitsch PetscCall(VecViennaCLGetArrayRead(x, &xarray)); 959566063dSJacob Faibussowitsch PetscCall(VecViennaCLGetArrayWrite(y, &yarray)); 9670baa948SKarl Rupp try { 9770baa948SKarl Rupp #if defined(PETSC_USE_COMPLEX) 9870baa948SKarl Rupp 9970baa948SKarl Rupp #else 10070baa948SKarl Rupp *yarray = *xarray; 10170baa948SKarl Rupp ilu->ROWSCALINGVIENNACL->apply(*yarray); 10270baa948SKarl Rupp #endif 103d71ae5a4SJacob Faibussowitsch } catch (char *ex) { 104d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "ViennaCL error: %s", ex); 105d71ae5a4SJacob Faibussowitsch } 1069566063dSJacob Faibussowitsch PetscCall(VecViennaCLRestoreArrayRead(x, &xarray)); 1079566063dSJacob Faibussowitsch PetscCall(VecViennaCLRestoreArrayWrite(y, &yarray)); 1089566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)y)); 1093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11070baa948SKarl Rupp } 111f1580f4eSBarry Smith 11270baa948SKarl Rupp /* 11370baa948SKarl Rupp PCDestroy_ROWSCALINGVIENNACL - Destroys the private context for the ROWSCALINGVIENNACL preconditioner 11470baa948SKarl Rupp that was created with PCCreate_ROWSCALINGVIENNACL(). 11570baa948SKarl Rupp 11670baa948SKarl Rupp Input Parameter: 11770baa948SKarl Rupp . pc - the preconditioner context 11870baa948SKarl Rupp 11970baa948SKarl Rupp Application Interface Routine: PCDestroy() 12070baa948SKarl Rupp */ 121d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_ROWSCALINGVIENNACL(PC pc) 122d71ae5a4SJacob Faibussowitsch { 12370baa948SKarl Rupp PC_ROWSCALINGVIENNACL *rowscaling = (PC_ROWSCALINGVIENNACL *)pc->data; 12470baa948SKarl Rupp 12570baa948SKarl Rupp PetscFunctionBegin; 12670baa948SKarl Rupp if (rowscaling->ROWSCALINGVIENNACL) { 12770baa948SKarl Rupp try { 12870baa948SKarl Rupp delete rowscaling->ROWSCALINGVIENNACL; 129d71ae5a4SJacob Faibussowitsch } catch (char *ex) { 130d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "ViennaCL error: %s", ex); 131d71ae5a4SJacob Faibussowitsch } 13270baa948SKarl Rupp } 13370baa948SKarl Rupp 13470baa948SKarl Rupp /* 13570baa948SKarl Rupp Free the private data structure that was hanging off the PC 13670baa948SKarl Rupp */ 1379566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13970baa948SKarl Rupp } 14070baa948SKarl Rupp 141ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_ROWSCALINGVIENNACL(PC pc, PetscOptionItems PetscOptionsObject) 142d71ae5a4SJacob Faibussowitsch { 14370baa948SKarl Rupp PetscFunctionBegin; 144d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "ROWSCALINGVIENNACL options"); 145d0609cedSBarry Smith PetscOptionsHeadEnd(); 1463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14770baa948SKarl Rupp } 14870baa948SKarl Rupp 14970baa948SKarl Rupp /*MC 15004c3f3b8SBarry Smith PCROWSCALINGVIENNACLL - A diagonal preconditioner (scaling rows of matrices by their norm) that can be used via the CUDA, 15104c3f3b8SBarry Smith OpenCL, and OpenMP backends of ViennaCL 15270baa948SKarl Rupp 15370baa948SKarl Rupp Level: advanced 15470baa948SKarl Rupp 155562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC` 15670baa948SKarl Rupp M*/ 157d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ROWSCALINGVIENNACL(PC pc) 158d71ae5a4SJacob Faibussowitsch { 15970baa948SKarl Rupp PC_ROWSCALINGVIENNACL *rowscaling; 16070baa948SKarl Rupp 16170baa948SKarl Rupp PetscFunctionBegin; 16270baa948SKarl Rupp /* 16370baa948SKarl Rupp Creates the private data structure for this preconditioner and 16470baa948SKarl Rupp attach it to the PC object. 16570baa948SKarl Rupp */ 1664dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&rowscaling)); 16770baa948SKarl Rupp pc->data = (void *)rowscaling; 16870baa948SKarl Rupp 16970baa948SKarl Rupp /* 17070baa948SKarl Rupp Initialize the pointer to zero 17170baa948SKarl Rupp Initialize number of v-cycles to default (1) 17270baa948SKarl Rupp */ 17370baa948SKarl Rupp rowscaling->ROWSCALINGVIENNACL = 0; 17470baa948SKarl Rupp 17570baa948SKarl Rupp /* 17670baa948SKarl Rupp Set the pointers for the functions that are provided above. 17770baa948SKarl Rupp Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 17870baa948SKarl Rupp are called, they will automatically call these functions. Note we 17970baa948SKarl Rupp choose not to provide a couple of these functions since they are 18070baa948SKarl Rupp not needed. 18170baa948SKarl Rupp */ 18270baa948SKarl Rupp pc->ops->apply = PCApply_ROWSCALINGVIENNACL; 18370baa948SKarl Rupp pc->ops->applytranspose = 0; 18470baa948SKarl Rupp pc->ops->setup = PCSetUp_ROWSCALINGVIENNACL; 18570baa948SKarl Rupp pc->ops->destroy = PCDestroy_ROWSCALINGVIENNACL; 18670baa948SKarl Rupp pc->ops->setfromoptions = PCSetFromOptions_ROWSCALINGVIENNACL; 18770baa948SKarl Rupp pc->ops->view = 0; 18870baa948SKarl Rupp pc->ops->applyrichardson = 0; 18970baa948SKarl Rupp pc->ops->applysymmetricleft = 0; 19070baa948SKarl Rupp pc->ops->applysymmetricright = 0; 1913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19270baa948SKarl Rupp } 193