xref: /petsc/src/dm/impls/sliced/sliced.c (revision 47c6ae997ffd1b2afd66b6474dff5950ae8613d1)
1 #define PETSCDM_DLL
2 
3 #include "petscda.h"     /*I      "petscda.h"     I*/
4 #include "petscmat.h"    /*I      "petscmat.h"    I*/
5 #include "private/dmimpl.h"    /*I      "petscmat.h"    I*/
6 
7 /* CSR storage of the nonzero structure of a bs*bs matrix */
8 typedef struct {
9   PetscInt bs,nz,*i,*j;
10 } SlicedBlockFills;
11 
12 typedef struct  {
13   Vec              globalvector;
14   PetscInt         bs,n,N,Nghosts,*ghosts;
15   PetscInt         d_nz,o_nz,*d_nnz,*o_nnz;
16   SlicedBlockFills *dfill,*ofill;
17 } DM_Sliced;
18 
19 #undef __FUNCT__
20 #define __FUNCT__ "SlicedGetMatrix"
21 /*@C
22     SlicedGetMatrix - Creates a matrix with the correct parallel layout required for
23       computing the Jacobian on a function defined using the informatin in Sliced.
24 
25     Collective on Sliced
26 
27     Input Parameter:
28 +   slice - the slice object
29 -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ,
30             or any type which inherits from one of these (such as MATAIJ, MATLUSOL, etc.).
31 
32     Output Parameters:
33 .   J  - matrix with the correct nonzero preallocation
34         (obviously without the correct Jacobian values)
35 
36     Level: advanced
37 
38     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
39        do not need to do it yourself.
40 
41 .seealso ISColoringView(), ISColoringGetIS(), MatFDColoringCreate(), DASetBlockFills()
42 
43 @*/
44 PetscErrorCode PETSCDM_DLLEXPORT SlicedGetMatrix(DM dm, const MatType mtype,Mat *J)
45 {
46   PetscErrorCode         ierr;
47   PetscInt               *globals,*sd_nnz,*so_nnz,rstart,bs,i;
48   ISLocalToGlobalMapping lmap,blmap;
49   void                   (*aij)(void) = PETSC_NULL;
50   DM_Sliced              *slice = (DM_Sliced*)dm->data;
51 
52   PetscFunctionBegin;
53   bs = slice->bs;
54   ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr);
55   ierr = MatSetSizes(*J,slice->n*bs,slice->n*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
56   ierr = MatSetType(*J,mtype);CHKERRQ(ierr);
57   ierr = MatSeqBAIJSetPreallocation(*J,bs,slice->d_nz,slice->d_nnz);CHKERRQ(ierr);
58   ierr = MatMPIBAIJSetPreallocation(*J,bs,slice->d_nz,slice->d_nnz,slice->o_nz,slice->o_nnz);CHKERRQ(ierr);
59   /* In general, we have to do extra work to preallocate for scalar (AIJ) matrices so we check whether it will do any
60   * good before going on with it. */
61   ierr = PetscObjectQueryFunction((PetscObject)*J,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
62   if (!aij) {
63     ierr = PetscObjectQueryFunction((PetscObject)*J,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
64   }
65   if (aij) {
66     if (bs == 1) {
67       ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz);CHKERRQ(ierr);
68       ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz,slice->o_nz,slice->o_nnz);CHKERRQ(ierr);
69     } else if (!slice->d_nnz) {
70       ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz*bs,PETSC_NULL);CHKERRQ(ierr);
71       ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz*bs,PETSC_NULL,slice->o_nz*bs,PETSC_NULL);CHKERRQ(ierr);
72     } else {
73       /* The user has provided preallocation per block-row, convert it to per scalar-row respecting SlicedSetBlockFills() if applicable */
74       ierr = PetscMalloc2(slice->n*bs,PetscInt,&sd_nnz,(!!slice->o_nnz)*slice->n*bs,PetscInt,&so_nnz);CHKERRQ(ierr);
75       for (i=0; i<slice->n*bs; i++) {
76         sd_nnz[i] = (slice->d_nnz[i/bs]-1) * (slice->ofill ? slice->ofill->i[i%bs+1]-slice->ofill->i[i%bs] : bs)
77                                            + (slice->dfill ? slice->dfill->i[i%bs+1]-slice->dfill->i[i%bs] : bs);
78         if (so_nnz) {
79           so_nnz[i] = slice->o_nnz[i/bs] * (slice->ofill ? slice->ofill->i[i%bs+1]-slice->ofill->i[i%bs] : bs);
80         }
81       }
82       ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz*bs,sd_nnz);CHKERRQ(ierr);
83       ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz*bs,sd_nnz,slice->o_nz*bs,so_nnz);CHKERRQ(ierr);
84       ierr = PetscFree2(sd_nnz,so_nnz);CHKERRQ(ierr);
85     }
86   }
87 
88   ierr = MatSetBlockSize(*J,bs);CHKERRQ(ierr);
89 
90   /* Set up the local to global map.  For the scalar map, we have to translate to entry-wise indexing instead of block-wise. */
91   ierr = PetscMalloc((slice->n+slice->Nghosts)*bs*sizeof(PetscInt),&globals);CHKERRQ(ierr);
92   ierr = MatGetOwnershipRange(*J,&rstart,PETSC_NULL);CHKERRQ(ierr);
93   for (i=0; i<slice->n*bs; i++) {
94     globals[i] = rstart + i;
95   }
96   for (i=0; i<slice->Nghosts*bs; i++) {
97     globals[slice->n*bs+i] = slice->ghosts[i/bs]*bs + i%bs;
98   }
99   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,(slice->n+slice->Nghosts)*bs,globals,PETSC_OWN_POINTER,&lmap);CHKERRQ(ierr);
100   ierr = ISLocalToGlobalMappingBlock(lmap,bs,&blmap);CHKERRQ(ierr);
101   ierr = MatSetLocalToGlobalMapping(*J,lmap);CHKERRQ(ierr);
102   ierr = MatSetLocalToGlobalMappingBlock(*J,blmap);CHKERRQ(ierr);
103   ierr = ISLocalToGlobalMappingDestroy(lmap);CHKERRQ(ierr);
104   ierr = ISLocalToGlobalMappingDestroy(blmap);CHKERRQ(ierr);
105   PetscFunctionReturn(0);
106 }
107 
108 #undef __FUNCT__
109 #define __FUNCT__ "SlicedSetGhosts"
110 /*@C
111     SlicedSetGhosts - Sets the global indices of other processes elements that will
112       be ghosts on this process
113 
114     Not Collective
115 
116     Input Parameters:
117 +    slice - the Sliced object
118 .    bs - block size
119 .    nlocal - number of local (owned, non-ghost) blocks
120 .    Nghosts - number of ghost blocks on this process
121 -    ghosts - global indices of each ghost block
122 
123     Level: advanced
124 
125 .seealso SlicedDestroy(), SlicedCreateGlobalVector(), SlicedGetGlobalIndices()
126 
127 @*/
128 PetscErrorCode PETSCDM_DLLEXPORT SlicedSetGhosts(DM dm,PetscInt bs,PetscInt nlocal,PetscInt Nghosts,const PetscInt ghosts[])
129 {
130   PetscErrorCode ierr;
131   DM_Sliced      *slice = (DM_Sliced*)dm->data;
132 
133   PetscFunctionBegin;
134   PetscValidHeaderSpecific(slice,DM_CLASSID,1);
135   ierr = PetscFree(slice->ghosts);CHKERRQ(ierr);
136   ierr = PetscMalloc(Nghosts*sizeof(PetscInt),&slice->ghosts);CHKERRQ(ierr);
137   ierr = PetscMemcpy(slice->ghosts,ghosts,Nghosts*sizeof(PetscInt));CHKERRQ(ierr);
138   slice->bs      = bs;
139   slice->n       = nlocal;
140   slice->Nghosts = Nghosts;
141   PetscFunctionReturn(0);
142 }
143 
144 #undef __FUNCT__
145 #define __FUNCT__ "SlicedSetPreallocation"
146 /*@C
147     SlicedSetPreallocation - sets the matrix memory preallocation for matrices computed by Sliced
148 
149     Not Collective
150 
151     Input Parameters:
152 +    slice - the Sliced object
153 .    d_nz  - number of block nonzeros per block row in diagonal portion of local
154            submatrix  (same for all local rows)
155 .    d_nnz - array containing the number of block nonzeros in the various block rows
156            of the in diagonal portion of the local (possibly different for each block
157            row) or PETSC_NULL.  You must leave room for the diagonal entry even if it is zero.
158 .    o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
159            submatrix (same for all local rows).
160 -    o_nnz - array containing the number of nonzeros in the various block rows of the
161            off-diagonal portion of the local submatrix (possibly different for
162            each block row) or PETSC_NULL.
163 
164     Notes:
165     See MatMPIBAIJSetPreallocation() for more details on preallocation.  If a scalar matrix (AIJ) is
166     obtained with SlicedGetMatrix(), the correct preallocation will be set, respecting SlicedSetBlockFills().
167 
168     Level: advanced
169 
170 .seealso SlicedDestroy(), SlicedCreateGlobalVector(), SlicedGetGlobalIndices(), MatMPIAIJSetPreallocation(),
171          MatMPIBAIJSetPreallocation(), SlicedGetMatrix(), SlicedSetBlockFills()
172 
173 @*/
174 PetscErrorCode PETSCDM_DLLEXPORT SlicedSetPreallocation(DM dm,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
175 {
176   DM_Sliced *slice = (DM_Sliced*)dm->data;
177 
178   PetscFunctionBegin;
179   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
180   slice->d_nz  = d_nz;
181   slice->d_nnz = (PetscInt*)d_nnz;
182   slice->o_nz  = o_nz;
183   slice->o_nnz = (PetscInt*)o_nnz;
184   PetscFunctionReturn(0);
185 }
186 
187 #undef __FUNCT__
188 #define __FUNCT__ "SlicedSetBlockFills_Private"
189 static PetscErrorCode SlicedSetBlockFills_Private(PetscInt bs,const PetscInt *fill,SlicedBlockFills **inf)
190 {
191   PetscErrorCode   ierr;
192   PetscInt         i,j,nz,*fi,*fj;
193   SlicedBlockFills *f;
194 
195   PetscFunctionBegin;
196   PetscValidPointer(inf,3);
197   if (*inf) {ierr = PetscFree3((*inf)->i,(*inf)->j,*inf);CHKERRQ(ierr);}
198   if (!fill) PetscFunctionReturn(0);
199   for (i=0,nz=0; i<bs*bs; i++) if (fill[i]) nz++;
200   ierr = PetscMalloc3(1,SlicedBlockFills,&f,bs+1,PetscInt,&fi,nz,PetscInt,&fj);CHKERRQ(ierr);
201   f->bs = bs;
202   f->nz = nz;
203   f->i  = fi;
204   f->j  = fj;
205   for (i=0,nz=0; i<bs; i++) {
206     fi[i] = nz;
207     for (j=0; j<bs; j++) if (fill[i*bs+j]) fj[nz++] = j;
208   }
209   fi[i] = nz;
210   *inf = f;
211   PetscFunctionReturn(0);
212 }
213 
214 #undef __FUNCT__
215 #define __FUNCT__ "SlicedSetBlockFills"
216 /*@C
217     SlicedSetBlockFills - Sets the fill pattern in each block for a multi-component problem
218     of the matrix returned by SlicedGetMatrix().
219 
220     Logically Collective on Sliced
221 
222     Input Parameter:
223 +   sliced - the Sliced object
224 .   dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block)
225 -   ofill - the fill pattern in the off-diagonal blocks
226 
227     Notes:
228     This only makes sense for multicomponent problems using scalar matrix formats (AIJ).
229     See DASetBlockFills() for example usage.
230 
231     Level: advanced
232 
233 .seealso SlicedGetMatrix(), DASetBlockFills()
234 @*/
235 PetscErrorCode PETSCDM_DLLEXPORT SlicedSetBlockFills(DM dm,const PetscInt *dfill,const PetscInt *ofill)
236 {
237   DM_Sliced      *slice = (DM_Sliced*)dm->data;
238   PetscErrorCode ierr;
239 
240   PetscFunctionBegin;
241   PetscValidHeaderSpecific(slice,DM_CLASSID,1);
242   ierr = SlicedSetBlockFills_Private(slice->bs,dfill,&slice->dfill);CHKERRQ(ierr);
243   ierr = SlicedSetBlockFills_Private(slice->bs,ofill,&slice->ofill);CHKERRQ(ierr);
244   PetscFunctionReturn(0);
245 }
246 
247 #undef __FUNCT__
248 #define __FUNCT__ "SlicedCreate"
249 /*@C
250     SlicedCreate - Creates a DM object, used to manage data for a unstructured problem
251 
252     Collective on MPI_Comm
253 
254     Input Parameter:
255 .   comm - the processors that will share the global vector
256 
257     Output Parameters:
258 .   slice - the slice object
259 
260     Level: advanced
261 
262 .seealso SlicedDestroy(), SlicedCreateGlobalVector(), SlicedGetGlobalIndices()
263 
264 @*/
265 PetscErrorCode PETSCDM_DLLEXPORT SlicedCreate(MPI_Comm comm,DM *dm)
266 {
267   PetscErrorCode ierr;
268   DM             p;
269   DM_Sliced      *slice;
270 
271   PetscFunctionBegin;
272   PetscValidPointer(slice,2);
273 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
274   ierr = DMInitializePackage(PETSC_NULL);CHKERRQ(ierr);
275 #endif
276 
277   ierr = PetscHeaderCreate(p,_p_DM,struct _DMOps,DM_CLASSID,0,"DM",comm,SlicedDestroy,0);CHKERRQ(ierr);
278   ierr = PetscNewLog(p,DM_Sliced,&slice);CHKERRQ(ierr);
279   p->data = slice;
280 
281   ierr = PetscObjectChangeTypeName((PetscObject)p,"Sliced");CHKERRQ(ierr);
282   p->ops->createglobalvector = SlicedCreateGlobalVector;
283   p->ops->getmatrix          = SlicedGetMatrix;
284   p->ops->destroy            = SlicedDestroy;
285   *dm = p;
286   PetscFunctionReturn(0);
287 }
288 
289 extern PetscErrorCode DMDestroy_Private(DM,PetscBool *);
290 
291 #undef __FUNCT__
292 #define __FUNCT__ "SlicedDestroy"
293 /*@C
294     SlicedDestroy - Destroys a vector slice.
295 
296     Collective on Sliced
297 
298     Input Parameter:
299 .   slice - the slice object
300 
301     Level: advanced
302 
303 .seealso SlicedCreate(), SlicedCreateGlobalVector(), SlicedGetGlobalIndices()
304 
305 @*/
306 PetscErrorCode PETSCDM_DLLEXPORT SlicedDestroy(DM dm)
307 {
308   PetscErrorCode ierr;
309   PetscBool      done;
310   DM_Sliced      *slice = (DM_Sliced*)dm->data;
311 
312   PetscFunctionBegin;
313   ierr = DMDestroy_Private(dm,&done);CHKERRQ(ierr);
314   if (!done) PetscFunctionReturn(0);
315 
316   if (slice->globalvector) {ierr = VecDestroy(slice->globalvector);CHKERRQ(ierr);}
317   ierr = PetscFree(slice->ghosts);CHKERRQ(ierr);
318   if (slice->dfill) {ierr = PetscFree3(slice->dfill,slice->dfill->i,slice->dfill->j);CHKERRQ(ierr);}
319   if (slice->ofill) {ierr = PetscFree3(slice->ofill,slice->ofill->i,slice->ofill->j);CHKERRQ(ierr);}
320   ierr = PetscFree(dm->data);CHKERRQ(ierr);
321   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
322   PetscFunctionReturn(0);
323 }
324 
325 
326 #undef __FUNCT__
327 #define __FUNCT__ "SlicedCreateGlobalVector"
328 /*@C
329     SlicedCreateGlobalVector - Creates a vector of the correct size to be gathered into
330         by the slice.
331 
332     Collective on Sliced
333 
334     Input Parameter:
335 .    slice - the slice object
336 
337     Output Parameters:
338 .   gvec - the global vector
339 
340     Level: advanced
341 
342     Notes: Once this has been created you cannot add additional arrays or vectors to be packed.
343 
344 .seealso SlicedDestroy(), SlicedCreate(), SlicedGetGlobalIndices()
345 
346 @*/
347 PetscErrorCode PETSCDM_DLLEXPORT SlicedCreateGlobalVector(DM dm,Vec *gvec)
348 {
349   PetscErrorCode     ierr;
350   PetscInt           bs,cnt;
351   DM_Sliced          *slice = (DM_Sliced*)dm->data;
352 
353   PetscFunctionBegin;
354   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
355   PetscValidPointer(gvec,2);
356   *gvec = 0;
357   if (slice->globalvector) {
358     ierr = PetscObjectGetReference((PetscObject)slice->globalvector,&cnt);CHKERRQ(ierr);
359     if (cnt == 1) {             /* Nobody else has a reference so we can just reference it and give it away */
360       *gvec = slice->globalvector;
361       ierr = PetscObjectReference((PetscObject)*gvec);CHKERRQ(ierr);
362       ierr = VecZeroEntries(*gvec);CHKERRQ(ierr);
363     } else {                    /* Someone else has a reference so we duplicate the global vector */
364       ierr = VecDuplicate(slice->globalvector,gvec);CHKERRQ(ierr);
365     }
366   } else {
367     bs = slice->bs;
368     ierr = VecCreateGhostBlock(((PetscObject)dm)->comm,bs,slice->n*bs,PETSC_DETERMINE,slice->Nghosts,slice->ghosts,&slice->globalvector);CHKERRQ(ierr);
369     *gvec = slice->globalvector;
370     ierr = PetscObjectReference((PetscObject)*gvec);CHKERRQ(ierr);
371   }
372   PetscFunctionReturn(0);
373 }
374 
375 #undef __FUNCT__
376 #define __FUNCT__ "SlicedGetGlobalIndices"
377 /*@C
378     SlicedGetGlobalIndices - Gets the global indices for all the local entries
379 
380     Collective on Sliced
381 
382     Input Parameter:
383 .    slice - the slice object
384 
385     Output Parameters:
386 .    idx - the individual indices for each packed vector/array
387 
388     Level: advanced
389 
390     Notes:
391        The idx parameters should be freed by the calling routine with PetscFree()
392 
393 .seealso SlicedDestroy(), SlicedCreateGlobalVector(), SlicedCreate()
394 
395 @*/
396 PetscErrorCode PETSCDM_DLLEXPORT SlicedGetGlobalIndices(DM dm,PetscInt *idx[])
397 {
398   PetscFunctionReturn(0);
399 }
400 
401 
402 /* Explanation of the missing functions for DA-style handling of the local vector:
403 
404    SlicedCreateLocalVector()
405    SlicedGlobalToLocalBegin()
406    SlicedGlobalToLocalEnd()
407 
408  There is no way to get the global form from a local form, so SlicedCreateLocalVector() is a memory leak without
409  external accounting for the global vector.  Also, Sliced intends the user to work with the VecGhost interface since the
410  ghosts are already ordered after the owned entries.  Contrast this to a DA where the local vector has a special
411  ordering described by the structured grid, hence it cannot share memory with the global form.  For this reason, users
412  of Sliced should work with the global vector and use
413 
414    VecGhostGetLocalForm(), VecGhostRestoreLocalForm()
415    VecGhostUpdateBegin(), VecGhostUpdateEnd()
416 
417  rather than the missing DA-style functions.  This is conceptually simpler and offers better performance than is
418  possible with the DA-style interface.
419 */
420