Actual source code: mpibstream.c
petsc-3.7.3 2016-08-01
1: #define PETSCMAT_DLL
3: #include <../src/mat/impls/baij/mpi/mpibaij.h>
4: #include <../src/mat/impls/baij/seq/bstream/bstream.h>
6: extern PetscErrorCode MatMult_SeqBSTRM_4(Mat,Vec,Vec);
7: extern PetscErrorCode MatMult_SeqBSTRM_5(Mat,Vec,Vec);
8: extern PetscErrorCode MatMultAdd_SeqBSTRM_4(Mat,Vec,Vec,Vec);
9: extern PetscErrorCode MatMultAdd_SeqBSTRM_5(Mat,Vec,Vec,Vec);
10: extern PetscErrorCode MatSOR_SeqBSTRM_4(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11: extern PetscErrorCode MatSOR_SeqBSTRM_5(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
16: PetscErrorCode MatMPIBSTRM_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,&bstrmA);
39: a->A->spptr = (void*) bstrmA;
40: bstrmA = (Mat_SeqBSTRM*) a->A->spptr;
41: bstrmA->rbs = bstrmA->cbs = bs;
42: PetscMalloc1(bs2*blen, &bstrmA->as);
44: PetscMalloc1(rbs, &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++) asp[i][k*cbs+j] = aa[k*bs2+j*rbs+i];
51: }
52: }
54: switch (bs) {
55: case 4:
56: a->A->ops->mult = MatMult_SeqBSTRM_4;
57: a->A->ops->sor = MatSOR_SeqBSTRM_4;
58: break;
59: case 5:
60: a->A->ops->mult = MatMult_SeqBSTRM_5;
61: a->A->ops->sor = MatSOR_SeqBSTRM_5;
62: break;
63: default:
64: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"not supported for block size %D yet",bs);
65: }
66: PetscFree(asp);
68: /*.....*/
69: blen = bi[MROW]-bi[0];
70: slen = blen*bs;
71: PetscNewLog(a->B,&bstrmB);
73: a->B->spptr = (void*) bstrmB;
74: bstrmB = (Mat_SeqBSTRM*) a->B->spptr;
75: bstrmB->rbs = bstrmB->cbs = bs;
77: PetscMalloc1(bs2*blen, &bstrmB->as);
78: PetscMalloc1(rbs, &bsp);
80: for (i=0; i<rbs; i++) bsp[i] = bstrmB->as + i*slen;
82: for (k=0; k<blen; k++) {
83: for (j=0; j<cbs; j++) {
84: for (i=0; i<rbs; i++) bsp[i][k*cbs+j] = ba[k*bs2+j*rbs+i];
85: }
86: }
88: switch (bs) {
89: case 4:
90: a->B->ops->multadd = MatMultAdd_SeqBSTRM_4;
91: break;
92: case 5:
93: a->B->ops->multadd = MatMultAdd_SeqBSTRM_5;
94: break;
95: default:
96: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"not supported for block size %D yet",bs);
97: }
98: PetscFree(bsp);
99: return(0);
100: }
102: extern PetscErrorCode MatAssemblyEnd_MPIBAIJ(Mat,MatAssemblyType);
106: PetscErrorCode MatAssemblyEnd_MPIBSTRM(Mat A, MatAssemblyType mode)
107: {
111: /*
112: Aij->inode.use = PETSC_FALSE;
113: Bij->inode.use = PETSC_FALSE;
114: */
115: MatAssemblyEnd_MPIBAIJ(A,mode);
116: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
118: /* Now calculate the permutation and grouping information. */
119: MatMPIBSTRM_create_bstrm(A);
120: return(0);
121: }
126: 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)
127: {
129: PetscMPIInt size;
132: MatCreate(comm,A);
133: MatSetSizes(*A,m,n,M,N);
134: MPI_Comm_size(comm,&size);
135: if (size > 1) {
136: MatSetType(*A,MATMPIBSTRM);
137: MatMPIBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);
138: } else {
139: MatSetType(*A,MATSEQBSTRM);
140: MatSeqBAIJSetPreallocation(*A,bs,d_nz,d_nnz);
141: }
142: return(0);
143: }
145: PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqBSTRM(Mat,MatType,MatReuse,Mat*);
146: extern PetscErrorCode MatMPIBAIJSetPreallocation_MPIBAIJ(Mat,PetscInt,PetscInt,const PetscInt *,PetscInt,const PetscInt *);
150: PetscErrorCode MatMPIBAIJSetPreallocation_MPIBSTRM(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
151: {
155: MatMPIBAIJSetPreallocation_MPIBAIJ(B,bs,d_nz,d_nnz,o_nz,o_nnz);
156: return(0);
157: }
161: PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPIBSTRM(Mat A,MatType type,MatReuse reuse,Mat *newmat)
162: {
164: Mat B = *newmat;
165: Mat_SeqBSTRM *bstrm;
168: if (reuse == MAT_INITIAL_MATRIX) {
169: MatDuplicate(A,MAT_COPY_VALUES,&B);
170: }
172: PetscNewLog(B,&bstrm);
173: B->spptr = (void*) bstrm;
175: /* Set function pointers for methods that we inherit from AIJ but override.
176: B->ops->duplicate = MatDuplicate_BSTRM;
177: B->ops->mult = MatMult_BSTRM;
178: B->ops->destroy = MatDestroy_MPIBSTRM;
179: */
180: B->ops->assemblyend = MatAssemblyEnd_MPIBSTRM;
182: /* If A has already been assembled, compute the permutation. */
183: if (A->assembled) {
184: MatMPIBSTRM_create_bstrm(B);
185: }
187: PetscObjectChangeTypeName((PetscObject) B, MATMPIBSTRM);
188: PetscObjectComposeFunction((PetscObject)B,"MatMPIBAIJSetPreallocation_C",MatMPIBAIJSetPreallocation_MPIBSTRM);
189: *newmat = B;
190: return(0);
191: }
195: PETSC_EXTERN PetscErrorCode MatCreate_MPIBSTRM(Mat A)
196: {
200: MatSetType(A,MATMPIBAIJ);
201: MatConvert_MPIBAIJ_MPIBSTRM(A,MATMPIBSTRM,MAT_INPLACE_MATRIX,&A);
202: return(0);
203: }
207: PETSC_EXTERN PetscErrorCode MatCreate_BSTRM(Mat A)
208: {
210: PetscMPIInt size;
213: MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
214: if (size == 1) {
215: MatSetType(A,MATSEQBSTRM);
216: } else {
217: MatSetType(A,MATMPIBSTRM);
218: }
219: return(0);
220: }