Actual source code: ilu.h

  1: /*
  2:     Kernels used in sparse ILU (and LU) and in the resulting triangular
  3:  solves. These are for block algorithms where the block sizes are on 
  4:  the order of 2-6+.

  6:     There are TWO versions of the macros below. 
  7:     1) standard for MatScalar == PetscScalar use the standard BLAS for 
  8:        block size larger than 7 and
  9:     2) handcoded Fortran single precision for the matrices, since BLAS
 10:        does not have some arguments in single and some in double.

 12: */

 16:  #include petscblaslapack.h

 18: /*
 19:       These are C kernels,they are contained in 
 20:    src/mat/impls/baij/seq
 21: */

 23: EXTERN PetscErrorCode  LINPACKdgefa(MatScalar*,PetscInt,PetscInt*);
 24: EXTERN PetscErrorCode  LINPACKdgedi(MatScalar*,PetscInt,PetscInt*,MatScalar*);
 25: EXTERN PetscErrorCode  Kernel_A_gets_inverse_A_2(MatScalar*,PetscReal);
 26: EXTERN PetscErrorCode  Kernel_A_gets_inverse_A_3(MatScalar*,PetscReal);

 28: #define Kernel_A_gets_inverse_A_4_nopivot(mat) 0;\
 29: {\
 30:   MatScalar d, di;\
 31: \
 32:   di = mat[0];\
 33:   mat[0] = d = 1.0 / di;\
 34:   mat[4] *= -d;\
 35:   mat[8] *= -d;\
 36:   mat[12] *= -d;\
 37:   mat[1] *= d;\
 38:   mat[2] *= d;\
 39:   mat[3] *= d;\
 40:   mat[5] += mat[4] * mat[1] * di;\
 41:   mat[6] += mat[4] * mat[2] * di;\
 42:   mat[7] += mat[4] * mat[3] * di;\
 43:   mat[9] += mat[8] * mat[1] * di;\
 44:   mat[10] += mat[8] * mat[2] * di;\
 45:   mat[11] += mat[8] * mat[3] * di;\
 46:   mat[13] += mat[12] * mat[1] * di;\
 47:   mat[14] += mat[12] * mat[2] * di;\
 48:   mat[15] += mat[12] * mat[3] * di;\
 49:   di = mat[5];\
 50:   mat[5] = d = 1.0 / di;\
 51:   mat[1] *= -d;\
 52:   mat[9] *= -d;\
 53:   mat[13] *= -d;\
 54:   mat[4] *= d;\
 55:   mat[6] *= d;\
 56:   mat[7] *= d;\
 57:   mat[0] += mat[1] * mat[4] * di;\
 58:   mat[2] += mat[1] * mat[6] * di;\
 59:   mat[3] += mat[1] * mat[7] * di;\
 60:   mat[8] += mat[9] * mat[4] * di;\
 61:   mat[10] += mat[9] * mat[6] * di;\
 62:   mat[11] += mat[9] * mat[7] * di;\
 63:   mat[12] += mat[13] * mat[4] * di;\
 64:   mat[14] += mat[13] * mat[6] * di;\
 65:   mat[15] += mat[13] * mat[7] * di;\
 66:   di = mat[10];\
 67:   mat[10] = d = 1.0 / di;\
 68:   mat[2] *= -d;\
 69:   mat[6] *= -d;\
 70:   mat[14] *= -d;\
 71:   mat[8] *= d;\
 72:   mat[9] *= d;\
 73:   mat[11] *= d;\
 74:   mat[0] += mat[2] * mat[8] * di;\
 75:   mat[1] += mat[2] * mat[9] * di;\
 76:   mat[3] += mat[2] * mat[11] * di;\
 77:   mat[4] += mat[6] * mat[8] * di;\
 78:   mat[5] += mat[6] * mat[9] * di;\
 79:   mat[7] += mat[6] * mat[11] * di;\
 80:   mat[12] += mat[14] * mat[8] * di;\
 81:   mat[13] += mat[14] * mat[9] * di;\
 82:   mat[15] += mat[14] * mat[11] * di;\
 83:   di = mat[15];\
 84:   mat[15] = d = 1.0 / di;\
 85:   mat[3] *= -d;\
 86:   mat[7] *= -d;\
 87:   mat[11] *= -d;\
 88:   mat[12] *= d;\
 89:   mat[13] *= d;\
 90:   mat[14] *= d;\
 91:   mat[0] += mat[3] * mat[12] * di;\
 92:   mat[1] += mat[3] * mat[13] * di;\
 93:   mat[2] += mat[3] * mat[14] * di;\
 94:   mat[4] += mat[7] * mat[12] * di;\
 95:   mat[5] += mat[7] * mat[13] * di;\
 96:   mat[6] += mat[7] * mat[14] * di;\
 97:   mat[8] += mat[11] * mat[12] * di;\
 98:   mat[9] += mat[11] * mat[13] * di;\
 99:   mat[10] += mat[11] * mat[14] * di;\
100: }

