1 2 /* 3 Defines projective product routines where A is a SeqAIJ matrix 4 C = R * A * R^T 5 */ 6 7 #include <../src/mat/impls/aij/seq/aij.h> 8 #include <../src/mat/utils/freespace.h> 9 #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 10 #include <petsctime.h> 11 12 /* #define RART_PROFILE */ 13 14 #undef __FUNCT__ 15 #define __FUNCT__ "MatDestroy_SeqAIJ_RARt" 16 PetscErrorCode MatDestroy_SeqAIJ_RARt(Mat A) 17 { 18 PetscErrorCode ierr; 19 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 20 Mat_RARt *rart = a->rart; 21 22 PetscFunctionBegin; 23 ierr = MatTransposeColoringDestroy(&rart->matcoloring);CHKERRQ(ierr); 24 ierr = MatDestroy(&rart->Rt);CHKERRQ(ierr); 25 ierr = MatDestroy(&rart->RARt);CHKERRQ(ierr); 26 ierr = MatDestroy(&rart->ARt);CHKERRQ(ierr); 27 ierr = PetscFree(rart->work);CHKERRQ(ierr); 28 29 A->ops->destroy = rart->destroy; 30 if (A->ops->destroy) { 31 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 32 } 33 ierr = PetscFree(rart);CHKERRQ(ierr); 34 PetscFunctionReturn(0); 35 } 36 37 #undef __FUNCT__ 38 #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ" 39 PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C) 40 { 41 PetscErrorCode ierr; 42 Mat P; 43 PetscInt *rti,*rtj; 44 Mat_RARt *rart; 45 MatTransposeColoring matcoloring; 46 ISColoring iscoloring; 47 Mat Rt_dense,RARt_dense; 48 #if defined(RART_PROFILE) 49 PetscLogDouble GColor=0.0,MCCreate=0.0,MDenCreate=0.0,t0,tf,etime=0.0; 50 #endif 51 Mat_SeqAIJ *c; 52 53 PetscFunctionBegin; 54 #if defined(RART_PROFILE) 55 ierr = PetscTime(&t0);CHKERRQ(ierr); 56 #endif 57 /* create symbolic P=Rt */ 58 ierr = MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 59 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,NULL,&P);CHKERRQ(ierr); 60 61 /* get symbolic C=Pt*A*P */ 62 ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 63 (*C)->rmap->bs = R->rmap->bs; 64 (*C)->cmap->bs = R->rmap->bs; 65 66 /* create a supporting struct */ 67 ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 68 c = (Mat_SeqAIJ*)(*C)->data; 69 c->rart = rart; 70 #if defined(RART_PROFILE) 71 ierr = PetscTime(&tf);CHKERRQ(ierr); 72 etime += tf - t0; 73 #endif 74 75 /* ------ Use coloring ---------- */ 76 /* Create MatTransposeColoring from symbolic C=R*A*R^T */ 77 #if defined(RART_PROFILE) 78 ierr = PetscTime(&t0);CHKERRQ(ierr); 79 #endif 80 ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 81 #if defined(RART_PROFILE) 82 ierr = PetscTime(&tf);CHKERRQ(ierr); 83 GColor += tf - t0; 84 ierr = PetscTime(&t0);CHKERRQ(ierr); 85 #endif 86 ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 87 88 rart->matcoloring = matcoloring; 89 90 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 91 #if defined(RART_PROFILE) 92 ierr = PetscTime(&tf);CHKERRQ(ierr); 93 MCCreate += tf - t0; 94 ierr = PetscTime(&t0);CHKERRQ(ierr); 95 #endif 96 97 /* Create Rt_dense */ 98 ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr); 99 ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 100 ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr); 101 ierr = MatSeqDenseSetPreallocation(Rt_dense,NULL);CHKERRQ(ierr); 102 103 Rt_dense->assembled = PETSC_TRUE; 104 rart->Rt = Rt_dense; 105 106 /* Create RARt_dense = R*A*Rt_dense */ 107 ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr); 108 ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 109 ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr); 110 ierr = MatSeqDenseSetPreallocation(RARt_dense,NULL);CHKERRQ(ierr); 111 112 rart->RARt = RARt_dense; 113 114 /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 115 ierr = PetscMalloc(A->rmap->n*4*sizeof(PetscScalar),&rart->work);CHKERRQ(ierr); 116 117 #if defined(RART_PROFILE) 118 ierr = PetscTime(&tf);CHKERRQ(ierr); 119 MDenCreate += tf - t0; 120 #endif 121 122 rart->destroy = (*C)->ops->destroy; 123 (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 124 125 /* clean up */ 126 ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 127 ierr = MatDestroy(&P);CHKERRQ(ierr); 128 129 #if defined(PETSC_USE_INFO) 130 { 131 PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n); 132 ierr = PetscInfo6(*C,"RARt_den %D %D; Rt_den %D %D, (RARt->nz %D)/(m*ncolors)=%g\n",RARt_dense->rmap->n,RARt_dense->cmap->n,Rt_dense->rmap->n,Rt_dense->cmap->n,c->nz,density);CHKERRQ(ierr); 133 #if defined(RART_PROFILE) 134 ierr = PetscInfo5(*C,"Sym = GetColor %g + MColorCreate %g + MDenCreate %g + other %g = %g\n",GColor,MCCreate,MDenCreate,etime,GColor+MCCreate+MDenCreate+etime);CHKERRQ(ierr); 135 #endif 136 } 137 #endif 138 PetscFunctionReturn(0); 139 } 140 141 /* 142 RAB = R * A * B, R and A in seqaij format, B in dense format; 143 */ 144 #undef __FUNCT__ 145 #define __FUNCT__ "MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense" 146 PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work) 147 { 148 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data; 149 PetscErrorCode ierr; 150 PetscScalar *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 151 MatScalar *aa,*ra; 152 PetscInt cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n; 153 PetscInt am2=2*am,am3=3*am,bm4=4*bm; 154 PetscScalar *d,*c,*c2,*c3,*c4; 155 PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n; 156 PetscInt rm2=2*rm,rm3=3*rm,colrm; 157 #if defined(RART_PROFILE) 158 PetscLogDouble t0,tf,Mult_sp_den1=0.0,Mult_sp_den2=0.0; 159 #endif 160 161 PetscFunctionBegin; 162 if (!dm || !dn) PetscFunctionReturn(0); 163 if (bm != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in A %D not equal rows in B %D\n",A->cmap->n,bm); 164 if (am != R->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in R %D not equal rows in A %D\n",R->cmap->n,am); 165 if (R->rmap->n != RAB->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows in RAB %D not equal rows in R %D\n",RAB->rmap->n,R->rmap->n); 166 if (B->cmap->n != RAB->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in RAB %D not equal columns in B %D\n",RAB->cmap->n,B->cmap->n); 167 168 { /* 169 This approach is not as good as original ones (will be removed later), but it reveals that 170 AB_den=A*B takes almost all execution time in R*A*B for src/ksp/ksp/examples/tutorials/ex56.c 171 */ 172 PetscBool via_matmatmult=PETSC_FALSE; 173 ierr = PetscOptionsGetBool(NULL,"-matrart_via_matmatmult",&via_matmatmult,NULL);CHKERRQ(ierr); 174 if (via_matmatmult) { 175 Mat AB_den; 176 ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,0.0,&AB_den);CHKERRQ(ierr); 177 178 #if defined(RART_PROFILE) 179 ierr = PetscTime(&t0);CHKERRQ(ierr); 180 #endif 181 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,AB_den);CHKERRQ(ierr); 182 #if defined(RART_PROFILE) 183 ierr = PetscTime(&tf);CHKERRQ(ierr); 184 Mult_sp_den1 += tf - t0; 185 ierr = PetscTime(&t0);CHKERRQ(ierr); 186 #endif 187 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(R,AB_den,RAB);CHKERRQ(ierr); 188 #if defined(RART_PROFILE) 189 ierr = PetscTime(&tf);CHKERRQ(ierr); 190 Mult_sp_den2 += tf - t0; 191 printf(" MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense...time = %g + %g = %g; ncolor: %d\n",Mult_sp_den1, Mult_sp_den2, Mult_sp_den1+ Mult_sp_den2,B->cmap->n); 192 #endif 193 ierr = MatDestroy(&AB_den);CHKERRQ(ierr); 194 PetscFunctionReturn(0); 195 } 196 } 197 198 ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 199 ierr = MatDenseGetArray(RAB,&d);CHKERRQ(ierr); 200 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 201 c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am; 202 for (col=0; col<cn-4; col += 4) { /* over columns of C */ 203 for (i=0; i<am; i++) { /* over rows of A in those columns */ 204 r1 = r2 = r3 = r4 = 0.0; 205 n = ai[i+1] - ai[i]; 206 aj = a->j + ai[i]; 207 aa = a->a + ai[i]; 208 for (j=0; j<n; j++) { 209 r1 += (*aa)*b1[*aj]; 210 r2 += (*aa)*b2[*aj]; 211 r3 += (*aa)*b3[*aj]; 212 r4 += (*aa++)*b4[*aj++]; 213 } 214 c[i] = r1; 215 c[am + i] = r2; 216 c[am2 + i] = r3; 217 c[am3 + i] = r4; 218 } 219 b1 += bm4; 220 b2 += bm4; 221 b3 += bm4; 222 b4 += bm4; 223 224 /* RAB[:,col] = R*C[:,col] */ 225 colrm = col*rm; 226 for (i=0; i<rm; i++) { /* over rows of R in those columns */ 227 r1 = r2 = r3 = r4 = 0.0; 228 n = r->i[i+1] - r->i[i]; 229 rj = r->j + r->i[i]; 230 ra = r->a + r->i[i]; 231 for (j=0; j<n; j++) { 232 r1 += (*ra)*c[*rj]; 233 r2 += (*ra)*c2[*rj]; 234 r3 += (*ra)*c3[*rj]; 235 r4 += (*ra++)*c4[*rj++]; 236 } 237 d[colrm + i] = r1; 238 d[colrm + rm + i] = r2; 239 d[colrm + rm2 + i] = r3; 240 d[colrm + rm3 + i] = r4; 241 } 242 } 243 for (; col<cn; col++) { /* over extra columns of C */ 244 for (i=0; i<am; i++) { /* over rows of A in those columns */ 245 r1 = 0.0; 246 n = a->i[i+1] - a->i[i]; 247 aj = a->j + a->i[i]; 248 aa = a->a + a->i[i]; 249 for (j=0; j<n; j++) { 250 r1 += (*aa++)*b1[*aj++]; 251 } 252 c[i] = r1; 253 } 254 b1 += bm; 255 256 for (i=0; i<rm; i++) { /* over rows of R in those columns */ 257 r1 = 0.0; 258 n = r->i[i+1] - r->i[i]; 259 rj = r->j + r->i[i]; 260 ra = r->a + r->i[i]; 261 for (j=0; j<n; j++) { 262 r1 += (*ra++)*c[*rj++]; 263 } 264 d[col*rm + i] = r1; 265 } 266 } 267 ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr); 268 269 ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 270 ierr = MatDenseRestoreArray(RAB,&d);CHKERRQ(ierr); 271 ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 272 ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 273 PetscFunctionReturn(0); 274 } 275 276 #undef __FUNCT__ 277 #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ" 278 PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C) 279 { 280 PetscErrorCode ierr; 281 Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 282 Mat_RARt *rart=c->rart; 283 MatTransposeColoring matcoloring; 284 Mat Rt,RARt; 285 #if defined(RART_PROFILE) 286 PetscLogDouble Mult_sp_den=0.0,app1=0.0,app2=0.0,t0,tf; 287 #endif 288 289 PetscFunctionBegin; 290 /* Get dense Rt by Apply MatTransposeColoring to R */ 291 matcoloring = rart->matcoloring; 292 Rt = rart->Rt; 293 294 #if defined(RART_PROFILE) 295 ierr = PetscTime(&t0);CHKERRQ(ierr); 296 #endif 297 ierr = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr); 298 #if defined(RART_PROFILE) 299 ierr = PetscTime(&tf);CHKERRQ(ierr); 300 app1 += tf - t0; 301 #endif 302 303 /* Get dense RARt = R*A*Rt -- dominates! */ 304 #if defined(RART_PROFILE) 305 ierr = PetscTime(&t0);CHKERRQ(ierr); 306 #endif 307 RARt = rart->RARt; 308 ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr); 309 #if defined(RART_PROFILE) 310 ierr = PetscTime(&tf);CHKERRQ(ierr); 311 Mult_sp_den += tf - t0; 312 #endif 313 314 /* Recover C from C_dense */ 315 #if defined(RART_PROFILE) 316 ierr = PetscTime(&t0);CHKERRQ(ierr); 317 #endif 318 ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr); 319 #if defined(RART_PROFILE) 320 ierr = PetscTime(&tf);CHKERRQ(ierr); 321 app2 += tf - t0; 322 #endif 323 324 #if defined(RART_PROFILE) 325 printf("MatRARtNumeric_SeqAIJ_SeqAIJ: ColorApp %g + %g + Mult_sp_den %g = %g\n",app1,app2,Mult_sp_den,app1+app2+Mult_sp_den);CHKERRQ(ierr); 326 #endif 327 PetscFunctionReturn(0); 328 } 329 330 #undef __FUNCT__ 331 #define __FUNCT__ "MatRARt_SeqAIJ_SeqAIJ" 332 PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C) 333 { 334 PetscErrorCode ierr; 335 PetscBool usecoloring = PETSC_FALSE,color_art; 336 337 PetscFunctionBegin; 338 ierr = PetscOptionsGetBool(NULL,"-matrart_color",&usecoloring,NULL);CHKERRQ(ierr); 339 340 if (!usecoloring) { /* Rt=R^T, C=R*A*Rt - avoid inefficient sparse inner products */ 341 Mat Rt; 342 Mat_SeqAIJ *c; 343 Mat_RARt *rart; 344 345 if (scall == MAT_INITIAL_MATRIX) { 346 ierr = MatTranspose(R,MAT_INITIAL_MATRIX,&Rt);CHKERRQ(ierr); 347 ierr = MatMatMatMult(R,A,Rt,scall,fill,C);CHKERRQ(ierr); 348 349 ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 350 rart->Rt = Rt; 351 c = (Mat_SeqAIJ*)(*C)->data; 352 c->rart = rart; 353 rart->destroy = (*C)->ops->destroy; 354 (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 355 } else { 356 c = (Mat_SeqAIJ*)(*C)->data; 357 rart = c->rart; 358 Rt = rart->Rt; 359 ierr = MatTranspose(R,MAT_REUSE_MATRIX,&Rt);CHKERRQ(ierr); 360 ierr = MatMatMatMult(R,A,Rt,scall,fill,C);CHKERRQ(ierr); 361 } 362 PetscFunctionReturn(0); 363 } 364 365 /* use coloring */ 366 /*--------------*/ 367 ierr = PetscOptionsGetBool(NULL,"-matrart_color_art",&color_art,NULL);CHKERRQ(ierr); 368 if (color_art) { /* apply coloring to A*R^T, not C = R*A*R^T */ 369 Mat ARt,RARt; 370 #if defined(RART_PROFILE) 371 PetscLogDouble t0,t1,t2; 372 #endif 373 Mat_SeqAIJ *c; 374 Mat_RARt *rart; 375 376 if (scall == MAT_INITIAL_MATRIX) { 377 ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 378 /* must use '-mat_no_inode' with '-matmattransmult_color 1' - do not knwo why? */ 379 ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,R,fill,&ARt);CHKERRQ(ierr); 380 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,&RARt);CHKERRQ(ierr); 381 382 *C = RARt; 383 ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 384 c = (Mat_SeqAIJ*)(*C)->data; 385 c->rart = rart; 386 rart->ARt = ARt; 387 rart->destroy = RARt->ops->destroy; 388 RARt->ops->destroy = MatDestroy_SeqAIJ_RARt; 389 ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 390 } 391 ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 392 c = (Mat_SeqAIJ*)(*C)->data; 393 rart = c->rart; 394 ARt = rart->ARt; 395 #if defined(RART_PROFILE) 396 ierr = PetscTime(&t0);CHKERRQ(ierr); 397 #endif 398 ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,ARt);CHKERRQ(ierr); /* dominate! */ 399 #if defined(RART_PROFILE) 400 ierr = PetscTime(&t1);CHKERRQ(ierr); 401 #endif 402 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(R,ARt,*C);CHKERRQ(ierr); 403 #if defined(RART_PROFILE) 404 ierr = PetscTime(&t2);CHKERRQ(ierr); 405 printf(" matrart_color_art_num = %g + %g = %g\n",t1-t0,t2-t1,t2-t0); 406 #endif 407 ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 408 PetscFunctionReturn(0); 409 } 410 411 if (scall == MAT_INITIAL_MATRIX) { /* apply coloring to C=R*A*R^T */ 412 ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 413 ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr); 414 ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 415 } 416 ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 417 ierr = MatRARtNumeric_SeqAIJ_SeqAIJ(A,R,*C);CHKERRQ(ierr); 418 ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 419 PetscFunctionReturn(0); 420 } 421