1 2 /* 3 Defines matrix-matrix product routines for pairs of SeqAIJ matrices 4 C = P * A * P^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 /* 12 MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ - Forms the symbolic product of two SeqAIJ matrices 13 C = P * A * P^T; 14 15 Note: C is assumed to be uncreated. 16 If this is not the case, Destroy C before calling this routine. 17 */ 18 #undef __FUNCT__ 19 #define __FUNCT__ "MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ" 20 PetscErrorCode MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat *C) 21 { 22 /* Note: This code is virtually identical to that of MatApplyPtAP_SeqAIJ_Symbolic */ 23 /* and MatMatMult_SeqAIJ_SeqAIJ_Symbolic. Perhaps they could be merged nicely. */ 24 PetscErrorCode ierr; 25 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 26 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; 27 PetscInt *ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pti,*ptj,*ptjj; 28 PetscInt *ci,*cj,*paj,*padenserow,*pasparserow,*denserow,*sparserow; 29 PetscInt an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N; 30 PetscInt i,j,k,pnzi,arow,anzj,panzi,ptrow,ptnzj,cnzi; 31 MatScalar *ca; 32 33 PetscFunctionBegin; 34 /* some error checking which could be moved into interface layer */ 35 if (pn!=am) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pn,am); 36 if (am!=an) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",am, an); 37 38 /* Set up timers */ 39 ierr = PetscLogEventBegin(MAT_Applypapt_symbolic,A,P,0,0);CHKERRQ(ierr); 40 41 /* Create ij structure of P^T */ 42 ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 43 44 /* Allocate ci array, arrays for fill computation and */ 45 /* free space for accumulating nonzero column info */ 46 ierr = PetscMalloc(((pm+1)*1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 47 ci[0] = 0; 48 49 ierr = PetscMalloc4(an,PetscInt,&padenserow,an,PetscInt,&pasparserow,pm,PetscInt,&denserow,pm,PetscInt,&sparserow);CHKERRQ(ierr); 50 ierr = PetscMemzero(padenserow,an*sizeof(PetscInt));CHKERRQ(ierr); 51 ierr = PetscMemzero(pasparserow,an*sizeof(PetscInt));CHKERRQ(ierr); 52 ierr = PetscMemzero(denserow,pm*sizeof(PetscInt));CHKERRQ(ierr); 53 ierr = PetscMemzero(sparserow,pm*sizeof(PetscInt));CHKERRQ(ierr); 54 55 /* Set initial free space to be nnz(A) scaled by aspect ratio of Pt. */ 56 /* This should be reasonable if sparsity of PAPt is similar to that of A. */ 57 ierr = PetscFreeSpaceGet((ai[am]/pn)*pm,&free_space); 58 current_space = free_space; 59 60 /* Determine fill for each row of C: */ 61 for (i=0;i<pm;i++) { 62 pnzi = pi[i+1] - pi[i]; 63 panzi = 0; 64 /* Get symbolic sparse row of PA: */ 65 for (j=0;j<pnzi;j++) { 66 arow = *pj++; 67 anzj = ai[arow+1] - ai[arow]; 68 ajj = aj + ai[arow]; 69 for (k=0;k<anzj;k++) { 70 if (!padenserow[ajj[k]]) { 71 padenserow[ajj[k]] = -1; 72 pasparserow[panzi++] = ajj[k]; 73 } 74 } 75 } 76 /* Using symbolic row of PA, determine symbolic row of C: */ 77 paj = pasparserow; 78 cnzi = 0; 79 for (j=0;j<panzi;j++) { 80 ptrow = *paj++; 81 ptnzj = pti[ptrow+1] - pti[ptrow]; 82 ptjj = ptj + pti[ptrow]; 83 for (k=0;k<ptnzj;k++) { 84 if (!denserow[ptjj[k]]) { 85 denserow[ptjj[k]] = -1; 86 sparserow[cnzi++] = ptjj[k]; 87 } 88 } 89 } 90 91 /* sort sparse representation */ 92 ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr); 93 94 /* If free space is not available, make more free space */ 95 /* Double the amount of total space in the list */ 96 if (current_space->local_remaining<cnzi) { 97 ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 98 } 99 100 /* Copy data into free space, and zero out dense row */ 101 ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));CHKERRQ(ierr); 102 current_space->array += cnzi; 103 current_space->local_used += cnzi; 104 current_space->local_remaining -= cnzi; 105 106 for (j=0;j<panzi;j++) { 107 padenserow[pasparserow[j]] = 0; 108 } 109 for (j=0;j<cnzi;j++) { 110 denserow[sparserow[j]] = 0; 111 } 112 ci[i+1] = ci[i] + cnzi; 113 } 114 /* column indices are in the list of free space */ 115 /* Allocate space for cj, initialize cj, and */ 116 /* destroy list of free space and other temporary array(s) */ 117 ierr = PetscMalloc((ci[pm]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 118 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 119 ierr = PetscFree4(padenserow,pasparserow,denserow,sparserow);CHKERRQ(ierr); 120 121 /* Allocate space for ca */ 122 ierr = PetscMalloc((ci[pm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 123 ierr = PetscMemzero(ca,(ci[pm]+1)*sizeof(MatScalar));CHKERRQ(ierr); 124 125 /* put together the new matrix */ 126 ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pm,pm,ci,cj,ca,C);CHKERRQ(ierr); 127 (*C)->rmap->bs = P->cmap->bs; 128 (*C)->cmap->bs = P->cmap->bs; 129 130 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 131 /* Since these are PETSc arrays, change flags to free them as necessary. */ 132 c = (Mat_SeqAIJ *)((*C)->data); 133 c->free_a = PETSC_TRUE; 134 c->free_ij = PETSC_TRUE; 135 c->nonew = 0; 136 137 /* Clean up. */ 138 ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 139 140 ierr = PetscLogEventEnd(MAT_Applypapt_symbolic,A,P,0,0);CHKERRQ(ierr); 141 PetscFunctionReturn(0); 142 } 143 144 /* 145 MatApplyPAPt_Numeric_SeqAIJ - Forms the numeric product of two SeqAIJ matrices 146 C = P * A * P^T; 147 Note: C must have been created by calling MatApplyPAPt_Symbolic_SeqAIJ. 148 */ 149 #undef __FUNCT__ 150 #define __FUNCT__ "MatApplyPAPt_Numeric_SeqAIJ_SeqAIJ" 151 PetscErrorCode MatApplyPAPt_Numeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 152 { 153 PetscErrorCode ierr; 154 PetscInt flops=0; 155 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 156 Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 157 Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 158 PetscInt *ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj=p->j,*paj,*pajdense,*ptj; 159 PetscInt *ci=c->i,*cj=c->j; 160 PetscInt an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N,cn=C->cmap->N,cm=C->rmap->N; 161 PetscInt i,j,k,k1,k2,pnzi,anzj,panzj,arow,ptcol,ptnzj,cnzi; 162 MatScalar *aa=a->a,*pa=p->a,*pta=p->a,*ptaj,*paa,*aaj,*ca=c->a,sum; 163 164 PetscFunctionBegin; 165 /* This error checking should be unnecessary if the symbolic was performed */ 166 if (pm!=cm) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pm,cm); 167 if (pn!=am) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pn,am); 168 if (am!=an) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",am, an); 169 if (pm!=cn) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pm, cn); 170 171 /* Set up timers */ 172 ierr = PetscLogEventBegin(MAT_Applypapt_numeric,A,P,C,0);CHKERRQ(ierr); 173 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 174 175 ierr = PetscMalloc3(an,MatScalar,&paa,an,PetscInt,&paj,an,PetscInt,&pajdense);CHKERRQ(ierr); 176 ierr = PetscMemzero(paa,an*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr); 177 178 for (i=0;i<pm;i++) { 179 /* Form sparse row of P*A */ 180 pnzi = pi[i+1] - pi[i]; 181 panzj = 0; 182 for (j=0;j<pnzi;j++) { 183 arow = *pj++; 184 anzj = ai[arow+1] - ai[arow]; 185 ajj = aj + ai[arow]; 186 aaj = aa + ai[arow]; 187 for (k=0;k<anzj;k++) { 188 if (!pajdense[ajj[k]]) { 189 pajdense[ajj[k]] = -1; 190 paj[panzj++] = ajj[k]; 191 } 192 paa[ajj[k]] += (*pa)*aaj[k]; 193 } 194 flops += 2*anzj; 195 pa++; 196 } 197 198 /* Sort the j index array for quick sparse axpy. */ 199 ierr = PetscSortInt(panzj,paj);CHKERRQ(ierr); 200 201 /* Compute P*A*P^T using sparse inner products. */ 202 /* Take advantage of pre-computed (i,j) of C for locations of non-zeros. */ 203 cnzi = ci[i+1] - ci[i]; 204 for (j=0;j<cnzi;j++) { 205 /* Form sparse inner product of current row of P*A with (*cj++) col of P^T. */ 206 ptcol = *cj++; 207 ptnzj = pi[ptcol+1] - pi[ptcol]; 208 ptj = pjj + pi[ptcol]; 209 ptaj = pta + pi[ptcol]; 210 sum = 0.; 211 k1 = 0; 212 k2 = 0; 213 while ((k1<panzj) && (k2<ptnzj)) { 214 if (paj[k1]==ptj[k2]) { 215 sum += paa[paj[k1++]]*ptaj[k2++]; 216 } else if (paj[k1] < ptj[k2]) { 217 k1++; 218 } else /* if (paj[k1] > ptj[k2]) */ { 219 k2++; 220 } 221 } 222 *ca++ = sum; 223 } 224 225 /* Zero the current row info for P*A */ 226 for (j=0;j<panzj;j++) { 227 paa[paj[j]] = 0.; 228 pajdense[paj[j]] = 0; 229 } 230 } 231 232 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 233 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 234 ierr = PetscFree3(paa,paj,pajdense);CHKERRQ(ierr); 235 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 236 ierr = PetscLogEventEnd(MAT_Applypapt_numeric,A,P,C,0);CHKERRQ(ierr); 237 PetscFunctionReturn(0); 238 } 239 240 #undef __FUNCT__ 241 #define __FUNCT__ "MatApplyPAPt_SeqAIJ_SeqAIJ" 242 PetscErrorCode MatApplyPAPt_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat *C) 243 { 244 PetscErrorCode ierr; 245 246 PetscFunctionBegin; 247 ierr = PetscLogEventBegin(MAT_Applypapt,A,P,0,0);CHKERRQ(ierr); 248 ierr = MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ(A,P,C);CHKERRQ(ierr); 249 ierr = MatApplyPAPt_Numeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr); 250 ierr = PetscLogEventEnd(MAT_Applypapt,A,P,0,0);CHKERRQ(ierr); 251 PetscFunctionReturn(0); 252 } 253 254 /*--------------------------------------------------*/ 255 /* 256 Defines projective product routines where A is a SeqAIJ matrix 257 C = R * A * R^T 258 */ 259 260 #undef __FUNCT__ 261 #define __FUNCT__ "PetscContainerDestroy_Mat_RARt" 262 PetscErrorCode PetscContainerDestroy_Mat_RARt(void *ptr) 263 { 264 PetscErrorCode ierr; 265 Mat_RARt *rart=(Mat_RARt*)ptr; 266 267 PetscFunctionBegin; 268 ierr = MatTransposeColoringDestroy(&rart->matcoloring);CHKERRQ(ierr); 269 ierr = MatDestroy(&rart->Rt);CHKERRQ(ierr); 270 ierr = MatDestroy(&rart->RARt);CHKERRQ(ierr); 271 ierr = PetscFree(rart->work);CHKERRQ(ierr); 272 ierr = PetscFree(rart);CHKERRQ(ierr); 273 PetscFunctionReturn(0); 274 } 275 276 #undef __FUNCT__ 277 #define __FUNCT__ "MatDestroy_SeqAIJ_RARt" 278 PetscErrorCode MatDestroy_SeqAIJ_RARt(Mat A) 279 { 280 PetscErrorCode ierr; 281 PetscContainer container; 282 Mat_RARt *rart=PETSC_NULL; 283 284 PetscFunctionBegin; 285 ierr = PetscObjectQuery((PetscObject)A,"Mat_RARt",(PetscObject *)&container);CHKERRQ(ierr); 286 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 287 ierr = PetscContainerGetPointer(container,(void **)&rart);CHKERRQ(ierr); 288 A->ops->destroy = rart->destroy; 289 if (A->ops->destroy) { 290 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 291 } 292 ierr = PetscObjectCompose((PetscObject)A,"Mat_RARt",0);CHKERRQ(ierr); 293 PetscFunctionReturn(0); 294 } 295 296 #undef __FUNCT__ 297 #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ" 298 PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C) 299 { 300 PetscErrorCode ierr; 301 Mat P; 302 PetscInt *rti,*rtj; 303 Mat_RARt *rart; 304 PetscContainer container; 305 MatTransposeColoring matcoloring; 306 ISColoring iscoloring; 307 Mat Rt_dense,RARt_dense; 308 PetscLogDouble GColor=0.0,MCCreate=0.0,MDenCreate=0.0,t0,tf,etime=0.0; 309 Mat_SeqAIJ *c; 310 311 PetscFunctionBegin; 312 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 313 /* create symbolic P=Rt */ 314 ierr = MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 315 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,PETSC_NULL,&P);CHKERRQ(ierr); 316 317 /* get symbolic C=Pt*A*P */ 318 ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 319 (*C)->rmap->bs = R->rmap->bs; 320 (*C)->cmap->bs = R->rmap->bs; 321 322 /* create a supporting struct */ 323 ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 324 325 /* attach the supporting struct to C */ 326 ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 327 ierr = PetscContainerSetPointer(container,rart);CHKERRQ(ierr); 328 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_RARt);CHKERRQ(ierr); 329 ierr = PetscObjectCompose((PetscObject)(*C),"Mat_RARt",(PetscObject)container);CHKERRQ(ierr); 330 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 331 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 332 etime += tf - t0; 333 334 /* Create MatTransposeColoring from symbolic C=R*A*R^T */ 335 c=(Mat_SeqAIJ*)(*C)->data; 336 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 337 ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 338 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 339 GColor += tf - t0; 340 341 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 342 ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 343 rart->matcoloring = matcoloring; 344 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 345 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 346 MCCreate += tf - t0; 347 348 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 349 /* Create Rt_dense */ 350 ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr); 351 ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 352 ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr); 353 ierr = MatSeqDenseSetPreallocation(Rt_dense,PETSC_NULL);CHKERRQ(ierr); 354 Rt_dense->assembled = PETSC_TRUE; 355 rart->Rt = Rt_dense; 356 357 /* Create RARt_dense = R*A*Rt_dense */ 358 ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr); 359 ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 360 ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr); 361 ierr = MatSeqDenseSetPreallocation(RARt_dense,PETSC_NULL);CHKERRQ(ierr); 362 rart->RARt = RARt_dense; 363 364 /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 365 ierr = PetscMalloc(A->rmap->n*4*sizeof(PetscScalar),&rart->work);CHKERRQ(ierr); 366 367 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 368 MDenCreate += tf - t0; 369 370 rart->destroy = (*C)->ops->destroy; 371 (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 372 373 /* clean up */ 374 ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 375 ierr = MatDestroy(&P);CHKERRQ(ierr); 376 377 #if defined(PETSC_USE_INFO) 378 { 379 PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n); 380 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); 381 ierr = PetscInfo5(*C,"Sym = GetColor %g + MColorCreate %g + MDenCreate %g + other %g = %g\n",GColor,MCCreate,MDenCreate,etime,GColor+MCCreate+MDenCreate+etime); 382 } 383 #endif 384 PetscFunctionReturn(0); 385 } 386 387 /* 388 RAB = R * A * B, R and A in seqaij format, B in dense format; 389 */ 390 #undef __FUNCT__ 391 #define __FUNCT__ "MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense" 392 PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work) 393 { 394 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data; 395 PetscErrorCode ierr; 396 PetscScalar *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 397 MatScalar *aa,*ra; 398 PetscInt cn=B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n; 399 PetscInt am2=2*am,am3=3*am,bm4=4*bm; 400 PetscScalar *d,*c,*c2,*c3,*c4; 401 PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n; 402 PetscInt rm2=2*rm,rm3=3*rm,colrm; 403 404 PetscFunctionBegin; 405 if (!dm || !dn) PetscFunctionReturn(0); 406 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); 407 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); 408 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); 409 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); 410 411 ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 412 ierr = MatDenseGetArray(RAB,&d);CHKERRQ(ierr); 413 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 414 c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am; 415 for (col=0; col<cn-4; col += 4){ /* over columns of C */ 416 for (i=0; i<am; i++) { /* over rows of A in those columns */ 417 r1 = r2 = r3 = r4 = 0.0; 418 n = ai[i+1] - ai[i]; 419 aj = a->j + ai[i]; 420 aa = a->a + ai[i]; 421 for (j=0; j<n; j++) { 422 r1 += (*aa)*b1[*aj]; 423 r2 += (*aa)*b2[*aj]; 424 r3 += (*aa)*b3[*aj]; 425 r4 += (*aa++)*b4[*aj++]; 426 } 427 c[i] = r1; 428 c[am + i] = r2; 429 c[am2 + i] = r3; 430 c[am3 + i] = r4; 431 } 432 b1 += bm4; 433 b2 += bm4; 434 b3 += bm4; 435 b4 += bm4; 436 437 /* RAB[:,col] = R*C[:,col] */ 438 colrm = col*rm; 439 for (i=0; i<rm; i++) { /* over rows of R in those columns */ 440 r1 = r2 = r3 = r4 = 0.0; 441 n = r->i[i+1] - r->i[i]; 442 rj = r->j + r->i[i]; 443 ra = r->a + r->i[i]; 444 for (j=0; j<n; j++) { 445 r1 += (*ra)*c[*rj]; 446 r2 += (*ra)*c2[*rj]; 447 r3 += (*ra)*c3[*rj]; 448 r4 += (*ra++)*c4[*rj++]; 449 } 450 d[colrm + i] = r1; 451 d[colrm + rm + i] = r2; 452 d[colrm + rm2 + i] = r3; 453 d[colrm + rm3 + i] = r4; 454 } 455 } 456 for (;col<cn; col++){ /* over extra columns of C */ 457 for (i=0; i<am; i++) { /* over rows of A in those columns */ 458 r1 = 0.0; 459 n = a->i[i+1] - a->i[i]; 460 aj = a->j + a->i[i]; 461 aa = a->a + a->i[i]; 462 for (j=0; j<n; j++) { 463 r1 += (*aa++)*b1[*aj++]; 464 } 465 c[i] = r1; 466 } 467 b1 += bm; 468 469 for (i=0; i<rm; i++) { /* over rows of R in those columns */ 470 r1 = 0.0; 471 n = r->i[i+1] - r->i[i]; 472 rj = r->j + r->i[i]; 473 ra = r->a + r->i[i]; 474 for (j=0; j<n; j++) { 475 r1 += (*ra++)*c[*rj++]; 476 } 477 d[col*rm + i] = r1; 478 } 479 } 480 ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr); 481 482 ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 483 ierr = MatDenseRestoreArray(RAB,&d);CHKERRQ(ierr); 484 ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 485 ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 486 PetscFunctionReturn(0); 487 } 488 489 #undef __FUNCT__ 490 #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ" 491 PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C) 492 { 493 PetscErrorCode ierr; 494 Mat_RARt *rart; 495 PetscContainer container; 496 MatTransposeColoring matcoloring; 497 Mat Rt,RARt; 498 PetscLogDouble Mult_sp_den=0.0,app1=0.0,app2=0.0,t0,tf; 499 500 PetscFunctionBegin; 501 ierr = PetscObjectQuery((PetscObject)C,"Mat_RARt",(PetscObject *)&container);CHKERRQ(ierr); 502 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 503 ierr = PetscContainerGetPointer(container,(void **)&rart);CHKERRQ(ierr); 504 505 /* Get dense Rt by Apply MatTransposeColoring to R */ 506 matcoloring = rart->matcoloring; 507 Rt = rart->Rt; 508 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 509 ierr = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr); 510 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 511 app1 += tf - t0; 512 513 /* Get dense RARt = R*A*Rt */ 514 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 515 RARt = rart->RARt; 516 ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr); 517 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 518 Mult_sp_den += tf - t0; 519 520 /* Recover C from C_dense */ 521 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 522 ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr); 523 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 524 app2 += tf - t0; 525 526 #if defined(PETSC_USE_INFO) 527 ierr = PetscInfo4(C,"Num = ColorApp %g + %g + Mult_sp_den %g = %g\n",app1,app2,Mult_sp_den,app1+app2+Mult_sp_den); 528 #endif 529 PetscFunctionReturn(0); 530 } 531 532 EXTERN_C_BEGIN 533 #undef __FUNCT__ 534 #define __FUNCT__ "MatRARt_SeqAIJ_SeqAIJ" 535 PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C) 536 { 537 PetscErrorCode ierr; 538 539 PetscFunctionBegin; 540 if (scall == MAT_INITIAL_MATRIX){ 541 ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr); 542 } 543 ierr = MatRARtNumeric_SeqAIJ_SeqAIJ(A,R,*C);CHKERRQ(ierr); 544 PetscFunctionReturn(0); 545 } 546 EXTERN_C_END 547