102: EXTERN PetscErrorCode Kernel_A_gets_inverse_A_4(MatScalar *,PetscReal);
103: # if defined(PETSC_HAVE_SSE)
104: EXTERN PetscErrorCode Kernel_A_gets_inverse_A_4_SSE(MatScalar *);
105: # endif
106: EXTERN PetscErrorCode Kernel_A_gets_inverse_A_5(MatScalar *,PetscReal);
107: EXTERN PetscErrorCode Kernel_A_gets_inverse_A_6(MatScalar *,PetscReal);
108: EXTERN PetscErrorCode Kernel_A_gets_inverse_A_7(MatScalar *,PetscReal);
109: EXTERN PetscErrorCode Kernel_A_gets_inverse_A_9(MatScalar *,PetscReal);

111: /*
112:     A = inv(A)    A_gets_inverse_A

114:    A      - square bs by bs array stored in column major order
115:    pivots - integer work array of length bs
116:    W      -  bs by 1 work array
117: */
118: #define Kernel_A_gets_inverse_A(bs,A,pivots,W) (LINPACKdgefa((A),(bs),(pivots)) || LINPACKdgedi((A),(bs),(pivots),(W)))

120: /* -----------------------------------------------------------------------*/

122: #if !defined(PETSC_USE_MAT_SINGLE)
123: /*
124:         Version that calls the BLAS directly
125: */
126: /*
127:       A = A * B   A_gets_A_times_B

129:    A, B - square bs by bs arrays stored in column major order
130:    W    - square bs by bs work array

132: */
133: #define Kernel_A_gets_A_times_B(bs,A,B,W) \
134: { \
135:   PetscBLASInt   _bbs;                 \
136:   PetscScalar    _one = 1.0,_zero = 0.0; \
137:   PetscErrorCode _ierr; \
138:   _bbs = PetscBLASIntCast(bs);        \
139:   _PetscMemcpy((W),(A),(bs)*(bs)*sizeof(MatScalar));CHKERRQ(_ierr); \
140:   BLASgemm_("N","N",&(_bbs),&(_bbs),&(_bbs),&_one,(W),&(_bbs),(B),&(_bbs),&_zero,(A),&(_bbs));\
141: }

143: /*

145:     A = A - B * C  A_gets_A_minus_B_times_C 

147:    A, B, C - square bs by bs arrays stored in column major order
148: */
149: #define Kernel_A_gets_A_minus_B_times_C(bs,A,B,C) \
150: { \
151:   PetscScalar  _mone = -1.0,_one = 1.0; \
152:   PetscBLASInt _bbs = PetscBLASIntCast(bs);        \
153:   BLASgemm_("N","N",&(_bbs),&(_bbs),&(_bbs),&_mone,(B),&(_bbs),(C),&(_bbs),&_one,(A),&(_bbs));\
154: }

156: /*

158:     A = A + B^T * C  A_gets_A_plus_Btranspose_times_C 

160:    A, B, C - square bs by bs arrays stored in column major order
161: */
162: #define Kernel_A_gets_A_plus_Btranspose_times_C(bs,A,B,C) \
163: { \
164:   PetscScalar  _one = 1.0; \
165:   PetscBLASInt _bbs = PetscBLASIntCast(bs);        \
166:   BLASgemm_("T","N",&(_bbs),&(_bbs),&(_bbs),&_one,(B),&(_bbs),(C),&(_bbs),&_one,(A),&(_bbs));\
167: }

169: /*
170:     v = v + A^T w  v_gets_v_plus_Atranspose_times_w

172:    v - array of length bs
173:    A - square bs by bs array
174:    w - array of length bs
175: */
176: #define  Kernel_v_gets_v_plus_Atranspose_times_w(bs,v,A,w) \
177: {  \
178:   PetscScalar  _one = 1.0; \
179:   PetscBLASInt _ione = 1, _bbs = PetscBLASIntCast(bs);                        \
180:   BLASgemv_("T",&(_bbs),&(_bbs),&_one,A,&(_bbs),w,&_ione,&_one,v,&_ione); \
181: } 

