xref: /petsc/src/mat/impls/localref/mlocalref.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1 
2 #include <petsc/private/matimpl.h>          /*I "petscmat.h" I*/
3 
4 typedef struct {
5   Mat Top;
6   PetscBool rowisblock;
7   PetscBool colisblock;
8   PetscErrorCode (*SetValues)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
9   PetscErrorCode (*SetValuesBlocked)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
10 } Mat_LocalRef;
11 
12 /* These need to be macros because they use sizeof */
13 #define IndexSpaceGet(buf,nrow,ncol,irowm,icolm) do {                   \
14     if (nrow + ncol > (PetscInt)(sizeof(buf)/sizeof(buf[0]))) {         \
15       CHKERRQ(PetscMalloc2(nrow,&irowm,ncol,&icolm)); \
16     } else {                                                            \
17       irowm = &buf[0];                                                  \
18       icolm = &buf[nrow];                                               \
19     }                                                                   \
20 } while (0)
21 
22 #define IndexSpaceRestore(buf,nrow,ncol,irowm,icolm) do {       \
23     if (nrow + ncol > (PetscInt)(sizeof(buf)/sizeof(buf[0]))) { \
24       CHKERRQ(PetscFree2(irowm,icolm));             \
25     }                                                           \
26 } while (0)
27 
28 static void BlockIndicesExpand(PetscInt n,const PetscInt idx[],PetscInt bs,PetscInt idxm[])
29 {
30   PetscInt i,j;
31   for (i=0; i<n; i++) {
32     for (j=0; j<bs; j++) {
33       idxm[i*bs+j] = idx[i]*bs + j;
34     }
35   }
36 }
37 
38 static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Block(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
39 {
40   Mat_LocalRef   *lr = (Mat_LocalRef*)A->data;
41   PetscInt       buf[4096],*irowm=NULL,*icolm; /* suppress maybe-uninitialized warning */
42 
43   PetscFunctionBegin;
44   if (!nrow || !ncol) PetscFunctionReturn(0);
45   IndexSpaceGet(buf,nrow,ncol,irowm,icolm);
46   CHKERRQ(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,nrow,irow,irowm));
47   CHKERRQ(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,ncol,icol,icolm));
48   CHKERRQ((*lr->SetValuesBlocked)(lr->Top,nrow,irowm,ncol,icolm,y,addv));
49   IndexSpaceRestore(buf,nrow,ncol,irowm,icolm);
50   PetscFunctionReturn(0);
51 }
52 
53 static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Scalar(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
54 {
55   Mat_LocalRef   *lr = (Mat_LocalRef*)A->data;
56   PetscInt       rbs,cbs,buf[4096],*irowm,*icolm;
57 
58   PetscFunctionBegin;
59   CHKERRQ(MatGetBlockSizes(A,&rbs,&cbs));
60   IndexSpaceGet(buf,nrow*rbs,ncol*cbs,irowm,icolm);
61   BlockIndicesExpand(nrow,irow,rbs,irowm);
62   BlockIndicesExpand(ncol,icol,cbs,icolm);
63   CHKERRQ(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,nrow*rbs,irowm,irowm));
64   CHKERRQ(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,ncol*cbs,icolm,icolm));
65   CHKERRQ((*lr->SetValues)(lr->Top,nrow*rbs,irowm,ncol*cbs,icolm,y,addv));
66   IndexSpaceRestore(buf,nrow*rbs,ncol*cbs,irowm,icolm);
67   PetscFunctionReturn(0);
68 }
69 
70 static PetscErrorCode MatSetValuesLocal_LocalRef_Scalar(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
71 {
72   Mat_LocalRef   *lr = (Mat_LocalRef*)A->data;
73   PetscInt       buf[4096],*irowm,*icolm;
74 
75   PetscFunctionBegin;
76   IndexSpaceGet(buf,nrow,ncol,irowm,icolm);
77   /* If the row IS defining this submatrix was an ISBLOCK, then the unblocked LGMapApply is the right one to use.  If
78    * instead it was (say) an ISSTRIDE with a block size > 1, then we need to use LGMapApplyBlock */
79   if (lr->rowisblock) {
80     CHKERRQ(ISLocalToGlobalMappingApply(A->rmap->mapping,nrow,irow,irowm));
81   } else {
82     CHKERRQ(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,nrow,irow,irowm));
83   }
84   /* As above, but for the column IS. */
85   if (lr->colisblock) {
86     CHKERRQ(ISLocalToGlobalMappingApply(A->cmap->mapping,ncol,icol,icolm));
87   } else {
88     CHKERRQ(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,ncol,icol,icolm));
89   }
90   CHKERRQ((*lr->SetValues)(lr->Top,nrow,irowm,ncol,icolm,y,addv));
91   IndexSpaceRestore(buf,nrow,ncol,irowm,icolm);
92   PetscFunctionReturn(0);
93 }
94 
95 /* Compose an IS with an ISLocalToGlobalMapping to map from IS source indices to global indices */
96 static PetscErrorCode ISL2GCompose(IS is,ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping *cltog)
97 {
98   const PetscInt *idx;
99   PetscInt       m,*idxm;
100   PetscInt       bs;
101   PetscBool      isblock;
102 
103   PetscFunctionBegin;
104   PetscValidHeaderSpecific(is,IS_CLASSID,1);
105   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,2);
106   PetscValidPointer(cltog,3);
107   CHKERRQ(PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock));
108   if (isblock) {
109     PetscInt lbs;
110 
111     CHKERRQ(ISGetBlockSize(is,&bs));
112     CHKERRQ(ISLocalToGlobalMappingGetBlockSize(ltog,&lbs));
113     if (bs == lbs) {
114       CHKERRQ(ISGetLocalSize(is,&m));
115       m    = m/bs;
116       CHKERRQ(ISBlockGetIndices(is,&idx));
117       CHKERRQ(PetscMalloc1(m,&idxm));
118       CHKERRQ(ISLocalToGlobalMappingApplyBlock(ltog,m,idx,idxm));
119       CHKERRQ(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),bs,m,idxm,PETSC_OWN_POINTER,cltog));
120       CHKERRQ(ISBlockRestoreIndices(is,&idx));
121       PetscFunctionReturn(0);
122     }
123   }
124   CHKERRQ(ISGetLocalSize(is,&m));
125   CHKERRQ(ISGetIndices(is,&idx));
126   CHKERRQ(ISGetBlockSize(is,&bs));
127   CHKERRQ(PetscMalloc1(m,&idxm));
128   if (ltog) {
129     CHKERRQ(ISLocalToGlobalMappingApply(ltog,m,idx,idxm));
130   } else {
131     CHKERRQ(PetscArraycpy(idxm,idx,m));
132   }
133   CHKERRQ(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),bs,m,idxm,PETSC_OWN_POINTER,cltog));
134   CHKERRQ(ISRestoreIndices(is,&idx));
135   PetscFunctionReturn(0);
136 }
137 
138 static PetscErrorCode ISL2GComposeBlock(IS is,ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping *cltog)
139 {
140   const PetscInt *idx;
141   PetscInt       m,*idxm,bs;
142 
143   PetscFunctionBegin;
144   PetscValidHeaderSpecific(is,IS_CLASSID,1);
145   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,2);
146   PetscValidPointer(cltog,3);
147   CHKERRQ(ISBlockGetLocalSize(is,&m));
148   CHKERRQ(ISBlockGetIndices(is,&idx));
149   CHKERRQ(ISLocalToGlobalMappingGetBlockSize(ltog,&bs));
150   CHKERRQ(PetscMalloc1(m,&idxm));
151   if (ltog) {
152     CHKERRQ(ISLocalToGlobalMappingApplyBlock(ltog,m,idx,idxm));
153   } else {
154     CHKERRQ(PetscArraycpy(idxm,idx,m));
155   }
156   CHKERRQ(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),bs,m,idxm,PETSC_OWN_POINTER,cltog));
157   CHKERRQ(ISBlockRestoreIndices(is,&idx));
158   PetscFunctionReturn(0);
159 }
160 
161 static PetscErrorCode MatDestroy_LocalRef(Mat B)
162 {
163   PetscFunctionBegin;
164   CHKERRQ(PetscFree(B->data));
165   PetscFunctionReturn(0);
166 }
167 
168 /*@
169    MatCreateLocalRef - Gets a logical reference to a local submatrix, for use in assembly
170 
171    Not Collective
172 
173    Input Parameters:
174 + A - Full matrix, generally parallel
175 . isrow - Local index set for the rows
176 - iscol - Local index set for the columns
177 
178    Output Parameter:
179 . newmat - New serial Mat
180 
181    Level: developer
182 
183    Notes:
184    Most will use MatGetLocalSubMatrix() which returns a real matrix corresponding to the local
185    block if it available, such as with matrix formats that store these blocks separately.
186 
187    The new matrix forwards MatSetValuesLocal() and MatSetValuesBlockedLocal() to the global system.
188    In general, it does not define MatMult() or any other functions.  Local submatrices can be nested.
189 
190 .seealso: MatSetValuesLocal(), MatSetValuesBlockedLocal(), MatGetLocalSubMatrix(), MatCreateSubMatrix()
191 @*/
192 PetscErrorCode  MatCreateLocalRef(Mat A,IS isrow,IS iscol,Mat *newmat)
193 {
194   Mat_LocalRef   *lr;
195   Mat            B;
196   PetscInt       m,n;
197   PetscBool      islr;
198 
199   PetscFunctionBegin;
200   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
201   PetscValidHeaderSpecific(isrow,IS_CLASSID,2);
202   PetscValidHeaderSpecific(iscol,IS_CLASSID,3);
203   PetscValidPointer(newmat,4);
204   PetscCheckFalse(!A->rmap->mapping,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Matrix must have local to global mapping provided before this call");
205   *newmat = NULL;
206 
207   CHKERRQ(MatCreate(PETSC_COMM_SELF,&B));
208   CHKERRQ(ISGetLocalSize(isrow,&m));
209   CHKERRQ(ISGetLocalSize(iscol,&n));
210   CHKERRQ(MatSetSizes(B,m,n,m,n));
211   CHKERRQ(PetscObjectChangeTypeName((PetscObject)B,MATLOCALREF));
212   CHKERRQ(MatSetUp(B));
213 
214   B->ops->destroy = MatDestroy_LocalRef;
215 
216   CHKERRQ(PetscNewLog(B,&lr));
217   B->data = (void*)lr;
218 
219   CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATLOCALREF,&islr));
220   if (islr) {
221     Mat_LocalRef *alr = (Mat_LocalRef*)A->data;
222     lr->Top = alr->Top;
223   } else {
224     /* This does not increase the reference count because MatLocalRef is not allowed to live longer than its parent */
225     lr->Top = A;
226   }
227   {
228     ISLocalToGlobalMapping rltog,cltog;
229     PetscInt               arbs,acbs,rbs,cbs;
230 
231     /* We will translate directly to global indices for the top level */
232     lr->SetValues        = MatSetValues;
233     lr->SetValuesBlocked = MatSetValuesBlocked;
234 
235     B->ops->setvalueslocal = MatSetValuesLocal_LocalRef_Scalar;
236 
237     CHKERRQ(ISL2GCompose(isrow,A->rmap->mapping,&rltog));
238     if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) {
239       CHKERRQ(PetscObjectReference((PetscObject)rltog));
240       cltog = rltog;
241     } else {
242       CHKERRQ(ISL2GCompose(iscol,A->cmap->mapping,&cltog));
243     }
244     /* Remember if the ISes we used to pull out the submatrix are of type ISBLOCK.  Will be used later in
245      * MatSetValuesLocal_LocalRef_Scalar. */
246     CHKERRQ(PetscObjectTypeCompare((PetscObject)isrow,ISBLOCK,&lr->rowisblock));
247     CHKERRQ(PetscObjectTypeCompare((PetscObject)iscol,ISBLOCK,&lr->colisblock));
248     CHKERRQ(MatSetLocalToGlobalMapping(B,rltog,cltog));
249     CHKERRQ(ISLocalToGlobalMappingDestroy(&rltog));
250     CHKERRQ(ISLocalToGlobalMappingDestroy(&cltog));
251 
252     CHKERRQ(MatGetBlockSizes(A,&arbs,&acbs));
253     CHKERRQ(ISGetBlockSize(isrow,&rbs));
254     CHKERRQ(ISGetBlockSize(iscol,&cbs));
255     /* Always support block interface insertion on submatrix */
256     CHKERRQ(PetscLayoutSetBlockSize(B->rmap,rbs));
257     CHKERRQ(PetscLayoutSetBlockSize(B->cmap,cbs));
258     if (arbs != rbs || acbs != cbs || (arbs == 1 && acbs == 1)) {
259       /* Top-level matrix has different block size, so we have to call its scalar insertion interface */
260       B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Scalar;
261     } else {
262       /* Block sizes match so we can forward values to the top level using the block interface */
263       B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Block;
264 
265       CHKERRQ(ISL2GComposeBlock(isrow,A->rmap->mapping,&rltog));
266       if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) {
267         CHKERRQ(PetscObjectReference((PetscObject)rltog));
268         cltog = rltog;
269       } else {
270         CHKERRQ(ISL2GComposeBlock(iscol,A->cmap->mapping,&cltog));
271       }
272       CHKERRQ(MatSetLocalToGlobalMapping(B,rltog,cltog));
273       CHKERRQ(ISLocalToGlobalMappingDestroy(&rltog));
274       CHKERRQ(ISLocalToGlobalMappingDestroy(&cltog));
275     }
276   }
277   *newmat = B;
278   PetscFunctionReturn(0);
279 }
280