Actual source code: pbjacobi.c
2: /*
3: Include files needed for the PBJacobi preconditioner:
4: pcimpl.h - private include file intended for use by all preconditioners
5: */
7: #include <private/matimpl.h>
8: #include <private/pcimpl.h> /*I "petscpc.h" I*/
10: /*
11: Private context (data structure) for the PBJacobi preconditioner.
12: */
13: typedef struct {
14: MatScalar *diag;
15: PetscInt bs,mbs;
16: } PC_PBJacobi;
21: static PetscErrorCode PCApply_PBJacobi_1(PC pc,Vec x,Vec y)
22: {
23: PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
24: PetscErrorCode ierr;
25: PetscInt i,m = jac->mbs;
26: const MatScalar *diag = jac->diag;
27: const PetscScalar *xx;
28: PetscScalar *yy;
31: VecGetArrayRead(x,&xx);
32: VecGetArray(y,&yy);
33: for (i=0; i<m; i++) {
34: yy[i] = diag[i]*xx[i];
35: }
36: VecRestoreArrayRead(x,&xx);
37: VecRestoreArray(y,&yy);
38: PetscLogFlops(2.0*m);
39: return(0);
40: }
44: static PetscErrorCode PCApply_PBJacobi_2(PC pc,Vec x,Vec y)
45: {
46: PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
47: PetscErrorCode ierr;
48: PetscInt i,m = jac->mbs;
49: const MatScalar *diag = jac->diag;
50: PetscScalar x0,x1,*xx,*yy;
51:
53: VecGetArray(x,&xx);
54: VecGetArray(y,&yy);
55: for (i=0; i<m; i++) {
56: x0 = xx[2*i]; x1 = xx[2*i+1];
57: yy[2*i] = diag[0]*x0 + diag[2]*x1;
58: yy[2*i+1] = diag[1]*x0 + diag[3]*x1;
59: diag += 4;
60: }
61: VecRestoreArray(x,&xx);
62: VecRestoreArray(y,&yy);
63: PetscLogFlops(6.0*m);
64: return(0);
65: }
68: static PetscErrorCode PCApply_PBJacobi_3(PC pc,Vec x,Vec y)
69: {
70: PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
71: PetscErrorCode ierr;
72: PetscInt i,m = jac->mbs;
73: const MatScalar *diag = jac->diag;
74: PetscScalar x0,x1,x2,*xx,*yy;
75:
77: VecGetArray(x,&xx);
78: VecGetArray(y,&yy);
79: for (i=0; i<m; i++) {
80: x0 = xx[3*i]; x1 = xx[3*i+1]; x2 = xx[3*i+2];
81: yy[3*i] = diag[0]*x0 + diag[3]*x1 + diag[6]*x2;
82: yy[3*i+1] = diag[1]*x0 + diag[4]*x1 + diag[7]*x2;
83: yy[3*i+2] = diag[2]*x0 + diag[5]*x1 + diag[8]*x2;
84: diag += 9;
85: }
86: VecRestoreArray(x,&xx);
87: VecRestoreArray(y,&yy);
88: PetscLogFlops(15.0*m);
89: return(0);
90: }
93: static PetscErrorCode PCApply_PBJacobi_4(PC pc,Vec x,Vec y)
94: {
95: PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
96: PetscErrorCode ierr;
97: PetscInt i,m = jac->mbs;
98: const MatScalar *diag = jac->diag;
99: PetscScalar x0,x1,x2,x3,*xx,*yy;
100:
102: VecGetArray(x,&xx);
103: VecGetArray(y,&yy);
104: for (i=0; i<m; i++) {
105: x0 = xx[4*i]; x1 = xx[4*i+1]; x2 = xx[4*i+2]; x3 = xx[4*i+3];
106: yy[4*i] = diag[0]*x0 + diag[4]*x1 + diag[8]*x2 + diag[12]*x3;
107: yy[4*i+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2 + diag[13]*x3;
108: yy[4*i+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3;
109: yy[4*i+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3;
110: diag += 16;
111: }
112: VecRestoreArray(x,&xx);
113: VecRestoreArray(y,&yy);
114: PetscLogFlops(28.0*m);
115: return(0);
116: }
119: static PetscErrorCode PCApply_PBJacobi_5(PC pc,Vec x,Vec y)
120: {
121: PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
122: PetscErrorCode ierr;
123: PetscInt i,m = jac->mbs;
124: const MatScalar *diag = jac->diag;
125: PetscScalar x0,x1,x2,x3,x4,*xx,*yy;
126:
128: VecGetArray(x,&xx);
129: VecGetArray(y,&yy);
130: for (i=0; i<m; i++) {
131: x0 = xx[5*i]; x1 = xx[5*i+1]; x2 = xx[5*i+2]; x3 = xx[5*i+3]; x4 = xx[5*i+4];
132: yy[5*i] = diag[0]*x0 + diag[5]*x1 + diag[10]*x2 + diag[15]*x3 + diag[20]*x4;
133: yy[5*i+1] = diag[1]*x0 + diag[6]*x1 + diag[11]*x2 + diag[16]*x3 + diag[21]*x4;
134: yy[5*i+2] = diag[2]*x0 + diag[7]*x1 + diag[12]*x2 + diag[17]*x3 + diag[22]*x4;
135: yy[5*i+3] = diag[3]*x0 + diag[8]*x1 + diag[13]*x2 + diag[18]*x3 + diag[23]*x4;
136: yy[5*i+4] = diag[4]*x0 + diag[9]*x1 + diag[14]*x2 + diag[19]*x3 + diag[24]*x4;
137: diag += 25;
138: }
139: VecRestoreArray(x,&xx);
140: VecRestoreArray(y,&yy);
141: PetscLogFlops(45.0*m);
142: return(0);
143: }
146: static PetscErrorCode PCApply_PBJacobi_6(PC pc,Vec x,Vec y)
147: {
148: PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
149: PetscErrorCode ierr;
150: PetscInt i,m = jac->mbs;
151: const MatScalar *diag = jac->diag;
152: PetscScalar x0,x1,x2,x3,x4,x5,*xx,*yy;
153:
155: VecGetArray(x,&xx);
156: VecGetArray(y,&yy);
157: for (i=0; i<m; i++) {
158: x0 = xx[6*i]; x1 = xx[6*i+1]; x2 = xx[6*i+2]; x3 = xx[6*i+3]; x4 = xx[6*i+4]; x5 = xx[6*i+5];
159: yy[6*i] = diag[0]*x0 + diag[6]*x1 + diag[12]*x2 + diag[18]*x3 + diag[24]*x4 + diag[30]*x5;
160: yy[6*i+1] = diag[1]*x0 + diag[7]*x1 + diag[13]*x2 + diag[19]*x3 + diag[25]*x4 + diag[31]*x5;
161: yy[6*i+2] = diag[2]*x0 + diag[8]*x1 + diag[14]*x2 + diag[20]*x3 + diag[26]*x4 + diag[32]*x5;
162: yy[6*i+3] = diag[3]*x0 + diag[9]*x1 + diag[15]*x2 + diag[21]*x3 + diag[27]*x4 + diag[33]*x5;
163: yy[6*i+4] = diag[4]*x0 + diag[10]*x1 + diag[16]*x2 + diag[22]*x3 + diag[28]*x4 + diag[34]*x5;
164: yy[6*i+5] = diag[5]*x0 + diag[11]*x1 + diag[17]*x2 + diag[23]*x3 + diag[29]*x4 + diag[35]*x5;
165: diag += 36;
166: }
167: VecRestoreArray(x,&xx);
168: VecRestoreArray(y,&yy);
169: PetscLogFlops(66.0*m);
170: return(0);
171: }
172: /* -------------------------------------------------------------------------- */
175: static PetscErrorCode PCSetUp_PBJacobi(PC pc)
176: {
177: PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
179: Mat A = pc->pmat;
182: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supported only for square matrices and square storage");
184: MatInvertBlockDiagonal(A,&jac->diag);
185: jac->bs = A->rmap->bs;
186: jac->mbs = A->rmap->n/A->rmap->bs;
187: switch (jac->bs){
188: case 1:
189: pc->ops->apply = PCApply_PBJacobi_1;
190: break;
191: case 2:
192: pc->ops->apply = PCApply_PBJacobi_2;
193: break;
194: case 3:
195: pc->ops->apply = PCApply_PBJacobi_3;
196: break;
197: case 4:
198: pc->ops->apply = PCApply_PBJacobi_4;
199: break;
200: case 5:
201: pc->ops->apply = PCApply_PBJacobi_5;
202: break;
203: case 6:
204: pc->ops->apply = PCApply_PBJacobi_6;
205: break;
206: default:
207: SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"not supported for block size %D",jac->bs);
208: }
210: return(0);
211: }
212: /* -------------------------------------------------------------------------- */
215: static PetscErrorCode PCDestroy_PBJacobi(PC pc)
216: {
220: /*
221: Free the private data structure that was hanging off the PC
222: */
223: PetscFree(pc->data);
224: return(0);
225: }
226: /* -------------------------------------------------------------------------- */
227: /*MC
228: PCPBJACOBI - Point block Jacobi
230: Level: beginner
232: Concepts: point block Jacobi
235: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC
237: M*/
242: PetscErrorCode PCCreate_PBJacobi(PC pc)
243: {
244: PC_PBJacobi *jac;
249: /*
250: Creates the private data structure for this preconditioner and
251: attach it to the PC object.
252: */
253: PetscNewLog(pc,PC_PBJacobi,&jac);
254: pc->data = (void*)jac;
256: /*
257: Initialize the pointers to vectors to ZERO; these will be used to store
258: diagonal entries of the matrix for fast preconditioner application.
259: */
260: jac->diag = 0;
262: /*
263: Set the pointers for the functions that are provided above.
264: Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
265: are called, they will automatically call these functions. Note we
266: choose not to provide a couple of these functions since they are
267: not needed.
268: */
269: pc->ops->apply = 0; /*set depending on the block size */
270: pc->ops->applytranspose = 0;
271: pc->ops->setup = PCSetUp_PBJacobi;
272: pc->ops->destroy = PCDestroy_PBJacobi;
273: pc->ops->setfromoptions = 0;
274: pc->ops->view = 0;
275: pc->ops->applyrichardson = 0;
276: pc->ops->applysymmetricleft = 0;
277: pc->ops->applysymmetricright = 0;
278: return(0);
279: }