Actual source code: pbjacobi.c
1: /*
2: Include files needed for the PBJacobi preconditioner:
3: pcimpl.h - private include file intended for use by all preconditioners
4: */
6: #include src/ksp/pc/pcimpl.h
8: /*
9: Private context (data structure) for the PBJacobi preconditioner.
10: */
11: typedef struct {
12: PetscScalar *diag;
13: PetscInt bs,mbs;
14: } PC_PBJacobi;
16: /*
17: Currently only implemented for baij matrices and directly access baij
18: data structures.
19: */
20: #include src/mat/impls/baij/mpi/mpibaij.h
21: #include src/inline/ilu.h
25: static PetscErrorCode PCApply_PBJacobi_2(PC pc,Vec x,Vec y)
26: {
27: PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
29: PetscInt i,m = jac->mbs;
30: PetscScalar *diag = jac->diag,x0,x1,*xx,*yy;
31:
33: VecGetArray(x,&xx);
34: VecGetArray(y,&yy);
35: for (i=0; i<m; i++) {
36: x0 = xx[2*i]; x1 = xx[2*i+1];
37: yy[2*i] = diag[0]*x0 + diag[2]*x1;
38: yy[2*i+1] = diag[1]*x0 + diag[3]*x1;
39: diag += 4;
40: }
41: VecRestoreArray(x,&xx);
42: VecRestoreArray(y,&yy);
43: PetscLogFlops(6*m);
44: return(0);
45: }
48: static PetscErrorCode PCApply_PBJacobi_3(PC pc,Vec x,Vec y)
49: {
50: PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
52: PetscInt i,m = jac->mbs;
53: PetscScalar *diag = jac->diag,x0,x1,x2,*xx,*yy;
54:
56: VecGetArray(x,&xx);
57: VecGetArray(y,&yy);
58: for (i=0; i<m; i++) {
59: x0 = xx[3*i]; x1 = xx[3*i+1]; x2 = xx[3*i+2];
60: yy[3*i] = diag[0]*x0 + diag[3]*x1 + diag[6]*x2;
61: yy[3*i+1] = diag[1]*x0 + diag[4]*x1 + diag[7]*x2;
62: yy[3*i+2] = diag[2]*x0 + diag[5]*x1 + diag[8]*x2;
63: diag += 9;
64: }
65: VecRestoreArray(x,&xx);
66: VecRestoreArray(y,&yy);
67: PetscLogFlops(15*m);
68: return(0);
69: }
72: static PetscErrorCode PCApply_PBJacobi_4(PC pc,Vec x,Vec y)
73: {
74: PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
76: PetscInt i,m = jac->mbs;
77: PetscScalar *diag = jac->diag,x0,x1,x2,x3,*xx,*yy;
78:
80: VecGetArray(x,&xx);
81: VecGetArray(y,&yy);
82: for (i=0; i<m; i++) {
83: x0 = xx[4*i]; x1 = xx[4*i+1]; x2 = xx[4*i+2]; x3 = xx[4*i+3];
84: yy[4*i] = diag[0]*x0 + diag[4]*x1 + diag[8]*x2 + diag[12]*x3;
85: yy[4*i+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2 + diag[13]*x3;
86: yy[4*i+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3;
87: yy[4*i+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3;
88: diag += 16;
89: }
90: VecRestoreArray(x,&xx);
91: VecRestoreArray(y,&yy);
92: PetscLogFlops(28*m);
93: return(0);
94: }
97: static PetscErrorCode PCApply_PBJacobi_5(PC pc,Vec x,Vec y)
98: {
99: PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
101: PetscInt i,m = jac->mbs;
102: PetscScalar *diag = jac->diag,x0,x1,x2,x3,x4,*xx,*yy;
103:
105: VecGetArray(x,&xx);
106: VecGetArray(y,&yy);
107: for (i=0; i<m; i++) {
108: x0 = xx[5*i]; x1 = xx[5*i+1]; x2 = xx[5*i+2]; x3 = xx[5*i+3]; x4 = xx[5*i+4];
109: yy[5*i] = diag[0]*x0 + diag[5]*x1 + diag[10]*x2 + diag[15]*x3 + diag[20]*x4;
110: yy[5*i+1] = diag[1]*x0 + diag[6]*x1 + diag[11]*x2 + diag[16]*x3 + diag[21]*x4;
111: yy[5*i+2] = diag[2]*x0 + diag[7]*x1 + diag[12]*x2 + diag[17]*x3 + diag[22]*x4;
112: yy[5*i+3] = diag[3]*x0 + diag[8]*x1 + diag[13]*x2 + diag[18]*x3 + diag[23]*x4;
113: yy[5*i+4] = diag[4]*x0 + diag[9]*x1 + diag[14]*x2 + diag[19]*x3 + diag[24]*x4;
114: diag += 25;
115: }
116: VecRestoreArray(x,&xx);
117: VecRestoreArray(y,&yy);
118: PetscLogFlops(45*m);
119: return(0);
120: }
121: /* -------------------------------------------------------------------------- */
122: EXTERN PetscErrorCode MatInvertBlockDiagonal_SeqBAIJ(Mat);
125: static PetscErrorCode PCSetUp_PBJacobi(PC pc)
126: {
127: PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
129: PetscMPIInt size;
130: PetscTruth seqbaij,mpibaij,baij;
131: Mat A = pc->pmat;
132: Mat_SeqBAIJ *a;
135: PetscTypeCompare((PetscObject)pc->pmat,MATSEQBAIJ,&seqbaij);
136: PetscTypeCompare((PetscObject)pc->pmat,MATMPIBAIJ,&mpibaij);
137: PetscTypeCompare((PetscObject)pc->pmat,MATBAIJ,&baij);
138: if (!seqbaij && !mpibaij && !baij) {
139: SETERRQ(PETSC_ERR_SUP,"Currently only supports BAIJ matrices");
140: }
141: MPI_Comm_size(pc->comm,&size);
142: if (mpibaij || (baij && (size > 1))) A = ((Mat_MPIBAIJ*)A->data)->A;
143: if (A->m != A->n) SETERRQ(PETSC_ERR_SUP,"Supported only for square matrices and square storage");
145: MatInvertBlockDiagonal_SeqBAIJ(A);
146: a = (Mat_SeqBAIJ*)A->data;
147: jac->diag = a->idiag;
148: jac->bs = A->bs;
149: jac->mbs = a->mbs;
150: switch (jac->bs){
151: case 2:
152: pc->ops->apply = PCApply_PBJacobi_2;
153: break;
154: case 3:
155: pc->ops->apply = PCApply_PBJacobi_3;
156: break;
157: case 4:
158: pc->ops->apply = PCApply_PBJacobi_4;
159: break;
160: case 5:
161: pc->ops->apply = PCApply_PBJacobi_5;
162: break;
163: default:
164: SETERRQ1(PETSC_ERR_SUP,"not supported for block size %D",jac->bs);
165: }
167: return(0);
168: }
169: /* -------------------------------------------------------------------------- */
172: static PetscErrorCode PCDestroy_PBJacobi(PC pc)
173: {
174: PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
178: /*
179: Free the private data structure that was hanging off the PC
180: */
181: PetscFree(jac);
182: return(0);
183: }
184: /* -------------------------------------------------------------------------- */
185: /*MC
186: PCPBJACOBI - Point block Jacobi
188: Level: beginner
190: Concepts: point block Jacobi
192: Notes: Only implemented for the BAIJ matrix formats.
194: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC
196: M*/
201: PetscErrorCode PCCreate_PBJacobi(PC pc)
202: {
203: PC_PBJacobi *jac;
208: /*
209: Creates the private data structure for this preconditioner and
210: attach it to the PC object.
211: */
212: PetscNew(PC_PBJacobi,&jac);
213: pc->data = (void*)jac;
215: /*
216: Logs the memory usage; this is not needed but allows PETSc to
217: monitor how much memory is being used for various purposes.
218: */
219: PetscLogObjectMemory(pc,sizeof(PC_PBJacobi));
221: /*
222: Initialize the pointers to vectors to ZERO; these will be used to store
223: diagonal entries of the matrix for fast preconditioner application.
224: */
225: jac->diag = 0;
227: /*
228: Set the pointers for the functions that are provided above.
229: Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
230: are called, they will automatically call these functions. Note we
231: choose not to provide a couple of these functions since they are
232: not needed.
233: */
234: pc->ops->apply = 0; /*set depending on the block size */
235: pc->ops->applytranspose = 0;
236: pc->ops->setup = PCSetUp_PBJacobi;
237: pc->ops->destroy = PCDestroy_PBJacobi;
238: pc->ops->setfromoptions = 0;
239: pc->ops->view = 0;
240: pc->ops->applyrichardson = 0;
241: pc->ops->applysymmetricleft = 0;
242: pc->ops->applysymmetricright = 0;
243: return(0);
244: }