xref: /petsc/src/ksp/pc/impls/vpbjacobi/vpbjacobi.c (revision 4877389942e948c907dd43b3a2d019cdca14ea44)
1f1be3500SJunchao Zhang #include <../src/ksp/pc/impls/vpbjacobi/vpbjacobi.h>
20da83c2eSBarry Smith 
3d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_VPBJacobi(PC pc, Vec x, Vec y)
4d71ae5a4SJacob Faibussowitsch {
50da83c2eSBarry Smith   PC_VPBJacobi      *jac = (PC_VPBJacobi *)pc->data;
60da83c2eSBarry Smith   PetscInt           i, ncnt = 0;
70da83c2eSBarry Smith   const MatScalar   *diag = jac->diag;
80da83c2eSBarry Smith   PetscInt           ib, jb, bs;
90da83c2eSBarry Smith   const PetscScalar *xx;
100da83c2eSBarry Smith   PetscScalar       *yy, x0, x1, x2, x3, x4, x5, x6;
110da83c2eSBarry Smith   PetscInt           nblocks;
120da83c2eSBarry Smith   const PetscInt    *bsizes;
130da83c2eSBarry Smith 
140da83c2eSBarry Smith   PetscFunctionBegin;
159566063dSJacob Faibussowitsch   PetscCall(MatGetVariableBlockSizes(pc->pmat, &nblocks, &bsizes));
169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
179566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
180da83c2eSBarry Smith   for (i = 0; i < nblocks; i++) {
190da83c2eSBarry Smith     bs = bsizes[i];
200da83c2eSBarry Smith     switch (bs) {
21d71ae5a4SJacob Faibussowitsch     case 1:
22d71ae5a4SJacob Faibussowitsch       yy[ncnt] = *diag * xx[ncnt];
23d71ae5a4SJacob Faibussowitsch       break;
240da83c2eSBarry Smith     case 2:
259371c9d4SSatish Balay       x0           = xx[ncnt];
269371c9d4SSatish Balay       x1           = xx[ncnt + 1];
270da83c2eSBarry Smith       yy[ncnt]     = diag[0] * x0 + diag[2] * x1;
280da83c2eSBarry Smith       yy[ncnt + 1] = diag[1] * x0 + diag[3] * x1;
290da83c2eSBarry Smith       break;
300da83c2eSBarry Smith     case 3:
319371c9d4SSatish Balay       x0           = xx[ncnt];
329371c9d4SSatish Balay       x1           = xx[ncnt + 1];
339371c9d4SSatish Balay       x2           = xx[ncnt + 2];
340da83c2eSBarry Smith       yy[ncnt]     = diag[0] * x0 + diag[3] * x1 + diag[6] * x2;
350da83c2eSBarry Smith       yy[ncnt + 1] = diag[1] * x0 + diag[4] * x1 + diag[7] * x2;
360da83c2eSBarry Smith       yy[ncnt + 2] = diag[2] * x0 + diag[5] * x1 + diag[8] * x2;
370da83c2eSBarry Smith       break;
380da83c2eSBarry Smith     case 4:
399371c9d4SSatish Balay       x0           = xx[ncnt];
409371c9d4SSatish Balay       x1           = xx[ncnt + 1];
419371c9d4SSatish Balay       x2           = xx[ncnt + 2];
429371c9d4SSatish Balay       x3           = xx[ncnt + 3];
430da83c2eSBarry Smith       yy[ncnt]     = diag[0] * x0 + diag[4] * x1 + diag[8] * x2 + diag[12] * x3;
440da83c2eSBarry Smith       yy[ncnt + 1] = diag[1] * x0 + diag[5] * x1 + diag[9] * x2 + diag[13] * x3;
450da83c2eSBarry Smith       yy[ncnt + 2] = diag[2] * x0 + diag[6] * x1 + diag[10] * x2 + diag[14] * x3;
460da83c2eSBarry Smith       yy[ncnt + 3] = diag[3] * x0 + diag[7] * x1 + diag[11] * x2 + diag[15] * x3;
470da83c2eSBarry Smith       break;
480da83c2eSBarry Smith     case 5:
499371c9d4SSatish Balay       x0           = xx[ncnt];
509371c9d4SSatish Balay       x1           = xx[ncnt + 1];
519371c9d4SSatish Balay       x2           = xx[ncnt + 2];
529371c9d4SSatish Balay       x3           = xx[ncnt + 3];
539371c9d4SSatish Balay       x4           = xx[ncnt + 4];
540da83c2eSBarry Smith       yy[ncnt]     = diag[0] * x0 + diag[5] * x1 + diag[10] * x2 + diag[15] * x3 + diag[20] * x4;
550da83c2eSBarry Smith       yy[ncnt + 1] = diag[1] * x0 + diag[6] * x1 + diag[11] * x2 + diag[16] * x3 + diag[21] * x4;
560da83c2eSBarry Smith       yy[ncnt + 2] = diag[2] * x0 + diag[7] * x1 + diag[12] * x2 + diag[17] * x3 + diag[22] * x4;
570da83c2eSBarry Smith       yy[ncnt + 3] = diag[3] * x0 + diag[8] * x1 + diag[13] * x2 + diag[18] * x3 + diag[23] * x4;
580da83c2eSBarry Smith       yy[ncnt + 4] = diag[4] * x0 + diag[9] * x1 + diag[14] * x2 + diag[19] * x3 + diag[24] * x4;
590da83c2eSBarry Smith       break;
600da83c2eSBarry Smith     case 6:
619371c9d4SSatish Balay       x0           = xx[ncnt];
629371c9d4SSatish Balay       x1           = xx[ncnt + 1];
639371c9d4SSatish Balay       x2           = xx[ncnt + 2];
649371c9d4SSatish Balay       x3           = xx[ncnt + 3];
659371c9d4SSatish Balay       x4           = xx[ncnt + 4];
669371c9d4SSatish Balay       x5           = xx[ncnt + 5];
670da83c2eSBarry Smith       yy[ncnt]     = diag[0] * x0 + diag[6] * x1 + diag[12] * x2 + diag[18] * x3 + diag[24] * x4 + diag[30] * x5;
680da83c2eSBarry Smith       yy[ncnt + 1] = diag[1] * x0 + diag[7] * x1 + diag[13] * x2 + diag[19] * x3 + diag[25] * x4 + diag[31] * x5;
690da83c2eSBarry Smith       yy[ncnt + 2] = diag[2] * x0 + diag[8] * x1 + diag[14] * x2 + diag[20] * x3 + diag[26] * x4 + diag[32] * x5;
700da83c2eSBarry Smith       yy[ncnt + 3] = diag[3] * x0 + diag[9] * x1 + diag[15] * x2 + diag[21] * x3 + diag[27] * x4 + diag[33] * x5;
710da83c2eSBarry Smith       yy[ncnt + 4] = diag[4] * x0 + diag[10] * x1 + diag[16] * x2 + diag[22] * x3 + diag[28] * x4 + diag[34] * x5;
720da83c2eSBarry Smith       yy[ncnt + 5] = diag[5] * x0 + diag[11] * x1 + diag[17] * x2 + diag[23] * x3 + diag[29] * x4 + diag[35] * x5;
730da83c2eSBarry Smith       break;
740da83c2eSBarry Smith     case 7:
759371c9d4SSatish Balay       x0           = xx[ncnt];
769371c9d4SSatish Balay       x1           = xx[ncnt + 1];
779371c9d4SSatish Balay       x2           = xx[ncnt + 2];
789371c9d4SSatish Balay       x3           = xx[ncnt + 3];
799371c9d4SSatish Balay       x4           = xx[ncnt + 4];
809371c9d4SSatish Balay       x5           = xx[ncnt + 5];
819371c9d4SSatish Balay       x6           = xx[ncnt + 6];
820da83c2eSBarry Smith       yy[ncnt]     = diag[0] * x0 + diag[7] * x1 + diag[14] * x2 + diag[21] * x3 + diag[28] * x4 + diag[35] * x5 + diag[42] * x6;
830da83c2eSBarry Smith       yy[ncnt + 1] = diag[1] * x0 + diag[8] * x1 + diag[15] * x2 + diag[22] * x3 + diag[29] * x4 + diag[36] * x5 + diag[43] * x6;
840da83c2eSBarry Smith       yy[ncnt + 2] = diag[2] * x0 + diag[9] * x1 + diag[16] * x2 + diag[23] * x3 + diag[30] * x4 + diag[37] * x5 + diag[44] * x6;
850da83c2eSBarry Smith       yy[ncnt + 3] = diag[3] * x0 + diag[10] * x1 + diag[17] * x2 + diag[24] * x3 + diag[31] * x4 + diag[38] * x5 + diag[45] * x6;
860da83c2eSBarry Smith       yy[ncnt + 4] = diag[4] * x0 + diag[11] * x1 + diag[18] * x2 + diag[25] * x3 + diag[32] * x4 + diag[39] * x5 + diag[46] * x6;
870da83c2eSBarry Smith       yy[ncnt + 5] = diag[5] * x0 + diag[12] * x1 + diag[19] * x2 + diag[26] * x3 + diag[33] * x4 + diag[40] * x5 + diag[47] * x6;
880da83c2eSBarry Smith       yy[ncnt + 6] = diag[6] * x0 + diag[13] * x1 + diag[20] * x2 + diag[27] * x3 + diag[34] * x4 + diag[41] * x5 + diag[48] * x6;
890da83c2eSBarry Smith       break;
900da83c2eSBarry Smith     default:
910da83c2eSBarry Smith       for (ib = 0; ib < bs; ib++) {
920da83c2eSBarry Smith         PetscScalar rowsum = 0;
93ad540459SPierre Jolivet         for (jb = 0; jb < bs; jb++) rowsum += diag[ib + jb * bs] * xx[ncnt + jb];
940da83c2eSBarry Smith         yy[ncnt + ib] = rowsum;
950da83c2eSBarry Smith       }
960da83c2eSBarry Smith     }
970da83c2eSBarry Smith     ncnt += bsizes[i];
980da83c2eSBarry Smith     diag += bsizes[i] * bsizes[i];
990da83c2eSBarry Smith   }
1009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
1019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
1020da83c2eSBarry Smith   PetscFunctionReturn(0);
1030da83c2eSBarry Smith }
1040da83c2eSBarry Smith 
105d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCSetUp_VPBJacobi_Host(PC pc)
106d71ae5a4SJacob Faibussowitsch {
1070da83c2eSBarry Smith   PC_VPBJacobi   *jac = (PC_VPBJacobi *)pc->data;
1080da83c2eSBarry Smith   Mat             A   = pc->pmat;
1090da83c2eSBarry Smith   MatFactorError  err;
1100da83c2eSBarry Smith   PetscInt        i, nsize = 0, nlocal;
1110da83c2eSBarry Smith   PetscInt        nblocks;
1120da83c2eSBarry Smith   const PetscInt *bsizes;
1130da83c2eSBarry Smith 
1140da83c2eSBarry Smith   PetscFunctionBegin;
1159566063dSJacob Faibussowitsch   PetscCall(MatGetVariableBlockSizes(pc->pmat, &nblocks, &bsizes));
1169566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &nlocal, NULL));
11708401ef6SPierre Jolivet   PetscCheck(!nlocal || nblocks, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatSetVariableBlockSizes() before using PCVPBJACOBI");
1180da83c2eSBarry Smith   if (!jac->diag) {
1190da83c2eSBarry Smith     for (i = 0; i < nblocks; i++) nsize += bsizes[i] * bsizes[i];
1209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nsize, &jac->diag));
1210da83c2eSBarry Smith   }
1229566063dSJacob Faibussowitsch   PetscCall(MatInvertVariableBlockDiagonal(A, nblocks, bsizes, jac->diag));
1239566063dSJacob Faibussowitsch   PetscCall(MatFactorGetError(A, &err));
1240da83c2eSBarry Smith   if (err) pc->failedreason = (PCFailedReason)err;
1250da83c2eSBarry Smith   pc->ops->apply = PCApply_VPBJacobi;
1260da83c2eSBarry Smith   PetscFunctionReturn(0);
1270da83c2eSBarry Smith }
1288a9c020eSBarry Smith 
129d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_VPBJacobi(PC pc)
130d71ae5a4SJacob Faibussowitsch {
131f1be3500SJunchao Zhang   PetscFunctionBegin;
132f1be3500SJunchao Zhang   /* In PCCreate_VPBJacobi() pmat might have not been set, so we wait to the last minute to do the dispatch */
133f1be3500SJunchao Zhang #if defined(PETSC_HAVE_CUDA)
134f1be3500SJunchao Zhang   PetscBool isCuda;
135f1be3500SJunchao Zhang   PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isCuda, MATSEQAIJCUSPARSE, MATMPIAIJCUSPARSE, ""));
136f1be3500SJunchao Zhang #endif
137f1be3500SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS_KERNELS)
138f1be3500SJunchao Zhang   PetscBool isKok;
139f1be3500SJunchao Zhang   PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isKok, MATSEQAIJKOKKOS, MATMPIAIJKOKKOS, ""));
140f1be3500SJunchao Zhang #endif
141f1be3500SJunchao Zhang 
142f1be3500SJunchao Zhang #if defined(PETSC_HAVE_CUDA)
1439371c9d4SSatish Balay   if (isCuda) PetscCall(PCSetUp_VPBJacobi_CUDA(pc));
1449371c9d4SSatish Balay   else
145f1be3500SJunchao Zhang #endif
146f1be3500SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS_KERNELS)
1479371c9d4SSatish Balay     if (isKok)
1489371c9d4SSatish Balay     PetscCall(PCSetUp_VPBJacobi_Kokkos(pc));
1499371c9d4SSatish Balay   else
150f1be3500SJunchao Zhang #endif
1519371c9d4SSatish Balay   {
1529371c9d4SSatish Balay     PetscCall(PCSetUp_VPBJacobi_Host(pc));
1539371c9d4SSatish Balay   }
154f1be3500SJunchao Zhang   PetscFunctionReturn(0);
155f1be3500SJunchao Zhang }
156f1be3500SJunchao Zhang 
157d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCDestroy_VPBJacobi(PC pc)
158d71ae5a4SJacob Faibussowitsch {
1590da83c2eSBarry Smith   PC_VPBJacobi *jac = (PC_VPBJacobi *)pc->data;
1600da83c2eSBarry Smith 
1610da83c2eSBarry Smith   PetscFunctionBegin;
1620da83c2eSBarry Smith   /*
1630da83c2eSBarry Smith       Free the private data structure that was hanging off the PC
1640da83c2eSBarry Smith   */
1659566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->diag));
1669566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1670da83c2eSBarry Smith   PetscFunctionReturn(0);
1680da83c2eSBarry Smith }
1690da83c2eSBarry Smith 
1700da83c2eSBarry Smith /*MC
1710da83c2eSBarry Smith      PCVPBJACOBI - Variable size point block Jacobi preconditioner
1720da83c2eSBarry Smith 
173f1580f4eSBarry Smith    Level: beginner
1740da83c2eSBarry Smith 
175f1580f4eSBarry Smith    Notes:
176f1580f4eSBarry Smith      See `PCJACOBI` for point Jacobi preconditioning, `PCPBJACOBI` for fixed point block size, and `PCBJACOBI` for large size blocks
177f1580f4eSBarry Smith 
178f1580f4eSBarry Smith      This works for `MATAIJ` matrices
1790da83c2eSBarry Smith 
1800da83c2eSBarry Smith      Uses dense LU factorization with partial pivoting to invert the blocks; if a zero pivot
1810da83c2eSBarry Smith      is detected a PETSc error is generated.
1820da83c2eSBarry Smith 
183f1580f4eSBarry Smith      One must call `MatSetVariableBlockSizes()` to use this preconditioner
18490dfcc32SBarry Smith 
1850da83c2eSBarry Smith    Developer Notes:
186f1580f4eSBarry Smith      This should support the `PCSetErrorIfFailure()` flag set to `PETSC_TRUE` to allow
1870da83c2eSBarry Smith      the factorization to continue even after a zero pivot is found resulting in a Nan and hence
188*48773899SPierre Jolivet      terminating `KSP` with a `KSP_DIVERGED_NANORINF` allowing
1890da83c2eSBarry Smith      a nonlinear solver/ODE integrator to recover without stopping the program as currently happens.
1900da83c2eSBarry Smith 
1910da83c2eSBarry Smith      Perhaps should provide an option that allows generation of a valid preconditioner
192f1580f4eSBarry Smith      even if a block is singular as the `PCJACOBI` does.
1930da83c2eSBarry Smith 
194db781477SPatrick Sanan .seealso: `MatSetVariableBlockSizes()`, `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCJACOBI`, `PCPBJACOBI`, `PCBJACOBI`
1950da83c2eSBarry Smith M*/
1960da83c2eSBarry Smith 
197d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_VPBJacobi(PC pc)
198d71ae5a4SJacob Faibussowitsch {
1990da83c2eSBarry Smith   PC_VPBJacobi *jac;
2000da83c2eSBarry Smith 
2010da83c2eSBarry Smith   PetscFunctionBegin;
2020da83c2eSBarry Smith   /*
2030da83c2eSBarry Smith      Creates the private data structure for this preconditioner and
2040da83c2eSBarry Smith      attach it to the PC object.
2050da83c2eSBarry Smith   */
2064dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
2070da83c2eSBarry Smith   pc->data = (void *)jac;
2080da83c2eSBarry Smith 
2090da83c2eSBarry Smith   /*
2100da83c2eSBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
2110da83c2eSBarry Smith      diagonal entries of the matrix for fast preconditioner application.
2120da83c2eSBarry Smith   */
2130da83c2eSBarry Smith   jac->diag = NULL;
2140da83c2eSBarry Smith 
2150da83c2eSBarry Smith   /*
2160da83c2eSBarry Smith       Set the pointers for the functions that are provided above.
2170da83c2eSBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
2180da83c2eSBarry Smith       are called, they will automatically call these functions.  Note we
2190da83c2eSBarry Smith       choose not to provide a couple of these functions since they are
2200da83c2eSBarry Smith       not needed.
2210da83c2eSBarry Smith   */
2220da83c2eSBarry Smith   pc->ops->apply               = PCApply_VPBJacobi;
2230a545947SLisandro Dalcin   pc->ops->applytranspose      = NULL;
2240da83c2eSBarry Smith   pc->ops->setup               = PCSetUp_VPBJacobi;
2250da83c2eSBarry Smith   pc->ops->destroy             = PCDestroy_VPBJacobi;
2260a545947SLisandro Dalcin   pc->ops->setfromoptions      = NULL;
2270a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
2280a545947SLisandro Dalcin   pc->ops->applysymmetricleft  = NULL;
2290a545947SLisandro Dalcin   pc->ops->applysymmetricright = NULL;
2300da83c2eSBarry Smith   PetscFunctionReturn(0);
2310da83c2eSBarry Smith }
232