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