183: /*
184:     v = v - A w  v_gets_v_minus_A_times_w

186:    v - array of length bs
187:    A - square bs by bs array
188:    w - array of length bs
189: */
190: #define  Kernel_v_gets_v_minus_A_times_w(bs,v,A,w) \
191: {  \
192:   PetscScalar  _mone = -1.0,_one = 1.0; \
193:   PetscBLASInt  _ione = 1,_bbs = PetscBLASIntCast(bs);                        \
194:   BLASgemv_("N",&(_bbs),&(_bbs),&_mone,A,&(_bbs),w,&_ione,&_one,v,&_ione); \
195: }

197: /*
198:     v = v + A w  v_gets_v_plus_A_times_w

200:    v - array of length bs
201:    A - square bs by bs array
202:    w - array of length bs
203: */
204: #define  Kernel_v_gets_v_plus_A_times_w(bs,v,A,w) \
205: {  \
206:   PetscScalar  _one = 1.0; \
207:   PetscBLASInt _ione = 1,_bbs = PetscBLASIntCast(bs);                        \
208:   BLASgemv_("N",&(_bbs),&(_bbs),&_one,A,&(_bbs),w,&_ione,&_one,v,&_ione); \
209: }

211: /*
212:     v = v + A w  v_gets_v_plus_Ar_times_w

214:    v - array of length bs
215:    A - square bs by bs array
216:    w - array of length bs
217: */
218: #define  Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,v,A,w) \
219: {  \
220:   PetscScalar  _one = 1.0; \
221:   PetscBLASInt _ione = 1,_bbs,_bncols; \
222:   _bbs = PetscBLASIntCast(bs); _bncols = PetscBLASIntCast(ncols); \
223:   BLASgemv_("N",&(_bbs),&(_bncols),&_one,A,&(_bbs),v,&_ione,&_one,w,&_ione); \
224: }

226: /*
227:     w = A v   w_gets_A_times_v

229:    v - array of length bs
230:    A - square bs by bs array
231:    w - array of length bs
232: */
233: #define Kernel_w_gets_A_times_v(bs,v,A,w) \
234: {  \
235:   PetscScalar  _zero = 0.0,_one = 1.0; \
236:   PetscBLASInt _ione = 1,_bbs = PetscBLASIntCast(bs);                        \
237:   BLASgemv_("N",&(_bbs),&(_bbs),&_one,A,&(_bbs),v,&_ione,&_zero,w,&_ione); \
238: }

240: /*
241:         z = A*x
242: */
243: #define Kernel_w_gets_Ar_times_v(bs,ncols,x,A,z) \
244: { \
245:   PetscScalar _one = 1.0,_zero = 0.0; \
246:   PetscBLASInt _ione = 1,_bbs,_bncols; \
247:   _bbs = PetscBLASIntCast(bs); _bncols = PetscBLASIntCast(ncols); \
248:   BLASgemv_("N",&(_bbs),&_bncols,&_one,A,&(_bbs),x,&_ione,&_zero,z,&_ione); \
249: }

251: /*
252:         z = trans(A)*x
253: */
254: #define Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,A,z) \
255: { \
256:   PetscScalar  _one = 1.0; \
257:   PetscBLASInt _ione = 1,_bbs,_bncols;\
258:   _bbs = PetscBLASIntCast(bs); _bncols = PetscBLASIntCast(ncols); \
259:   BLASgemv_("T",&_bbs,&_bncols,&_one,A,&_bbs,x,&_ione,&_one,z,&_ione); \
260: }

262: #else 
263: /*
264:        Version that calls Fortran routines; can handle different precision
265:    of matrix (array) and vectors
266: */
267: /*
268:      These are Fortran kernels: They replace certain BLAS routines but
269:    have some arguments that may be single precision,rather than double
270:    These routines are provided in src/fortran/kernels/sgemv.F 
271:    They are pretty pitiful but get the job done. The intention is 
272:    that for important block sizes (currently 1,2,3,4,5,6,7) custom inlined 
273:    code is used.
274: */

276: /* BGL kernels */
277: #if defined(PETSC_USE_FORTRAN_KERNELS_BGL)
278: #define msgemv  msgemv_bgl
279: #define msgemvp msgemvp_bgl
280: #define msgemvm msgemvm_bgl
281: #define msgemvt msgemvt_bgl
282: #define msgemmi msgemmi_bgl
283: #define msgemm  msgemm_bgl
284: #endif

