xref: /petsc/src/ksp/pc/impls/saviennacl/saviennacl.cxx (revision 2b338477d770c45a88d8c20a3ee06ef156cc4d7c)
107598726SKarl Rupp /*
207598726SKarl Rupp    Include files needed for the ViennaCL Smoothed Aggregation preconditioner:
307598726SKarl Rupp      pcimpl.h - private include file intended for use by all preconditioners
407598726SKarl Rupp */
599acd6aaSStefano Zampini #define PETSC_SKIP_IMMINTRIN_H_CUDAWORKAROUND 1
607598726SKarl Rupp #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
707598726SKarl Rupp #include <../src/mat/impls/aij/seq/aij.h>
807598726SKarl Rupp #include <../src/vec/vec/impls/dvecimpl.h>
907598726SKarl Rupp #include <../src/mat/impls/aij/seq/seqviennacl/viennaclmatimpl.h>
1007598726SKarl Rupp #include <../src/vec/vec/impls/seq/seqviennacl/viennaclvecimpl.h>
1107598726SKarl Rupp #include <viennacl/linalg/amg.hpp>
1207598726SKarl Rupp 
1307598726SKarl Rupp /*
1407598726SKarl Rupp    Private context (data structure) for the SAVIENNACL preconditioner.
1507598726SKarl Rupp */
1607598726SKarl Rupp typedef struct {
1707598726SKarl Rupp   viennacl::linalg::amg_precond<viennacl::compressed_matrix<PetscScalar>> *SAVIENNACL;
1807598726SKarl Rupp } PC_SAVIENNACL;
1907598726SKarl Rupp 
2007598726SKarl Rupp /*
2107598726SKarl Rupp    PCSetUp_SAVIENNACL - Prepares for the use of the SAVIENNACL preconditioner
2207598726SKarl Rupp                         by setting data structures and options.
2307598726SKarl Rupp 
2407598726SKarl Rupp    Input Parameter:
2507598726SKarl Rupp .  pc - the preconditioner context
2607598726SKarl Rupp 
2707598726SKarl Rupp    Application Interface Routine: PCSetUp()
2807598726SKarl Rupp 
29f1580f4eSBarry Smith    Note:
3007598726SKarl Rupp    The interface routine PCSetUp() is not usually called directly by
3107598726SKarl Rupp    the user, but instead is called by PCApply() if necessary.
3207598726SKarl Rupp */
33d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_SAVIENNACL(PC pc)
34d71ae5a4SJacob Faibussowitsch {
3507598726SKarl Rupp   PC_SAVIENNACL      *sa  = (PC_SAVIENNACL *)pc->data;
3607598726SKarl Rupp   PetscBool           flg = PETSC_FALSE;
3707598726SKarl Rupp   Mat_SeqAIJViennaCL *gpustruct;
3807598726SKarl Rupp 
3907598726SKarl Rupp   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATSEQAIJVIENNACL, &flg));
4128b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Currently only handles ViennaCL matrices");
42*371d2eb7SMartin Diehl   if (pc->setupcalled) {
4307598726SKarl Rupp     try {
4407598726SKarl Rupp       delete sa->SAVIENNACL;
45d71ae5a4SJacob Faibussowitsch     } catch (char *ex) {
46d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "ViennaCL error: %s", ex);
47d71ae5a4SJacob Faibussowitsch     }
4807598726SKarl Rupp   }
4907598726SKarl Rupp   try {
5007598726SKarl Rupp #if defined(PETSC_USE_COMPLEX)
5107598726SKarl Rupp     gpustruct = NULL;
5207598726SKarl Rupp     SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "No support for complex arithmetic in SAVIENNACL preconditioner");
5307598726SKarl Rupp #else
549566063dSJacob Faibussowitsch     PetscCall(MatViennaCLCopyToGPU(pc->pmat));
55f4f49eeaSPierre Jolivet     gpustruct = (Mat_SeqAIJViennaCL *)pc->pmat->spptr;
5607598726SKarl Rupp 
5707598726SKarl Rupp     viennacl::linalg::amg_tag amg_tag_sa_pmis;
5807598726SKarl Rupp     amg_tag_sa_pmis.set_coarsening_method(viennacl::linalg::AMG_COARSENING_METHOD_MIS2_AGGREGATION);
5907598726SKarl Rupp     amg_tag_sa_pmis.set_interpolation_method(viennacl::linalg::AMG_INTERPOLATION_METHOD_SMOOTHED_AGGREGATION);
6007598726SKarl Rupp     ViennaCLAIJMatrix *mat = (ViennaCLAIJMatrix *)gpustruct->mat;
6107598726SKarl Rupp     sa->SAVIENNACL         = new viennacl::linalg::amg_precond<viennacl::compressed_matrix<PetscScalar>>(*mat, amg_tag_sa_pmis);
6207598726SKarl Rupp     sa->SAVIENNACL->setup();
6307598726SKarl Rupp #endif
64d71ae5a4SJacob Faibussowitsch   } catch (char *ex) {
65d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "ViennaCL error: %s", ex);
66d71ae5a4SJacob Faibussowitsch   }
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6807598726SKarl Rupp }
6907598726SKarl Rupp 
7007598726SKarl Rupp /*
7107598726SKarl Rupp    PCApply_SAVIENNACL - Applies the SAVIENNACL preconditioner to a vector.
7207598726SKarl Rupp 
7307598726SKarl Rupp    Input Parameters:
7407598726SKarl Rupp .  pc - the preconditioner context
7507598726SKarl Rupp .  x - input vector
7607598726SKarl Rupp 
7707598726SKarl Rupp    Output Parameter:
7807598726SKarl Rupp .  y - output vector
7907598726SKarl Rupp 
8007598726SKarl Rupp    Application Interface Routine: PCApply()
8107598726SKarl Rupp  */
82d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_SAVIENNACL(PC pc, Vec x, Vec y)
83d71ae5a4SJacob Faibussowitsch {
8407598726SKarl Rupp   PC_SAVIENNACL                       *sac = (PC_SAVIENNACL *)pc->data;
8507598726SKarl Rupp   PetscBool                            flg1, flg2;
8607598726SKarl Rupp   viennacl::vector<PetscScalar> const *xarray = NULL;
8707598726SKarl Rupp   viennacl::vector<PetscScalar>       *yarray = NULL;
8807598726SKarl Rupp 
8907598726SKarl Rupp   PetscFunctionBegin;
9007598726SKarl Rupp   /*how to apply a certain fixed number of iterations?*/
919566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)x, VECSEQVIENNACL, &flg1));
929566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)y, VECSEQVIENNACL, &flg2));
9357508eceSPierre Jolivet   PetscCheck(flg1 && flg2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Currently only handles ViennaCL vectors");
9448a46eb9SPierre Jolivet   if (!sac->SAVIENNACL) PetscCall(PCSetUp_SAVIENNACL(pc));
959566063dSJacob Faibussowitsch   PetscCall(VecViennaCLGetArrayRead(x, &xarray));
969566063dSJacob Faibussowitsch   PetscCall(VecViennaCLGetArrayWrite(y, &yarray));
9707598726SKarl Rupp   try {
98c4163675SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
9907598726SKarl Rupp     *yarray = *xarray;
10007598726SKarl Rupp     sac->SAVIENNACL->apply(*yarray);
10107598726SKarl Rupp #endif
102d71ae5a4SJacob Faibussowitsch   } catch (char *ex) {
103d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "ViennaCL error: %s", ex);
104d71ae5a4SJacob Faibussowitsch   }
1059566063dSJacob Faibussowitsch   PetscCall(VecViennaCLRestoreArrayRead(x, &xarray));
1069566063dSJacob Faibussowitsch   PetscCall(VecViennaCLRestoreArrayWrite(y, &yarray));
1079566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)y));
1083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10907598726SKarl Rupp }
110f1580f4eSBarry Smith 
11107598726SKarl Rupp /*
11207598726SKarl Rupp    PCDestroy_SAVIENNACL - Destroys the private context for the SAVIENNACL preconditioner
11307598726SKarl Rupp    that was created with PCCreate_SAVIENNACL().
11407598726SKarl Rupp 
11507598726SKarl Rupp    Input Parameter:
11607598726SKarl Rupp .  pc - the preconditioner context
11707598726SKarl Rupp 
11807598726SKarl Rupp    Application Interface Routine: PCDestroy()
11907598726SKarl Rupp */
120d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_SAVIENNACL(PC pc)
121d71ae5a4SJacob Faibussowitsch {
12207598726SKarl Rupp   PC_SAVIENNACL *sac = (PC_SAVIENNACL *)pc->data;
12307598726SKarl Rupp 
12407598726SKarl Rupp   PetscFunctionBegin;
12507598726SKarl Rupp   if (sac->SAVIENNACL) {
12607598726SKarl Rupp     try {
12707598726SKarl Rupp       delete sac->SAVIENNACL;
128d71ae5a4SJacob Faibussowitsch     } catch (char *ex) {
129d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "ViennaCL error: %s", ex);
130d71ae5a4SJacob Faibussowitsch     }
13107598726SKarl Rupp   }
13207598726SKarl Rupp 
13307598726SKarl Rupp   /*
13407598726SKarl Rupp       Free the private data structure that was hanging off the PC
13507598726SKarl Rupp   */
1369566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13807598726SKarl Rupp }
13907598726SKarl Rupp 
140ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_SAVIENNACL(PC pc, PetscOptionItems PetscOptionsObject)
141d71ae5a4SJacob Faibussowitsch {
14207598726SKarl Rupp   PetscFunctionBegin;
143d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SAVIENNACL options");
144d0609cedSBarry Smith   PetscOptionsHeadEnd();
1453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14607598726SKarl Rupp }
14707598726SKarl Rupp 
14807598726SKarl Rupp /*MC
14904c3f3b8SBarry Smith   PCSAVIENNACL  - A smoothed agglomeration algorithm that can be used via the CUDA, OpenCL, and OpenMP backends of ViennaCL
15007598726SKarl Rupp 
15107598726SKarl Rupp   Level: advanced
15207598726SKarl Rupp 
153feefa0e1SJacob Faibussowitsch   Developer Notes:
154f1580f4eSBarry Smith   This `PCType` does not appear to be registered
15507598726SKarl Rupp 
156562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`
15707598726SKarl Rupp M*/
15804c3f3b8SBarry Smith 
159d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_SAVIENNACL(PC pc)
160d71ae5a4SJacob Faibussowitsch {
16107598726SKarl Rupp   PC_SAVIENNACL *sac;
16207598726SKarl Rupp 
16307598726SKarl Rupp   PetscFunctionBegin;
16407598726SKarl Rupp   /*
16507598726SKarl Rupp      Creates the private data structure for this preconditioner and
16607598726SKarl Rupp      attach it to the PC object.
16707598726SKarl Rupp   */
1684dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&sac));
16907598726SKarl Rupp   pc->data = (void *)sac;
17007598726SKarl Rupp 
17107598726SKarl Rupp   /*
17207598726SKarl Rupp      Initialize the pointer to zero
17307598726SKarl Rupp      Initialize number of v-cycles to default (1)
17407598726SKarl Rupp   */
17507598726SKarl Rupp   sac->SAVIENNACL = 0;
17607598726SKarl Rupp 
17707598726SKarl Rupp   /*
17807598726SKarl Rupp       Set the pointers for the functions that are provided above.
17907598726SKarl Rupp       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
18007598726SKarl Rupp       are called, they will automatically call these functions.  Note we
18107598726SKarl Rupp       choose not to provide a couple of these functions since they are
18207598726SKarl Rupp       not needed.
18307598726SKarl Rupp   */
18407598726SKarl Rupp   pc->ops->apply               = PCApply_SAVIENNACL;
18507598726SKarl Rupp   pc->ops->applytranspose      = 0;
18607598726SKarl Rupp   pc->ops->setup               = PCSetUp_SAVIENNACL;
18707598726SKarl Rupp   pc->ops->destroy             = PCDestroy_SAVIENNACL;
18807598726SKarl Rupp   pc->ops->setfromoptions      = PCSetFromOptions_SAVIENNACL;
18907598726SKarl Rupp   pc->ops->view                = 0;
19007598726SKarl Rupp   pc->ops->applyrichardson     = 0;
19107598726SKarl Rupp   pc->ops->applysymmetricleft  = 0;
19207598726SKarl Rupp   pc->ops->applysymmetricright = 0;
1933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19407598726SKarl Rupp }
195