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