286: #ifdef PETSC_HAVE_FORTRAN_CAPS
287: #define msgemv_  MSGEMV
288: #define msgemvp_ MSGEMVP
289: #define msgemvm_ MSGEMVM
290: #define msgemvt_ MSGEMVT
291: #define msgemmi_ MSGEMMI
292: #define msgemm_  MSGEMM
293: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
294: #define msgemv_  msgemv
295: #define msgemvp_ msgemvp
296: #define msgemvm_ msgemvm
297: #define msgemvt_ msgemvt
298: #define msgemmi_ msgemmi
299: #define msgemm_  msgemm
300: #endif
302: EXTERN void msgemv_(PetscInt*,PetscInt *,MatScalar*,PetscScalar*,PetscScalar*);
303: EXTERN void msgemvp_(PetscInt*,PetscInt *,MatScalar*,PetscScalar*,PetscScalar*);
304: EXTERN void msgemvm_(PetscInt*,PetscInt *,MatScalar*,PetscScalar*,PetscScalar*);
305: EXTERN void msgemvt_(PetscInt*,PetscInt *,MatScalar*,PetscScalar*,PetscScalar*);
306: EXTERN void msgemmi_(PetscInt*,MatScalar*,MatScalar*,MatScalar*);
307: EXTERN void msgemm_(PetscInt*,MatScalar*,MatScalar*,MatScalar*);

310: /*
311:       A = A * B   A_gets_A_times_B

313:    A, B - square bs by bs arrays stored in column major order
314:    W    - square bs by bs work array

316: */
317: #define Kernel_A_gets_A_times_B(bs,A,B,W) \
318: { \
319:   PetscErrorCode _PetscMemcpy((W),(A),(bs)*(bs)*sizeof(MatScalar));CHKERRQ(_ierr); \
320:   msgemmi_(&bs,A,B,W); \
321: }

323: /*

325:     A = A - B * C  A_gets_A_minus_B_times_C 

327:    A, B, C - square bs by bs arrays stored in column major order
328: */
329: #define Kernel_A_gets_A_minus_B_times_C(bs,A,B,C) \
330: { \
331:   msgemm_(&bs,A,B,C); \
332: }

334: /*
335:     v = v - A w  v_gets_v_minus_A_times_w

337:    v - array of length bs
338:    A - square bs by bs array
339:    w - array of length bs
340: */
341: #define  Kernel_v_gets_v_minus_A_times_w(bs,v,A,w) \
342: {  \
343:   msgemvm_(&bs,&bs,A,w,v); \
344: }

346: /*
347:     v = v + A w  v_gets_v_plus_A_times_w

349:    v - array of length bs
350:    A - square bs by bs array
351:    w - array of length bs
352: */
353: #define  Kernel_v_gets_v_plus_A_times_w(bs,v,A,w) \
354: {  \
355:   msgemvp_(&bs,&bs,A,w,v);\
356: }

358: /*
359:     v = v + A w  v_gets_v_plus_Ar_times_w

361:    v - array of length bs
362:    A - square bs by bs array
363:    w - array of length bs
364: */
365: #define  Kernel_w_gets_w_plus_Ar_times_v(bs,ncol,v,A,w) \
366: {  \
367:   msgemvp_(&bs,&ncol,A,v,w);\
368: }

370: /*
371:     w = A v   w_gets_A_times_v

373:    v - array of length bs
374:    A - square bs by bs array
375:    w - array of length bs
376: */
377: #define Kernel_w_gets_A_times_v(bs,v,A,w) \
378: {  \
379:   msgemv_(&bs,&bs,A,v,w);\
380: }
381: 
382: /*
383:         z = A*x
384: */
385: #define Kernel_w_gets_Ar_times_v(bs,ncols,x,A,z) \
386: { \
387:   msgemv_(&bs,&ncols,A,x,z);\
388: }

390: /*
391:         z = trans(A)*x
392: */
393: #define Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,A,z) \
394: { \
395:   msgemvt_(&bs,&ncols,A,x,z);\
396: }

398: /* These do not work yet */
399: #define Kernel_A_gets_A_plus_Btranspose_times_C(bs,A,B,C) 
400: #define Kernel_v_gets_v_plus_Atranspose_times_w(bs,v,A,w) 


403: #endif

405: #endif