xref: /petsc/src/ksp/pc/impls/vpbjacobi/vpbjacobi.c (revision 0da83c2e51e31cc9f188e522ae9f307edf38c9c2)
1*0da83c2eSBarry Smith 
2*0da83c2eSBarry Smith /*
3*0da83c2eSBarry Smith    Include files needed for the variable size block PBJacobi preconditioner:
4*0da83c2eSBarry Smith      pcimpl.h - private include file intended for use by all preconditioners
5*0da83c2eSBarry Smith */
6*0da83c2eSBarry Smith 
7*0da83c2eSBarry Smith #include <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
8*0da83c2eSBarry Smith 
9*0da83c2eSBarry Smith /*
10*0da83c2eSBarry Smith    Private context (data structure) for the VPBJacobi preconditioner.
11*0da83c2eSBarry Smith */
12*0da83c2eSBarry Smith typedef struct {
13*0da83c2eSBarry Smith   MatScalar *diag;
14*0da83c2eSBarry Smith } PC_VPBJacobi;
15*0da83c2eSBarry Smith 
16*0da83c2eSBarry Smith 
17*0da83c2eSBarry Smith static PetscErrorCode PCApply_VPBJacobi(PC pc,Vec x,Vec y)
18*0da83c2eSBarry Smith {
19*0da83c2eSBarry Smith   PC_VPBJacobi      *jac = (PC_VPBJacobi*)pc->data;
20*0da83c2eSBarry Smith   PetscErrorCode    ierr;
21*0da83c2eSBarry Smith   PetscInt          i,ncnt = 0;
22*0da83c2eSBarry Smith   const MatScalar   *diag = jac->diag;
23*0da83c2eSBarry Smith   PetscInt          ib,jb,bs;
24*0da83c2eSBarry Smith   const PetscScalar *xx;
25*0da83c2eSBarry Smith   PetscScalar       *yy,x0,x1,x2,x3,x4,x5,x6;
26*0da83c2eSBarry Smith   PetscInt          nblocks;
27*0da83c2eSBarry Smith   const PetscInt    *bsizes;
28*0da83c2eSBarry Smith 
29*0da83c2eSBarry Smith   PetscFunctionBegin;
30*0da83c2eSBarry Smith   ierr = MatGetVariableBlockSizes(pc->pmat,&nblocks,&bsizes);CHKERRQ(ierr);
31*0da83c2eSBarry Smith   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
32*0da83c2eSBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
33*0da83c2eSBarry Smith   for (i=0; i<nblocks; i++) {
34*0da83c2eSBarry Smith     bs = bsizes[i];
35*0da83c2eSBarry Smith     switch (bs) {
36*0da83c2eSBarry Smith     case 1:
37*0da83c2eSBarry Smith       yy[ncnt] = *diag*xx[ncnt];
38*0da83c2eSBarry Smith       break;
39*0da83c2eSBarry Smith     case 2:
40*0da83c2eSBarry Smith       x0         = xx[ncnt]; x1 = xx[ncnt+1];
41*0da83c2eSBarry Smith       yy[ncnt]   = diag[0]*x0 + diag[2]*x1;
42*0da83c2eSBarry Smith       yy[ncnt+1] = diag[1]*x0 + diag[3]*x1;
43*0da83c2eSBarry Smith       break;
44*0da83c2eSBarry Smith     case 3:
45*0da83c2eSBarry Smith       x0 = xx[ncnt]; x1 = xx[ncnt+1]; x2 = xx[ncnt+2];
46*0da83c2eSBarry Smith       yy[ncnt]   = diag[0]*x0 + diag[3]*x1 + diag[6]*x2;
47*0da83c2eSBarry Smith       yy[ncnt+1] = diag[1]*x0 + diag[4]*x1 + diag[7]*x2;
48*0da83c2eSBarry Smith       yy[ncnt+2] = diag[2]*x0 + diag[5]*x1 + diag[8]*x2;
49*0da83c2eSBarry Smith       break;
50*0da83c2eSBarry Smith     case 4:
51*0da83c2eSBarry Smith       x0 = xx[ncnt]; x1 = xx[ncnt+1]; x2 = xx[ncnt+2]; x3 = xx[ncnt+3];
52*0da83c2eSBarry Smith       yy[ncnt]   = diag[0]*x0 + diag[4]*x1 + diag[8]*x2  + diag[12]*x3;
53*0da83c2eSBarry Smith       yy[ncnt+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2  + diag[13]*x3;
54*0da83c2eSBarry Smith       yy[ncnt+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3;
55*0da83c2eSBarry Smith       yy[ncnt+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3;
56*0da83c2eSBarry Smith       break;
57*0da83c2eSBarry Smith     case 5:
58*0da83c2eSBarry Smith       x0 = xx[ncnt]; x1 = xx[ncnt+1]; x2 = xx[ncnt+2]; x3 = xx[ncnt+3]; x4 = xx[ncnt+4];
59*0da83c2eSBarry Smith       yy[ncnt]   = diag[0]*x0 + diag[5]*x1 + diag[10]*x2  + diag[15]*x3 + diag[20]*x4;
60*0da83c2eSBarry Smith       yy[ncnt+1] = diag[1]*x0 + diag[6]*x1 + diag[11]*x2  + diag[16]*x3 + diag[21]*x4;
61*0da83c2eSBarry Smith       yy[ncnt+2] = diag[2]*x0 + diag[7]*x1 + diag[12]*x2 + diag[17]*x3 + diag[22]*x4;
62*0da83c2eSBarry Smith       yy[ncnt+3] = diag[3]*x0 + diag[8]*x1 + diag[13]*x2 + diag[18]*x3 + diag[23]*x4;
63*0da83c2eSBarry Smith       yy[ncnt+4] = diag[4]*x0 + diag[9]*x1 + diag[14]*x2 + diag[19]*x3 + diag[24]*x4;
64*0da83c2eSBarry Smith       break;
65*0da83c2eSBarry Smith     case 6:
66*0da83c2eSBarry Smith       x0 = xx[ncnt]; x1 = xx[ncnt+1]; x2 = xx[ncnt+2]; x3 = xx[ncnt+3]; x4 = xx[ncnt+4]; x5 = xx[ncnt+5];
67*0da83c2eSBarry Smith       yy[ncnt]   = diag[0]*x0 + diag[6]*x1  + diag[12]*x2  + diag[18]*x3 + diag[24]*x4 + diag[30]*x5;
68*0da83c2eSBarry Smith       yy[ncnt+1] = diag[1]*x0 + diag[7]*x1  + diag[13]*x2  + diag[19]*x3 + diag[25]*x4 + diag[31]*x5;
69*0da83c2eSBarry Smith       yy[ncnt+2] = diag[2]*x0 + diag[8]*x1  + diag[14]*x2  + diag[20]*x3 + diag[26]*x4 + diag[32]*x5;
70*0da83c2eSBarry Smith       yy[ncnt+3] = diag[3]*x0 + diag[9]*x1  + diag[15]*x2  + diag[21]*x3 + diag[27]*x4 + diag[33]*x5;
71*0da83c2eSBarry Smith       yy[ncnt+4] = diag[4]*x0 + diag[10]*x1 + diag[16]*x2  + diag[22]*x3 + diag[28]*x4 + diag[34]*x5;
72*0da83c2eSBarry Smith       yy[ncnt+5] = diag[5]*x0 + diag[11]*x1 + diag[17]*x2  + diag[23]*x3 + diag[29]*x4 + diag[35]*x5;
73*0da83c2eSBarry Smith       break;
74*0da83c2eSBarry Smith     case 7:
75*0da83c2eSBarry 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];
76*0da83c2eSBarry Smith       yy[ncnt]   = diag[0]*x0 + diag[7]*x1  + diag[14]*x2  + diag[21]*x3 + diag[28]*x4 + diag[35]*x5 + diag[42]*x6;
77*0da83c2eSBarry 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;
78*0da83c2eSBarry 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;
79*0da83c2eSBarry 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;
80*0da83c2eSBarry 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;
81*0da83c2eSBarry 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;
82*0da83c2eSBarry 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;
83*0da83c2eSBarry Smith       break;
84*0da83c2eSBarry Smith     default:
85*0da83c2eSBarry Smith       for (ib=0; ib<bs; ib++){
86*0da83c2eSBarry Smith         PetscScalar rowsum = 0;
87*0da83c2eSBarry Smith         for (jb=0; jb<bs; jb++){
88*0da83c2eSBarry Smith           rowsum += diag[ib+jb*bs] * xx[ncnt+jb];
89*0da83c2eSBarry Smith         }
90*0da83c2eSBarry Smith         yy[ncnt+ib] = rowsum;
91*0da83c2eSBarry Smith       }
92*0da83c2eSBarry Smith     }
93*0da83c2eSBarry Smith     ncnt += bsizes[i];
94*0da83c2eSBarry Smith     diag += bsizes[i]*bsizes[i];
95*0da83c2eSBarry Smith   }
96*0da83c2eSBarry Smith   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
97*0da83c2eSBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
98*0da83c2eSBarry Smith   PetscFunctionReturn(0);
99*0da83c2eSBarry Smith }
100*0da83c2eSBarry Smith 
101*0da83c2eSBarry Smith 
102*0da83c2eSBarry Smith 
103*0da83c2eSBarry Smith /* -------------------------------------------------------------------------- */
104*0da83c2eSBarry Smith static PetscErrorCode PCSetUp_VPBJacobi(PC pc)
105*0da83c2eSBarry Smith {
106*0da83c2eSBarry Smith   PC_VPBJacobi    *jac = (PC_VPBJacobi*)pc->data;
107*0da83c2eSBarry Smith   PetscErrorCode ierr;
108*0da83c2eSBarry Smith   Mat            A = pc->pmat;
109*0da83c2eSBarry Smith   MatFactorError err;
110*0da83c2eSBarry Smith   PetscInt       i,nsize = 0,nlocal;
111*0da83c2eSBarry Smith   PetscInt       nblocks;
112*0da83c2eSBarry Smith   const PetscInt *bsizes;
113*0da83c2eSBarry Smith 
114*0da83c2eSBarry Smith   PetscFunctionBegin;
115*0da83c2eSBarry Smith   ierr = MatGetVariableBlockSizes(pc->pmat,&nblocks,&bsizes);CHKERRQ(ierr);
116*0da83c2eSBarry Smith   ierr = MatGetLocalSize(pc->pmat,&nlocal,NULL);CHKERRQ(ierr);
117*0da83c2eSBarry Smith   if (nlocal && !nblocks) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatSetVariableBlockSizes() before using PCVPBJACOBI");
118*0da83c2eSBarry Smith   if (!jac->diag) {
119*0da83c2eSBarry Smith     for (i=0; i<nblocks; i++) nsize += bsizes[i]*bsizes[i];
120*0da83c2eSBarry Smith     ierr = PetscMalloc1(nsize,&jac->diag);CHKERRQ(ierr);
121*0da83c2eSBarry Smith   }
122*0da83c2eSBarry Smith   ierr = MatInvertVariableBlockDiagonal(A,nblocks,bsizes,jac->diag);CHKERRQ(ierr);
123*0da83c2eSBarry Smith   ierr = MatFactorGetError(A,&err);CHKERRQ(ierr);
124*0da83c2eSBarry Smith   if (err) pc->failedreason = (PCFailedReason)err;
125*0da83c2eSBarry Smith   pc->ops->apply = PCApply_VPBJacobi;
126*0da83c2eSBarry Smith   PetscFunctionReturn(0);
127*0da83c2eSBarry Smith }
128*0da83c2eSBarry Smith /* -------------------------------------------------------------------------- */
129*0da83c2eSBarry Smith static PetscErrorCode PCDestroy_VPBJacobi(PC pc)
130*0da83c2eSBarry Smith {
131*0da83c2eSBarry Smith   PC_VPBJacobi    *jac = (PC_VPBJacobi*)pc->data;
132*0da83c2eSBarry Smith   PetscErrorCode ierr;
133*0da83c2eSBarry Smith 
134*0da83c2eSBarry Smith   PetscFunctionBegin;
135*0da83c2eSBarry Smith   /*
136*0da83c2eSBarry Smith       Free the private data structure that was hanging off the PC
137*0da83c2eSBarry Smith   */
138*0da83c2eSBarry Smith   ierr = PetscFree(jac->diag);CHKERRQ(ierr);
139*0da83c2eSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
140*0da83c2eSBarry Smith   PetscFunctionReturn(0);
141*0da83c2eSBarry Smith }
142*0da83c2eSBarry Smith 
143*0da83c2eSBarry Smith /* -------------------------------------------------------------------------- */
144*0da83c2eSBarry Smith /*MC
145*0da83c2eSBarry Smith      PCVPBJACOBI - Variable size point block Jacobi preconditioner
146*0da83c2eSBarry Smith 
147*0da83c2eSBarry Smith 
148*0da83c2eSBarry Smith    Notes:
149*0da83c2eSBarry Smith     See PCJACOBI for point Jacobi preconditioning, PCPBJACOBI for fixed point block size, and PCBJACOBI for large size blocks
150*0da83c2eSBarry Smith 
151*0da83c2eSBarry Smith    This works for AIJ matrices
152*0da83c2eSBarry Smith 
153*0da83c2eSBarry Smith    Uses dense LU factorization with partial pivoting to invert the blocks; if a zero pivot
154*0da83c2eSBarry Smith    is detected a PETSc error is generated.
155*0da83c2eSBarry Smith 
156*0da83c2eSBarry Smith    One must call MatSetVariableBlockSizes() to use this preconditioner
157*0da83c2eSBarry Smith    Developer Notes:
158*0da83c2eSBarry Smith     This should support the PCSetErrorIfFailure() flag set to PETSC_TRUE to allow
159*0da83c2eSBarry Smith    the factorization to continue even after a zero pivot is found resulting in a Nan and hence
160*0da83c2eSBarry Smith    terminating KSP with a KSP_DIVERGED_NANORIF allowing
161*0da83c2eSBarry Smith    a nonlinear solver/ODE integrator to recover without stopping the program as currently happens.
162*0da83c2eSBarry Smith 
163*0da83c2eSBarry Smith    Perhaps should provide an option that allows generation of a valid preconditioner
164*0da83c2eSBarry Smith    even if a block is singular as the PCJACOBI does.
165*0da83c2eSBarry Smith 
166*0da83c2eSBarry Smith    Level: beginner
167*0da83c2eSBarry Smith 
168*0da83c2eSBarry Smith   Concepts: variable point block Jacobi
169*0da83c2eSBarry Smith 
170*0da83c2eSBarry Smith .seealso:  MatSetVariableBlockSizes(), PCCreate(), PCSetType(), PCType (for list of available types), PC, PCJACOBI
171*0da83c2eSBarry Smith 
172*0da83c2eSBarry Smith M*/
173*0da83c2eSBarry Smith 
174*0da83c2eSBarry Smith PETSC_EXTERN PetscErrorCode PCCreate_VPBJacobi(PC pc)
175*0da83c2eSBarry Smith {
176*0da83c2eSBarry Smith   PC_VPBJacobi   *jac;
177*0da83c2eSBarry Smith   PetscErrorCode ierr;
178*0da83c2eSBarry Smith 
179*0da83c2eSBarry Smith   PetscFunctionBegin;
180*0da83c2eSBarry Smith   /*
181*0da83c2eSBarry Smith      Creates the private data structure for this preconditioner and
182*0da83c2eSBarry Smith      attach it to the PC object.
183*0da83c2eSBarry Smith   */
184*0da83c2eSBarry Smith   ierr     = PetscNewLog(pc,&jac);CHKERRQ(ierr);
185*0da83c2eSBarry Smith   pc->data = (void*)jac;
186*0da83c2eSBarry Smith 
187*0da83c2eSBarry Smith   /*
188*0da83c2eSBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
189*0da83c2eSBarry Smith      diagonal entries of the matrix for fast preconditioner application.
190*0da83c2eSBarry Smith   */
191*0da83c2eSBarry Smith   jac->diag = NULL;
192*0da83c2eSBarry Smith 
193*0da83c2eSBarry Smith   /*
194*0da83c2eSBarry Smith       Set the pointers for the functions that are provided above.
195*0da83c2eSBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
196*0da83c2eSBarry Smith       are called, they will automatically call these functions.  Note we
197*0da83c2eSBarry Smith       choose not to provide a couple of these functions since they are
198*0da83c2eSBarry Smith       not needed.
199*0da83c2eSBarry Smith   */
200*0da83c2eSBarry Smith   pc->ops->apply               = PCApply_VPBJacobi;
201*0da83c2eSBarry Smith   pc->ops->applytranspose      = 0;
202*0da83c2eSBarry Smith   pc->ops->setup               = PCSetUp_VPBJacobi;
203*0da83c2eSBarry Smith   pc->ops->destroy             = PCDestroy_VPBJacobi;
204*0da83c2eSBarry Smith   pc->ops->setfromoptions      = 0;
205*0da83c2eSBarry Smith   pc->ops->applyrichardson     = 0;
206*0da83c2eSBarry Smith   pc->ops->applysymmetricleft  = 0;
207*0da83c2eSBarry Smith   pc->ops->applysymmetricright = 0;
208*0da83c2eSBarry Smith   PetscFunctionReturn(0);
209*0da83c2eSBarry Smith }
210*0da83c2eSBarry Smith 
211*0da83c2eSBarry Smith 
212