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