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