Actual source code: mpibstream.c

  1: #define PETSCMAT_DLL

 3:  #include ../src/mat/impls/baij/mpi/mpibaij.h
 4:  #include ../src/mat/impls/baij/seq/bstream/bstream.h



 16: PetscErrorCode MPIBSTRM_create_bstrm(Mat A)
 17: {
 18:   Mat_MPIBAIJ     *a = (Mat_MPIBAIJ *)A->data;
 19:   Mat_SeqBAIJ     *Aij = (Mat_SeqBAIJ*)(a->A->data), *Bij = (Mat_SeqBAIJ*)(a->B->data);
 20:   /* */
 21:   Mat_SeqBSTRM   *bstrmA, *bstrmB;
 22:   PetscInt       MROW = Aij->mbs, bs = a->A->rmap->bs;
 23:   PetscInt       *ai = Aij->i, *bi = Bij->i;
 24:   PetscInt       i,j,k;
 25:   PetscScalar    *aa = Aij->a,*ba = Bij->a;

 27:   PetscInt      bs2,  rbs, cbs, slen, blen;
 29:   PetscScalar **asp;
 30:   PetscScalar **bsp;

 33:   rbs = cbs = bs;
 34:   bs2 = bs*bs;
 35:   blen = ai[MROW]-ai[0];
 36:   slen = blen*bs;

 38:   PetscNewLog(a->A,Mat_SeqBSTRM,&bstrmA);
 39:   a->A->spptr = (void *) bstrmA;
 40:   bstrmA = (Mat_SeqBSTRM*) a->A->spptr;
 41:   bstrmA->rbs = bstrmA->cbs = bs;
 42:   PetscMalloc(bs2*blen*sizeof(PetscScalar), &bstrmA->as);

 44:   PetscMalloc(rbs*sizeof(PetscScalar *), &asp);

 46:   for(i=0;i<rbs;i++) asp[i] = bstrmA->as + i*slen;

 48:   for (k=0; k<blen; k++) {
 49:     for (j=0; j<cbs; j++)
 50:     for (i=0; i<rbs; i++)
 51:         asp[i][k*cbs+j] = aa[k*bs2+j*rbs+i];
 52:   }
 53:   switch (bs){
 54:     case 4:
 55:       a->A->ops->mult  = MatMult_SeqBSTRM_4;
 56:       a->A->ops->sor   = MatSOR_SeqBSTRM_4;
 57:       break;
 58:     case 5:
 59:       a->A->ops->mult  = MatMult_SeqBSTRM_5;
 60:       a->A->ops->sor   = MatSOR_SeqBSTRM_5;
 61:       break;
 62:     default:
 63:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"not supported for block size %D yet",bs);
 64:   }
 65:   PetscFree(asp);

 67: /*.....*/
 68:   blen = bi[MROW]-bi[0];
 69:   slen = blen*bs;
 70:   PetscNewLog(a->B,Mat_SeqBSTRM,&bstrmB);
 71:   a->B->spptr = (void *) bstrmB;
 72:   bstrmB = (Mat_SeqBSTRM*) a->B->spptr;
 73:   bstrmB->rbs = bstrmB->cbs = bs;
 74:   PetscMalloc(bs2*blen*sizeof(PetscScalar), &bstrmB->as);

 76:   PetscMalloc(rbs*sizeof(PetscScalar *), &bsp);

 78:   for(i=0;i<rbs;i++) bsp[i] = bstrmB->as + i*slen;

 80:   for (k=0; k<blen; k++) {
 81:     for (j=0; j<cbs; j++)
 82:     for (i=0; i<rbs; i++)
 83:         bsp[i][k*cbs+j] = ba[k*bs2+j*rbs+i];
 84:   }
 85:   switch (bs){
 86:     case 4:
 87:       a->B->ops->multadd = MatMultAdd_SeqBSTRM_4;
 88:       break;
 89:     case 5:
 90:       a->B->ops->multadd = MatMultAdd_SeqBSTRM_5;
 91:       break;
 92:     default:
 93:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"not supported for block size %D yet",bs);
 94:   }
 95:   PetscFree(bsp);

 97:   return(0);
 98: }


