xref: /petsc/src/ksp/pc/impls/vpbjacobi/vpbjacobi.c (revision 8a9c020ef9386a9599af7d232cee023eb5ec42c5)
10da83c2eSBarry Smith 
20da83c2eSBarry Smith /*
30da83c2eSBarry Smith    Include files needed for the variable size block PBJacobi preconditioner:
40da83c2eSBarry Smith      pcimpl.h - private include file intended for use by all preconditioners
50da83c2eSBarry Smith */
60da83c2eSBarry Smith 
70da83c2eSBarry Smith #include <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
80da83c2eSBarry Smith 
90da83c2eSBarry Smith /*
100da83c2eSBarry Smith    Private context (data structure) for the VPBJacobi preconditioner.
110da83c2eSBarry Smith */
120da83c2eSBarry Smith typedef struct {
130da83c2eSBarry Smith   MatScalar *diag;
140da83c2eSBarry Smith } PC_VPBJacobi;
150da83c2eSBarry Smith 
160da83c2eSBarry Smith static PetscErrorCode PCApply_VPBJacobi(PC pc,Vec x,Vec y)
170da83c2eSBarry Smith {
180da83c2eSBarry Smith   PC_VPBJacobi      *jac = (PC_VPBJacobi*)pc->data;
190da83c2eSBarry Smith   PetscInt          i,ncnt = 0;
200da83c2eSBarry Smith   const MatScalar   *diag = jac->diag;
210da83c2eSBarry Smith   PetscInt          ib,jb,bs;
220da83c2eSBarry Smith   const PetscScalar *xx;
230da83c2eSBarry Smith   PetscScalar       *yy,x0,x1,x2,x3,x4,x5,x6;
240da83c2eSBarry Smith   PetscInt          nblocks;
250da83c2eSBarry Smith   const PetscInt    *bsizes;
260da83c2eSBarry Smith 
270da83c2eSBarry Smith   PetscFunctionBegin;
289566063dSJacob Faibussowitsch   PetscCall(MatGetVariableBlockSizes(pc->pmat,&nblocks,&bsizes));
299566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&xx));
309566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y,&yy));
310da83c2eSBarry Smith   for (i=0; i<nblocks; i++) {
320da83c2eSBarry Smith     bs = bsizes[i];
330da83c2eSBarry Smith     switch (bs) {
340da83c2eSBarry Smith     case 1:
350da83c2eSBarry Smith       yy[ncnt] = *diag*xx[ncnt];
360da83c2eSBarry Smith       break;
370da83c2eSBarry Smith     case 2:
380da83c2eSBarry Smith       x0         = xx[ncnt]; x1 = xx[ncnt+1];
390da83c2eSBarry Smith       yy[ncnt]   = diag[0]*x0 + diag[2]*x1;
400da83c2eSBarry Smith       yy[ncnt+1] = diag[1]*x0 + diag[3]*x1;
410da83c2eSBarry Smith       break;
420da83c2eSBarry Smith     case 3:
430da83c2eSBarry Smith       x0 = xx[ncnt]; x1 = xx[ncnt+1]; x2 = xx[ncnt+2];
440da83c2eSBarry Smith       yy[ncnt]   = diag[0]*x0 + diag[3]*x1 + diag[6]*x2;
450da83c2eSBarry Smith       yy[ncnt+1] = diag[1]*x0 + diag[4]*x1 + diag[7]*x2;
460da83c2eSBarry Smith       yy[ncnt+2] = diag[2]*x0 + diag[5]*x1 + diag[8]*x2;
470da83c2eSBarry Smith       break;
480da83c2eSBarry Smith     case 4:
490da83c2eSBarry Smith       x0 = xx[ncnt]; x1 = xx[ncnt+1]; x2 = xx[ncnt+2]; x3 = xx[ncnt+3];
500da83c2eSBarry Smith       yy[ncnt]   = diag[0]*x0 + diag[4]*x1 + diag[8]*x2  + diag[12]*x3;
510da83c2eSBarry Smith       yy[ncnt+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2  + diag[13]*x3;
520da83c2eSBarry Smith       yy[ncnt+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3;
530da83c2eSBarry Smith       yy[ncnt+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3;
540da83c2eSBarry Smith       break;
550da83c2eSBarry Smith     case 5:
560da83c2eSBarry Smith       x0 = xx[ncnt]; x1 = xx[ncnt+1]; x2 = xx[ncnt+2]; x3 = xx[ncnt+3]; x4 = xx[ncnt+4];
570da83c2eSBarry Smith       yy[ncnt]   = diag[0]*x0 + diag[5]*x1 + diag[10]*x2  + diag[15]*x3 + diag[20]*x4;
580da83c2eSBarry Smith       yy[ncnt+1] = diag[1]*x0 + diag[6]*x1 + diag[11]*x2  + diag[16]*x3 + diag[21]*x4;
590da83c2eSBarry Smith       yy[ncnt+2] = diag[2]*x0 + diag[7]*x1 + diag[12]*x2 + diag[17]*x3 + diag[22]*x4;
600da83c2eSBarry Smith       yy[ncnt+3] = diag[3]*x0 + diag[8]*x1 + diag[13]*x2 + diag[18]*x3 + diag[23]*x4;
610da83c2eSBarry Smith       yy[ncnt+4] = diag[4]*x0 + diag[9]*x1 + diag[14]*x2 + diag[19]*x3 + diag[24]*x4;
620da83c2eSBarry Smith       break;
630da83c2eSBarry Smith     case 6:
640da83c2eSBarry Smith       x0 = xx[ncnt]; x1 = xx[ncnt+1]; x2 = xx[ncnt+2]; x3 = xx[ncnt+3]; x4 = xx[ncnt+4]; x5 = xx[ncnt+5];
650da83c2eSBarry Smith       yy[ncnt]   = diag[0]*x0 + diag[6]*x1  + diag[12]*x2  + diag[18]*x3 + diag[24]*x4 + diag[30]*x5;
660da83c2eSBarry Smith       yy[ncnt+1] = diag[1]*x0 + diag[7]*x1  + diag[13]*x2  + diag[19]*x3 + diag[25]*x4 + diag[31]*x5;
670da83c2eSBarry Smith       yy[ncnt+2] = diag[2]*x0 + diag[8]*x1  + diag[14]*x2  + diag[20]*x3 + diag[26]*x4 + diag[32]*x5;
680da83c2eSBarry Smith       yy[ncnt+3] = diag[3]*x0 + diag[9]*x1  + diag[15]*x2  + diag[21]*x3 + diag[27]*x4 + diag[33]*x5;
690da83c2eSBarry Smith       yy[ncnt+4] = diag[4]*x0 + diag[10]*x1 + diag[16]*x2  + diag[22]*x3 + diag[28]*x4 + diag[34]*x5;
700da83c2eSBarry Smith       yy[ncnt+5] = diag[5]*x0 + diag[11]*x1 + diag[17]*x2  + diag[23]*x3 + diag[29]*x4 + diag[35]*x5;
710da83c2eSBarry Smith       break;
720da83c2eSBarry Smith     case 7:
730da83c2eSBarry Smith       x0 = xx[ncnt]; x1 = xx[ncnt+1]; x2 = xx[ncnt+2]; x3 = xx[ncnt+3]; x4 = xx[ncnt+4]; x5 = xx[ncnt+5]; x6 = xx[ncnt+6];
740da83c2eSBarry Smith       yy[ncnt]   = diag[0]*x0 + diag[7]*x1  + diag[14]*x2  + diag[21]*x3 + diag[28]*x4 + diag[35]*x5 + diag[42]*x6;
750da83c2eSBarry 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;
760da83c2eSBarry 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;
770da83c2eSBarry 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;
780da83c2eSBarry 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;
790da83c2eSBarry 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;
800da83c2eSBarry 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;
810da83c2eSBarry Smith       break;
820da83c2eSBarry Smith     default:
830da83c2eSBarry Smith       for (ib=0; ib<bs; ib++) {
840da83c2eSBarry Smith         PetscScalar rowsum = 0;
850da83c2eSBarry Smith         for (jb=0; jb<bs; jb++) {
860da83c2eSBarry Smith           rowsum += diag[ib+jb*bs] * xx[ncnt+jb];
870da83c2eSBarry Smith         }
880da83c2eSBarry Smith         yy[ncnt+ib] = rowsum;
890da83c2eSBarry Smith       }
900da83c2eSBarry Smith     }
910da83c2eSBarry Smith     ncnt += bsizes[i];
920da83c2eSBarry Smith     diag += bsizes[i]*bsizes[i];
930da83c2eSBarry Smith   }
949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&xx));
959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y,&yy));
960da83c2eSBarry Smith   PetscFunctionReturn(0);
970da83c2eSBarry Smith }
980da83c2eSBarry Smith 
990da83c2eSBarry Smith static PetscErrorCode PCSetUp_VPBJacobi(PC pc)
1000da83c2eSBarry Smith {
1010da83c2eSBarry Smith   PC_VPBJacobi    *jac = (PC_VPBJacobi*)pc->data;
1020da83c2eSBarry Smith   Mat            A = pc->pmat;
1030da83c2eSBarry Smith   MatFactorError err;
1040da83c2eSBarry Smith   PetscInt       i,nsize = 0,nlocal;
1050da83c2eSBarry Smith   PetscInt       nblocks;
1060da83c2eSBarry Smith   const PetscInt *bsizes;
1070da83c2eSBarry Smith 
1080da83c2eSBarry Smith   PetscFunctionBegin;
1099566063dSJacob Faibussowitsch   PetscCall(MatGetVariableBlockSizes(pc->pmat,&nblocks,&bsizes));
1109566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat,&nlocal,NULL));
11108401ef6SPierre Jolivet   PetscCheck(!nlocal || nblocks,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatSetVariableBlockSizes() before using PCVPBJACOBI");
1120da83c2eSBarry Smith   if (!jac->diag) {
1130da83c2eSBarry Smith     for (i=0; i<nblocks; i++) nsize += bsizes[i]*bsizes[i];
1149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nsize,&jac->diag));
1150da83c2eSBarry Smith   }
1169566063dSJacob Faibussowitsch   PetscCall(MatInvertVariableBlockDiagonal(A,nblocks,bsizes,jac->diag));
1179566063dSJacob Faibussowitsch   PetscCall(MatFactorGetError(A,&err));
1180da83c2eSBarry Smith   if (err) pc->failedreason = (PCFailedReason)err;
1190da83c2eSBarry Smith   pc->ops->apply = PCApply_VPBJacobi;
1200da83c2eSBarry Smith   PetscFunctionReturn(0);
1210da83c2eSBarry Smith }
122*8a9c020eSBarry Smith 
1230da83c2eSBarry Smith static PetscErrorCode PCDestroy_VPBJacobi(PC pc)
1240da83c2eSBarry Smith {
1250da83c2eSBarry Smith   PC_VPBJacobi    *jac = (PC_VPBJacobi*)pc->data;
1260da83c2eSBarry Smith 
1270da83c2eSBarry Smith   PetscFunctionBegin;
1280da83c2eSBarry Smith   /*
1290da83c2eSBarry Smith       Free the private data structure that was hanging off the PC
1300da83c2eSBarry Smith   */
1319566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->diag));
1329566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1330da83c2eSBarry Smith   PetscFunctionReturn(0);
1340da83c2eSBarry Smith }
1350da83c2eSBarry Smith 
1360da83c2eSBarry Smith /*MC
1370da83c2eSBarry Smith      PCVPBJACOBI - Variable size point block Jacobi preconditioner
1380da83c2eSBarry Smith 
1390da83c2eSBarry Smith    Notes:
1400da83c2eSBarry Smith      See PCJACOBI for point Jacobi preconditioning, PCPBJACOBI for fixed point block size, and PCBJACOBI for large size blocks
1410da83c2eSBarry Smith 
1420da83c2eSBarry Smith      This works for AIJ matrices
1430da83c2eSBarry Smith 
1440da83c2eSBarry Smith      Uses dense LU factorization with partial pivoting to invert the blocks; if a zero pivot
1450da83c2eSBarry Smith      is detected a PETSc error is generated.
1460da83c2eSBarry Smith 
1470da83c2eSBarry Smith      One must call MatSetVariableBlockSizes() to use this preconditioner
14890dfcc32SBarry Smith 
1490da83c2eSBarry Smith    Developer Notes:
1500da83c2eSBarry Smith      This should support the PCSetErrorIfFailure() flag set to PETSC_TRUE to allow
1510da83c2eSBarry Smith      the factorization to continue even after a zero pivot is found resulting in a Nan and hence
1520da83c2eSBarry Smith      terminating KSP with a KSP_DIVERGED_NANORIF allowing
1530da83c2eSBarry Smith      a nonlinear solver/ODE integrator to recover without stopping the program as currently happens.
1540da83c2eSBarry Smith 
1550da83c2eSBarry Smith      Perhaps should provide an option that allows generation of a valid preconditioner
1560da83c2eSBarry Smith      even if a block is singular as the PCJACOBI does.
1570da83c2eSBarry Smith 
1580da83c2eSBarry Smith    Level: beginner
1590da83c2eSBarry Smith 
16090dfcc32SBarry Smith .seealso:  MatSetVariableBlockSizes(), PCCreate(), PCSetType(), PCType (for list of available types), PC, PCJACOBI, PCPBJACOBI, PCBJACOBI
1610da83c2eSBarry Smith 
1620da83c2eSBarry Smith M*/
1630da83c2eSBarry Smith 
1640da83c2eSBarry Smith PETSC_EXTERN PetscErrorCode PCCreate_VPBJacobi(PC pc)
1650da83c2eSBarry Smith {
1660da83c2eSBarry Smith   PC_VPBJacobi   *jac;
1670da83c2eSBarry Smith 
1680da83c2eSBarry Smith   PetscFunctionBegin;
1690da83c2eSBarry Smith   /*
1700da83c2eSBarry Smith      Creates the private data structure for this preconditioner and
1710da83c2eSBarry Smith      attach it to the PC object.
1720da83c2eSBarry Smith   */
1739566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&jac));
1740da83c2eSBarry Smith   pc->data = (void*)jac;
1750da83c2eSBarry Smith 
1760da83c2eSBarry Smith   /*
1770da83c2eSBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
1780da83c2eSBarry Smith      diagonal entries of the matrix for fast preconditioner application.
1790da83c2eSBarry Smith   */
1800da83c2eSBarry Smith   jac->diag = NULL;
1810da83c2eSBarry Smith 
1820da83c2eSBarry Smith   /*
1830da83c2eSBarry Smith       Set the pointers for the functions that are provided above.
1840da83c2eSBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
1850da83c2eSBarry Smith       are called, they will automatically call these functions.  Note we
1860da83c2eSBarry Smith       choose not to provide a couple of these functions since they are
1870da83c2eSBarry Smith       not needed.
1880da83c2eSBarry Smith   */
1890da83c2eSBarry Smith   pc->ops->apply               = PCApply_VPBJacobi;
1900a545947SLisandro Dalcin   pc->ops->applytranspose      = NULL;
1910da83c2eSBarry Smith   pc->ops->setup               = PCSetUp_VPBJacobi;
1920da83c2eSBarry Smith   pc->ops->destroy             = PCDestroy_VPBJacobi;
1930a545947SLisandro Dalcin   pc->ops->setfromoptions      = NULL;
1940a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
1950a545947SLisandro Dalcin   pc->ops->applysymmetricleft  = NULL;
1960a545947SLisandro Dalcin   pc->ops->applysymmetricright = NULL;
1970da83c2eSBarry Smith   PetscFunctionReturn(0);
1980da83c2eSBarry Smith }
199