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