104: PetscErrorCode MatAssemblyEnd_MPIBSTRM(Mat A, MatAssemblyType mode)
105: {

109:   /*
110:     Aij->inode.use = PETSC_FALSE;
111:     Bij->inode.use = PETSC_FALSE;
112:   */
113:   MatAssemblyEnd_MPIBAIJ(A,mode);
114:   if (mode == MAT_FLUSH_ASSEMBLY) return(0);

116:   /* Now calculate the permutation and grouping information. */
117:   MPIBSTRM_create_bstrm(A);
118:   return(0);
119: }


124: PetscErrorCode MatCreateMPIBSTRM(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)
125: {
127:   PetscMPIInt    size;

130:   MatCreate(comm,A);
131:   MatSetSizes(*A,m,n,M,N);
132:   MPI_Comm_size(comm,&size);
133:   if (size > 1) {
134:     MatSetType(*A,MATMPIBSTRM);
135:     MatMPIBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);
136:   } else {
137:     MatSetType(*A,MATSEQBSTRM);
138:     MatSeqBAIJSetPreallocation(*A,bs,d_nz,d_nnz);
139:   }
140:   return(0);
141: }


151: PetscErrorCode MatMPIBAIJSetPreallocation_MPIBSTRM(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
152: {

156:   MatMPIBAIJSetPreallocation_MPIBAIJ(B,bs,d_nz,d_nnz,o_nz,o_nnz);
157: /*
158:   MatConvert_SeqBAIJ_SeqBSTRM(b->A, MATSEQBSTRM, MAT_REUSE_MATRIX, &b->A);
159:   MatConvert_SeqBAIJ_SeqBSTRM(b->B, MATSEQBSTRM, MAT_REUSE_MATRIX, &b->B);
160: */
161:   return(0);
162: }

168: PetscErrorCode MatConvert_MPIBAIJ_MPIBSTRM(Mat A,const MatType type,MatReuse reuse,Mat *newmat)
169: {
171:   Mat            B = *newmat;
172:   Mat_SeqBSTRM  *bstrm;


176:   if (reuse == MAT_INITIAL_MATRIX) {
177:     MatDuplicate(A,MAT_COPY_VALUES,&B);
178:   }

180:   PetscNewLog(B,   Mat_SeqBSTRM,&bstrm);
181:   B->spptr    = (void *) bstrm;

183:   /* Set function pointers for methods that we inherit from AIJ but override.
184:      B->ops->duplicate   = MatDuplicate_BSTRM;
185:      B->ops->mult        = MatMult_BSTRM;
186:      B->ops->destroy     = MatDestroy_MPIBSTRM;
187:   */
188:   B->ops->assemblyend = MatAssemblyEnd_MPIBSTRM;

190:   /* If A has already been assembled, compute the permutation. */
191:   if (A->assembled) {
192:     MPIBSTRM_create_bstrm(B);
193:   }

195:   PetscObjectChangeTypeName( (PetscObject) B, MATMPIBSTRM);
196:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIBAIJSetPreallocation_C",
197:                                      "MatMPIBAIJSetPreallocation_MPIBSTRM",
198:                                      MatMPIBAIJSetPreallocation_MPIBSTRM);
199:   *newmat = B;
200:   return(0);
201: }

207: PetscErrorCode MatCreate_MPIBSTRM(Mat A)
208: {

212:   MatSetType(A,MATMPIBAIJ);
213:   MatConvert_MPIBAIJ_MPIBSTRM(A,MATMPIBSTRM,MAT_REUSE_MATRIX,&A);
214:   return(0);
215: }

221: PetscErrorCode MatCreate_BSTRM(Mat A)
222: {
224:   PetscMPIInt    size;

227:   MPI_Comm_size(((PetscObject)A)->comm,&size);
228:   if (size == 1) {
229:     MatSetType(A,MATSEQBSTRM);
230:   } else {
231:     MatSetType(A,MATMPIBSTRM);
232:   }
233:   return(0);
234: }