xref: /petsc/src/mat/impls/localref/mlocalref.c (revision 520db06cfec378d7c34fb86a68534e48a6477cdf)
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   ierr = ISGetLocalSize(is,&m);CHKERRQ(ierr);
101   ierr = ISGetIndices(is,&idx);CHKERRQ(ierr);
102 #if defined(PETSC_USE_DEBUG)
103   {
104     PetscInt i;
105     for (i=0; i<m; i++) {
106       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);
107     }
108   }
109 #endif
110   ierr = PetscMalloc(m*sizeof(PetscInt),&idxm);CHKERRQ(ierr);
111   if (ltog) {
112     ierr = ISLocalToGlobalMappingApply(ltog,m,idx,idxm);CHKERRQ(ierr);
113   } else {
114     ierr = PetscMemcpy(idxm,idx,m*sizeof(PetscInt));CHKERRQ(ierr);
115   }
116   ierr = ISLocalToGlobalMappingCreate(((PetscObject)is)->comm,m,idxm,PETSC_OWN_POINTER,cltog);CHKERRQ(ierr);
117   ierr = ISRestoreIndices(is,&idx);CHKERRQ(ierr);
118   PetscFunctionReturn(0);
119 }
120 
121 #undef __FUNCT__
122 #define __FUNCT__ "ISL2GComposeBlock"
123 static PetscErrorCode ISL2GComposeBlock(IS is,ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping *cltog)
124 {
125   PetscErrorCode ierr;
126   const PetscInt *idx;
127   PetscInt m,*idxm;
128 
129   PetscFunctionBegin;
130   ierr = ISBlockGetLocalSize(is,&m);CHKERRQ(ierr);
131   ierr = ISBlockGetIndices(is,&idx);CHKERRQ(ierr);
132 #if defined(PETSC_USE_DEBUG)
133   {
134     PetscInt i;
135     for (i=0; i<m; i++) {
136       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);
137     }
138   }
139 #endif
140   ierr = PetscMalloc(m*sizeof(PetscInt),&idxm);CHKERRQ(ierr);
141   if (ltog) {
142     ierr = ISLocalToGlobalMappingApply(ltog,m,idx,idxm);CHKERRQ(ierr);
143   } else {
144     ierr = PetscMemcpy(idxm,idx,m*sizeof(PetscInt));CHKERRQ(ierr);
145   }
146   ierr = ISLocalToGlobalMappingCreate(((PetscObject)is)->comm,m,idxm,PETSC_OWN_POINTER,cltog);CHKERRQ(ierr);
147   ierr = ISBlockRestoreIndices(is,&idx);CHKERRQ(ierr);
148   PetscFunctionReturn(0);
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "MatDestroy_LocalRef"
153 static PetscErrorCode MatDestroy_LocalRef(Mat B)
154 {
155   PetscErrorCode ierr;
156 
157   PetscFunctionBegin;
158   ierr = PetscFree(B->data);CHKERRQ(ierr);
159   PetscFunctionReturn(0);
160 }
161 
162 
163 #undef __FUNCT__
164 #define __FUNCT__ "MatCreateLocalRef"
165 /*@
166    MatCreateLocalRef - Gets a logical reference to a local submatrix, for use in assembly
167 
168    Not Collective
169 
170    Input Arguments:
171 + A - Full matrix, generally parallel
172 . isrow - Local index set for the rows
173 - iscol - Local index set for the columns
174 
175    Output Arguments:
176 . newmat - New serial Mat
177 
178    Level: developer
179 
180    Notes:
181    Most will use MatGetLocalSubMatrix() which returns a real matrix corresponding to the local
182    block if it available, such as with matrix formats that store these blocks separately.
183 
184    The new matrix forwards MatSetValuesLocal() and MatSetValuesBlockedLocal() to the global system.
185    In general, it does not define MatMult() or any other functions.  Local submatrices can be nested.
186 
187 .seealso: MatSetValuesLocal(), MatSetValuesBlockedLocal(), MatGetLocalSubMatrix(), MatCreateSubMatrix()
188 @*/
189 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLocalRef(Mat A,IS isrow,IS iscol,Mat *newmat)
190 {
191   PetscErrorCode ierr;
192   Mat_LocalRef   *lr;
193   Mat            B;
194   PetscInt       m,n;
195   PetscBool      islr;
196 
197   PetscFunctionBegin;
198   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
199   PetscValidHeaderSpecific(isrow,IS_CLASSID,2);
200   PetscValidHeaderSpecific(iscol,IS_CLASSID,3);
201   PetscValidPointer(newmat,4);
202   *newmat = 0;
203 
204   ierr = MatCreate(PETSC_COMM_SELF,&B);CHKERRQ(ierr);
205   ierr = ISGetLocalSize(isrow,&m);CHKERRQ(ierr);
206   ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
207   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
208   ierr = PetscObjectChangeTypeName((PetscObject)B,MATLOCALREF);CHKERRQ(ierr);
209 
210   B->ops->destroy = MatDestroy_LocalRef;
211 
212   ierr = PetscNewLog(B,Mat_LocalRef,&lr);CHKERRQ(ierr);
213   B->data = (void*)lr;
214 
215   ierr = PetscTypeCompare((PetscObject)A,MATLOCALREF,&islr);CHKERRQ(ierr);
216   if (islr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Nesting MatLocalRef not implemented yet");
217   else {
218     ISLocalToGlobalMapping rltog,cltog;
219     PetscInt abs,rbs,cbs;
220 
221     /* This does not increase the reference count because MatLocalRef is not allowed to live longer than its parent */
222     lr->Top = A;
223 
224     /* The immediate parent is the top level and we will translate directly to global indices */
225     lr->SetValues        = MatSetValues;
226     lr->SetValuesBlocked = MatSetValuesBlocked;
227 
228     B->ops->setvalueslocal = MatSetValuesLocal_LocalRef_Scalar;
229     ierr = ISL2GCompose(isrow,A->rmapping,&rltog);CHKERRQ(ierr);
230     if (isrow == iscol && A->rmapping == A->cmapping) {
231       ierr = PetscObjectReference((PetscObject)rltog);CHKERRQ(ierr);
232       cltog = rltog;
233     } else {
234       ierr = ISL2GCompose(iscol,A->cmapping,&cltog);CHKERRQ(ierr);
235     }
236     ierr = MatSetLocalToGlobalMapping(B,rltog,cltog);CHKERRQ(ierr);
237     ierr = ISLocalToGlobalMappingDestroy(rltog);CHKERRQ(ierr);
238     ierr = ISLocalToGlobalMappingDestroy(cltog);CHKERRQ(ierr);
239 
240     ierr = MatGetBlockSize(A,&abs);CHKERRQ(ierr);
241     ierr = ISGetBlockSize(isrow,&rbs);CHKERRQ(ierr);
242     ierr = ISGetBlockSize(iscol,&cbs);CHKERRQ(ierr);
243     if (rbs == cbs) {           /* submatrix has block structure, so user can insert values with blocked interface */
244       ierr = PetscLayoutSetBlockSize(A->rmap,abs);CHKERRQ(ierr);
245       ierr = PetscLayoutSetBlockSize(A->rmap,abs);CHKERRQ(ierr);
246       if (abs != rbs) {
247         /* Top-level matrix has different block size, so we have to call its scalar insertion interface */
248         B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Scalar;
249       } else {
250         /* Block sizes match so we can forward values to the top level using the block interface */
251         B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Block;
252         ierr = ISL2GComposeBlock(isrow,A->rbmapping,&rltog);CHKERRQ(ierr);
253         if (isrow == iscol && A->rbmapping == A->cbmapping) {
254           ierr =  PetscObjectReference((PetscObject)rltog);CHKERRQ(ierr);
255           cltog = rltog;
256         } else {
257           ierr = ISL2GComposeBlock(iscol,A->cbmapping,&cltog);CHKERRQ(ierr);
258         }
259         ierr = MatSetLocalToGlobalMappingBlock(B,rltog,cltog);CHKERRQ(ierr);
260         ierr = ISLocalToGlobalMappingDestroy(rltog);CHKERRQ(ierr);
261         ierr = ISLocalToGlobalMappingDestroy(cltog);CHKERRQ(ierr);
262       }
263     }
264   }
265   *newmat = B;
266   PetscFunctionReturn(0);
267 }
268