Actual source code: mpisbstream.c
1: #define PETSCMAT_DLL
3: #include ../src/mat/impls/sbaij/mpi/mpisbaij.h
4: #include ../src/mat/impls/sbaij/seq/sbstream/sbstream.h
17: PetscErrorCode MPISBSTRM_create_sbstrm(Mat A)
18: {
19: Mat_MPISBAIJ *a = (Mat_MPISBAIJ *)A->data;
20: Mat_SeqSBAIJ *Aij = (Mat_SeqSBAIJ*)(a->A->data), *Bij = (Mat_SeqSBAIJ*)(a->B->data);
21: /*
22: */
23: Mat_SeqSBSTRM *sbstrmA, *sbstrmB;
24: PetscInt MROW = Aij->mbs, bs = a->A->rmap->bs;
26: /* PetscInt m = A->rmap->n;*/ /* Number of rows in the matrix. */
27: /* PetscInt nd = a->A->cmap->n;*/ /* number of columns in diagonal portion */
28: PetscInt *ai = Aij->i, *bi = Bij->i; /* From the CSR representation; points to the beginning of each row. */
29: PetscInt i,j,k;
30: PetscScalar *aa = Aij->a,*ba = Bij->a;
32: PetscInt bs2, rbs, cbs, slen, blen;
34: PetscScalar **asp;
35: PetscScalar **bsp;
38: /* printf(" --- in MPISBSTRM_create_sbstrm, m=%d, nd=%d, bs=%d, MROW=%d\n", m,nd,bs,MROW); */
40: rbs = cbs = bs;
41: bs2 = bs*bs;
42: blen = ai[MROW]-ai[0];
43: slen = blen*bs;
44:
45: /* printf(" --- blen=%d, slen=%d\n", blen, slen); */
47: PetscNewLog(a->A,Mat_SeqSBSTRM,&sbstrmA);
48: a->A->spptr = (void *) sbstrmA;
49: sbstrmA = (Mat_SeqSBSTRM*) a->A->spptr;
50: sbstrmA->rbs = sbstrmA->cbs = bs;
51: PetscMalloc(bs2*blen*sizeof(PetscScalar), &sbstrmA->as);
53: PetscMalloc(rbs*sizeof(PetscScalar *), &asp);
55: for(i=0;i<rbs;i++) asp[i] = sbstrmA->as + i*slen;
57: for (k=0; k<blen; k++) {
58: for (j=0; j<cbs; j++)
59: for (i=0; i<rbs; i++)
60: asp[i][k*cbs+j] = aa[k*bs2+j*rbs+i];
61: }
62: switch (bs){
63: case 4:
64: a->A->ops->mult = MatMult_SeqSBSTRM_4;
65: a->A->ops->multtranspose = MatMult_SeqSBSTRM_4;
66: /* a->A->ops->sor = MatSOR_SeqSBSTRM_4; */
67: break;
68: case 5:
69: a->A->ops->mult = MatMult_SeqSBSTRM_5;
70: a->A->ops->multtranspose = MatMult_SeqSBSTRM_5;
71: /* a->A->ops->sor = MatSOR_SeqSBSTRM_5; */
72: break;
73: default:
74: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"not supported for block size %D yet",bs);
75: }
76: PetscFree(asp);
79: /*.....*/
80: blen = bi[MROW]-bi[0];
81: slen = blen*bs;
82: PetscNewLog(a->B,Mat_SeqSBSTRM,&sbstrmB);
83: a->B->spptr = (void *) sbstrmB;
84: sbstrmB = (Mat_SeqSBSTRM*) a->B->spptr;
85: sbstrmB->rbs = sbstrmB->cbs = bs;
86: PetscMalloc(bs2*blen*sizeof(PetscScalar), &sbstrmB->as);
88: PetscMalloc(rbs*sizeof(PetscScalar *), &bsp);
90: for(i=0;i<rbs;i++) bsp[i] = sbstrmB->as + i*slen;
92: for (k=0; k<blen; k++) {
93: for (j=0; j<cbs; j++)
94: for (i=0; i<rbs; i++)
95: bsp[i][k*cbs+j] = ba[k*bs2+j*rbs+i];
96: }
97: switch (bs){
98: case 4:
99: /* a->B->ops->mult = MatMult_SeqSBSTRM_4; */
100: a->B->ops->multtranspose = MatMultTranspose_SeqBSTRM_4;
101: a->B->ops->multadd = MatMultAdd_SeqBSTRM_4;
102: /* a->B->ops->multtransposeadd = MatMultAdd_SeqSBSTRM_4; */
103: break;
104: case 5:
105: /* a->B->ops->mult = MatMult_SeqSBSTRM_5; */
106: a->B->ops->multtranspose = MatMultTranspose_SeqBSTRM_5;
107: a->B->ops->multadd = MatMultAdd_SeqBSTRM_5;
108: /* a->B->ops->multtransposeadd = MatMultAdd_SeqSBSTRM_5; */
109: break;
110: default:
111: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"not supported for block size %D yet",bs);
112: }
113: PetscFree(bsp);
115: return(0);
116: }
122: PetscErrorCode MatAssemblyEnd_MPISBSTRM(Mat A, MatAssemblyType mode)
123: {
127: /*
128: Aij->inode.use = PETSC_FALSE;
129: Bij->inode.use = PETSC_FALSE;
130: */
131: MatAssemblyEnd_MPISBAIJ(A,mode);
132: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
134: /* Now calculate the permutation and grouping information. */
135: MPISBSTRM_create_sbstrm(A);
136: return(0);
137: }
142: PetscErrorCode MatCreateMPISBSTRM(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[],Mat *A)
143: {
145: PetscMPIInt size;
148: MatCreate(comm,A);
149: MatSetSizes(*A,m,n,M,N);
150: MPI_Comm_size(comm,&size);
151: if (size > 1) {
152: MatSetType(*A,MATMPISBSTRM);
153: MatMPISBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);
154: } else {
155: MatSetType(*A,MATSEQSBSTRM);
156: MatSeqSBAIJSetPreallocation(*A,bs,d_nz,d_nnz);
157: }
158: return(0);
159: }
169: PetscErrorCode MatMPISBAIJSetPreallocation_MPISBSTRM(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
170: {
174: MatMPISBAIJSetPreallocation_MPISBAIJ(B,bs,d_nz,d_nnz,o_nz,o_nnz);
175: /*
176: MatConvert_SeqSBAIJ_SeqSBSTRM(b->A, MATSEQSBSTRM, MAT_REUSE_MATRIX, &b->A);
177: MatConvert_SeqSBAIJ_SeqSBSTRM(b->B, MATSEQSBSTRM, MAT_REUSE_MATRIX, &b->B);
178: */
179: return(0);
180: }
186: PetscErrorCode MatConvert_MPISBAIJ_MPISBSTRM(Mat A,const MatType type,MatReuse reuse,Mat *newmat)
187: {
189: Mat B = *newmat;
190: Mat_SeqSBSTRM *sbstrm;
194: if (reuse == MAT_INITIAL_MATRIX) {
195: MatDuplicate(A,MAT_COPY_VALUES,&B);
196: }
197: /* printf(" --- in MatConvert_MPISBAIJ_MPISBSTRM -- 1 \n"); */
199: PetscNewLog(B, Mat_SeqSBSTRM,&sbstrm);
200: B->spptr = (void *) sbstrm;
202: /* Set function pointers for methods that we inherit from AIJ but override.
203: B->ops->duplicate = MatDuplicate_SBSTRM;
204: B->ops->mult = MatMult_SBSTRM;
205: B->ops->destroy = MatDestroy_MPISBSTRM;
206: */
207: B->ops->assemblyend = MatAssemblyEnd_MPISBSTRM;
209: /* If A has already been assembled, compute the permutation. */
210: if (A->assembled) {
211: MPISBSTRM_create_sbstrm(B);
212: }
214: PetscObjectChangeTypeName( (PetscObject) B, MATMPISBSTRM);
215: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
216: "MatMPISBAIJSetPreallocation_MPISBSTRM",
217: MatMPISBAIJSetPreallocation_MPISBSTRM);
218: *newmat = B;
219: return(0);
220: }
226: PetscErrorCode MatCreate_MPISBSTRM(Mat A)
227: {
231: MatSetType(A,MATMPISBAIJ);
232: MatConvert_MPISBAIJ_MPISBSTRM(A,MATMPISBSTRM,MAT_REUSE_MATRIX,&A);
233: return(0);
234: }
241: PetscErrorCode MatCreate_SBSTRM(Mat A)
242: {
244: PetscMPIInt size;
247: MPI_Comm_size(((PetscObject)A)->comm,&size);
248: if (size == 1) {
249: MatSetType(A,MATSEQSBSTRM);
250: } else {
251: MatSetType(A,MATMPISBSTRM);
252: }
253: return(0);
254: }