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