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