1 2 #include <../src/mat/impls/aij/seq/aij.h> 3 #include <../src/mat/utils/freespace.h> 4 #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 5 #include <petsctime.h> 6 7 /* 8 Defines projective product routines where A is a SeqAIJ matrix 9 C = R * A * R^T 10 */ 11 12 #undef __FUNCT__ 13 #define __FUNCT__ "PetscContainerDestroy_Mat_RARt" 14 PetscErrorCode PetscContainerDestroy_Mat_RARt(void *ptr) 15 { 16 PetscErrorCode ierr; 17 Mat_RARt *rart=(Mat_RARt*)ptr; 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 = PetscFree(rart->work);CHKERRQ(ierr); 24 ierr = PetscFree(rart);CHKERRQ(ierr); 25 PetscFunctionReturn(0); 26 } 27 28 #undef __FUNCT__ 29 #define __FUNCT__ "MatDestroy_SeqAIJ_RARt" 30 PetscErrorCode MatDestroy_SeqAIJ_RARt(Mat A) 31 { 32 PetscErrorCode ierr; 33 PetscContainer container; 34 Mat_RARt *rart=NULL; 35 36 PetscFunctionBegin; 37 ierr = PetscObjectQuery((PetscObject)A,"Mat_RARt",(PetscObject*)&container);CHKERRQ(ierr); 38 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 39 ierr = PetscContainerGetPointer(container,(void**)&rart);CHKERRQ(ierr); 40 A->ops->destroy = rart->destroy; 41 if (A->ops->destroy) { 42 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 43 } 44 ierr = PetscObjectCompose((PetscObject)A,"Mat_RARt",0);CHKERRQ(ierr); 45 PetscFunctionReturn(0); 46 } 47 48 #undef __FUNCT__ 49 #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ" 50 PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C) 51 { 52 PetscErrorCode ierr; 53 Mat P; 54 PetscInt *rti,*rtj; 55 Mat_RARt *rart; 56 PetscContainer container; 57 MatTransposeColoring matcoloring; 58 ISColoring iscoloring; 59 Mat Rt_dense,RARt_dense; 60 PetscLogDouble GColor=0.0,MCCreate=0.0,MDenCreate=0.0,t0,tf,etime=0.0; 61 Mat_SeqAIJ *c; 62 63 PetscFunctionBegin; 64 ierr = PetscTime(&t0);CHKERRQ(ierr); 65 /* create symbolic P=Rt */ 66 ierr = MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 67 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,NULL,&P);CHKERRQ(ierr); 68 69 /* get symbolic C=Pt*A*P */ 70 ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 71 (*C)->rmap->bs = R->rmap->bs; 72 (*C)->cmap->bs = R->rmap->bs; 73 74 /* create a supporting struct */ 75 ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 76 77 /* attach the supporting struct to C */ 78 ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 79 ierr = PetscContainerSetPointer(container,rart);CHKERRQ(ierr); 80 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_RARt);CHKERRQ(ierr); 81 ierr = PetscObjectCompose((PetscObject)(*C),"Mat_RARt",(PetscObject)container);CHKERRQ(ierr); 82 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 83 ierr = PetscTime(&tf);CHKERRQ(ierr); 84 etime += tf - t0; 85 86 /* ------ Use coloring ---------- */ 87 /* Create MatTransposeColoring from symbolic C=R*A*R^T */ 88 c = (Mat_SeqAIJ*)(*C)->data; 89 ierr = PetscTime(&t0);CHKERRQ(ierr); 90 ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 91 ierr = PetscTime(&tf);CHKERRQ(ierr); 92 GColor += tf - t0; 93 94 ierr = PetscTime(&t0);CHKERRQ(ierr); 95 ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 96 97 rart->matcoloring = matcoloring; 98 99 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 100 ierr = PetscTime(&tf);CHKERRQ(ierr); 101 MCCreate += tf - t0; 102 103 ierr = PetscTime(&t0);CHKERRQ(ierr); 104 /* Create Rt_dense */ 105 ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr); 106 ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 107 ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr); 108 ierr = MatSeqDenseSetPreallocation(Rt_dense,NULL);CHKERRQ(ierr); 109 110 Rt_dense->assembled = PETSC_TRUE; 111 rart->Rt = Rt_dense; 112 113 /* Create RARt_dense = R*A*Rt_dense */ 114 ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr); 115 ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 116 ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr); 117 ierr = MatSeqDenseSetPreallocation(RARt_dense,NULL);CHKERRQ(ierr); 118 119 rart->RARt = RARt_dense; 120 121 /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 122 ierr = PetscMalloc(A->rmap->n*4*sizeof(PetscScalar),&rart->work);CHKERRQ(ierr); 123 124 ierr = PetscTime(&tf);CHKERRQ(ierr); 125 MDenCreate += tf - t0; 126 127 rart->destroy = (*C)->ops->destroy; 128 (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 129 130 /* clean up */ 131 ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 132 ierr = MatDestroy(&P);CHKERRQ(ierr); 133 134 #if defined(PETSC_USE_INFO) 135 { 136 PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n); 137 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); 138 ierr = PetscInfo5(*C,"Sym = GetColor %g + MColorCreate %g + MDenCreate %g + other %g = %g\n",GColor,MCCreate,MDenCreate,etime,GColor+MCCreate+MDenCreate+etime);CHKERRQ(ierr); 139 } 140 #endif 141 PetscFunctionReturn(0); 142 } 143 144 /* 145 RAB = R * A * B, R and A in seqaij format, B in dense format; 146 */ 147 #undef __FUNCT__ 148 #define __FUNCT__ "MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense" 149 PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work) 150 { 151 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data; 152 PetscErrorCode ierr; 153 PetscScalar *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 154 MatScalar *aa,*ra; 155 PetscInt cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n; 156 PetscInt am2=2*am,am3=3*am,bm4=4*bm; 157 PetscScalar *d,*c,*c2,*c3,*c4; 158 PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n; 159 PetscInt rm2=2*rm,rm3=3*rm,colrm; 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 ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 169 ierr = MatDenseGetArray(RAB,&d);CHKERRQ(ierr); 170 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 171 c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am; 172 for (col=0; col<cn-4; col += 4) { /* over columns of C */ 173 for (i=0; i<am; i++) { /* over rows of A in those columns */ 174 r1 = r2 = r3 = r4 = 0.0; 175 n = ai[i+1] - ai[i]; 176 aj = a->j + ai[i]; 177 aa = a->a + ai[i]; 178 for (j=0; j<n; j++) { 179 r1 += (*aa)*b1[*aj]; 180 r2 += (*aa)*b2[*aj]; 181 r3 += (*aa)*b3[*aj]; 182 r4 += (*aa++)*b4[*aj++]; 183 } 184 c[i] = r1; 185 c[am + i] = r2; 186 c[am2 + i] = r3; 187 c[am3 + i] = r4; 188 } 189 b1 += bm4; 190 b2 += bm4; 191 b3 += bm4; 192 b4 += bm4; 193 194 /* RAB[:,col] = R*C[:,col] */ 195 colrm = col*rm; 196 for (i=0; i<rm; i++) { /* over rows of R in those columns */ 197 r1 = r2 = r3 = r4 = 0.0; 198 n = r->i[i+1] - r->i[i]; 199 rj = r->j + r->i[i]; 200 ra = r->a + r->i[i]; 201 for (j=0; j<n; j++) { 202 r1 += (*ra)*c[*rj]; 203 r2 += (*ra)*c2[*rj]; 204 r3 += (*ra)*c3[*rj]; 205 r4 += (*ra++)*c4[*rj++]; 206 } 207 d[colrm + i] = r1; 208 d[colrm + rm + i] = r2; 209 d[colrm + rm2 + i] = r3; 210 d[colrm + rm3 + i] = r4; 211 } 212 } 213 for (; col<cn; col++) { /* over extra columns of C */ 214 for (i=0; i<am; i++) { /* over rows of A in those columns */ 215 r1 = 0.0; 216 n = a->i[i+1] - a->i[i]; 217 aj = a->j + a->i[i]; 218 aa = a->a + a->i[i]; 219 for (j=0; j<n; j++) { 220 r1 += (*aa++)*b1[*aj++]; 221 } 222 c[i] = r1; 223 } 224 b1 += bm; 225 226 for (i=0; i<rm; i++) { /* over rows of R in those columns */ 227 r1 = 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 } 234 d[col*rm + i] = r1; 235 } 236 } 237 ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr); 238 239 ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 240 ierr = MatDenseRestoreArray(RAB,&d);CHKERRQ(ierr); 241 ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 242 ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 243 PetscFunctionReturn(0); 244 } 245 246 #undef __FUNCT__ 247 #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ" 248 PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C) 249 { 250 PetscErrorCode ierr; 251 Mat_RARt *rart; 252 PetscContainer container; 253 MatTransposeColoring matcoloring; 254 Mat Rt,RARt; 255 PetscLogDouble Mult_sp_den=0.0,app1=0.0,app2=0.0,t0,tf; 256 257 PetscFunctionBegin; 258 ierr = PetscObjectQuery((PetscObject)C,"Mat_RARt",(PetscObject*)&container);CHKERRQ(ierr); 259 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 260 ierr = PetscContainerGetPointer(container,(void**)&rart);CHKERRQ(ierr); 261 262 /* Get dense Rt by Apply MatTransposeColoring to R */ 263 matcoloring = rart->matcoloring; 264 Rt = rart->Rt; 265 266 ierr = PetscTime(&t0);CHKERRQ(ierr); 267 ierr = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr); 268 ierr = PetscTime(&tf);CHKERRQ(ierr); 269 app1 += tf - t0; 270 271 /* Get dense RARt = R*A*Rt -- dominates! */ 272 ierr = PetscTime(&t0);CHKERRQ(ierr); 273 RARt = rart->RARt; 274 ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr); 275 ierr = PetscTime(&tf);CHKERRQ(ierr); 276 277 Mult_sp_den += tf - t0; 278 279 /* Recover C from C_dense */ 280 ierr = PetscTime(&t0);CHKERRQ(ierr); 281 ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr); 282 ierr = PetscTime(&tf);CHKERRQ(ierr); 283 284 app2 += tf - t0; 285 286 #if defined(PETSC_USE_INFO) 287 ierr = PetscInfo4(C,"Num = ColorApp %g + %g + Mult_sp_den %g = %g\n",app1,app2,Mult_sp_den,app1+app2+Mult_sp_den);CHKERRQ(ierr); 288 #endif 289 PetscFunctionReturn(0); 290 } 291 292 #undef __FUNCT__ 293 #define __FUNCT__ "MatRARt_SeqAIJ_SeqAIJ" 294 PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C) 295 { 296 PetscErrorCode ierr; 297 PetscBool usecoloring = PETSC_TRUE; 298 299 PetscFunctionBegin; 300 ierr = PetscOptionsGetBool(NULL,"-matrart_color",&usecoloring,NULL);CHKERRQ(ierr); 301 302 if (!usecoloring) { 303 Mat Rt; 304 ierr = MatTranspose(R,MAT_INITIAL_MATRIX,&Rt);CHKERRQ(ierr); /* replace MAT_INITIAL_MATRIX with scall if !usecoloring is better */ 305 ierr = MatMatMatMult(R,A,Rt,scall,fill,C);CHKERRQ(ierr); 306 ierr = MatDestroy(&Rt);CHKERRQ(ierr); 307 PetscFunctionReturn(0); 308 } 309 310 /* use coloring */ 311 PetscBool newimpl=PETSC_FALSE; 312 ierr = PetscOptionsGetBool(NULL,"-matrart_new",&newimpl,NULL);CHKERRQ(ierr); 313 if (newimpl) { 314 Mat ARt,RARt; 315 printf(" MatRARt_new ....\n"); 316 if (scall == MAT_INITIAL_MATRIX) { 317 ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 318 ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,R,fill,&ARt);CHKERRQ(ierr); 319 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,&RARt);CHKERRQ(ierr); 320 ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 321 } 322 ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 323 ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,ARt);CHKERRQ(ierr); 324 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(R,ARt,RARt);CHKERRQ(ierr); 325 *C = RARt; 326 ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 327 328 ierr = MatDestroy(&ARt);CHKERRQ(ierr); 329 PetscFunctionReturn(0); 330 } 331 if (scall == MAT_INITIAL_MATRIX) { 332 ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 333 ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr); 334 ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 335 } 336 ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 337 ierr = MatRARtNumeric_SeqAIJ_SeqAIJ(A,R,*C);CHKERRQ(ierr); 338 ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 339 PetscFunctionReturn(0); 340 } 341