Actual source code: sliced.c
1: #define PETSCDM_DLL
2:
3: #include petscda.h
4: #include petscmat.h
5: #include private/dmimpl.h
8: typedef struct _SlicedOps *SlicedOps;
9: struct _SlicedOps {
10: DMOPS(Sliced)
11: };
13: struct _p_Sliced {
14: PETSCHEADER(struct _SlicedOps);
15: DMHEADER
16: Vec globalvector;
17: PetscInt bs,n,N,Nghosts,*ghosts;
18: PetscInt d_nz,o_nz,*d_nnz,*o_nnz;
19: };
23: /*@C
24: SlicedGetMatrix - Creates a matrix with the correct parallel layout required for
25: computing the Jacobian on a function defined using the informatin in Sliced.
27: Collective on Sliced
29: Input Parameter:
30: + slice - the slice object
31: - mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ,
32: or any type which inherits from one of these (such as MATAIJ, MATLUSOL, etc.).
34: Output Parameters:
35: . J - matrix with the correct nonzero preallocation
36: (obviously without the correct Jacobian values)
38: Level: advanced
40: Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
41: do not need to do it yourself.
43: .seealso ISColoringView(), ISColoringGetIS(), MatFDColoringCreate(), DASetBlockFills()
45: @*/
46: PetscErrorCode SlicedGetMatrix(Sliced slice, const MatType mtype,Mat *J)
47: {
48: PetscErrorCode ierr;
49: PetscInt *globals,rstart,i;
50: ISLocalToGlobalMapping lmap;
53: MatCreate(((PetscObject)slice)->comm,J);
54: MatSetSizes(*J,slice->n,slice->n,PETSC_DETERMINE,PETSC_DETERMINE);
55: MatSetType(*J,mtype);
56: MatSetBlockSize(*J,slice->bs);
57: MatSeqAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz);
58: MatMPIAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz,slice->o_nz,slice->o_nnz);
59: MatSeqBAIJSetPreallocation(*J,slice->bs,slice->d_nz,slice->d_nnz);
60: MatMPIBAIJSetPreallocation(*J,slice->bs,slice->d_nz,slice->d_nnz,slice->o_nz,slice->o_nnz);
62: PetscMalloc((slice->n+slice->Nghosts+1)*sizeof(PetscInt),&globals);
63: MatGetOwnershipRange(*J,&rstart,PETSC_NULL);
64: for (i=0; i<slice->n; i++) {
65: globals[i] = rstart + i;
66: }
67: PetscMemcpy(globals+slice->n,slice->ghosts,slice->Nghosts*sizeof(PetscInt));
68: ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,slice->n+slice->Nghosts,globals,&lmap);
69: PetscFree(globals);
70: MatSetLocalToGlobalMapping(*J,lmap);
71: ISLocalToGlobalMappingDestroy(lmap);
72: return(0);
73: }
77: /*@C
78: SlicedSetGhosts - Sets the global indices of other processes elements that will
79: be ghosts on this process
81: Not Collective
83: Input Parameters:
84: + slice - the Sliced object
85: . bs - block size
86: . nlocal - number of local (non-ghost) entries
87: . Nghosts - number of ghosts on this process
88: - ghosts - indices of all the ghost points
90: Level: advanced
92: .seealso SlicedDestroy(), SlicedCreateGlobalVector(), SlicedGetGlobalIndices()
94: @*/
95: PetscErrorCode SlicedSetGhosts(Sliced slice,PetscInt bs,PetscInt nlocal,PetscInt Nghosts,const PetscInt ghosts[])
96: {
101: PetscFree(slice->ghosts);
102: PetscMalloc((1+Nghosts)*sizeof(PetscInt),&slice->ghosts);
103: PetscMemcpy(slice->ghosts,ghosts,Nghosts*sizeof(PetscInt));
104: slice->bs = bs;
105: slice->n = nlocal;
106: slice->Nghosts = Nghosts;
107: return(0);
108: }
112: /*@C
113: SlicedSetPreallocation - sets the matrix memory preallocation for matrices computed by Sliced
115: Not Collective
117: Input Parameters:
118: + slice - the Sliced object
119: . d_nz - maximum number of nonzeros in any row of diagonal block
120: . d_nnz - number of nonzeros in each row of diagonal block
121: . o_nz - maximum number of nonzeros in any row of off-diagonal block
122: . o_nnz - number of nonzeros in each row of off-diagonal block
125: Level: advanced
127: .seealso SlicedDestroy(), SlicedCreateGlobalVector(), SlicedGetGlobalIndices(), MatMPIAIJSetPreallocation(),
128: MatMPIBAIJSetPreallocation()
130: @*/
131: PetscErrorCode SlicedSetPreallocation(Sliced slice,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
132: {
135: slice->d_nz = d_nz;
136: slice->d_nnz = (PetscInt*)d_nnz;
137: slice->o_nz = o_nz;
138: slice->o_nnz = (PetscInt*)o_nnz;
139: return(0);
140: }
144: /*@C
145: SlicedCreate - Creates a DM object, used to manage data for a unstructured problem
147: Collective on MPI_Comm
149: Input Parameter:
150: . comm - the processors that will share the global vector
152: Output Parameters:
153: . slice - the slice object
155: Level: advanced
157: .seealso SlicedDestroy(), SlicedCreateGlobalVector(), SlicedGetGlobalIndices()
159: @*/
160: PetscErrorCode SlicedCreate(MPI_Comm comm,Sliced *slice)
161: {
163: Sliced p;
167: *slice = PETSC_NULL;
168: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
169: DMInitializePackage(PETSC_NULL);
170: #endif
172: PetscHeaderCreate(p,_p_Sliced,struct _SlicedOps,DM_COOKIE,0,"DM",comm,SlicedDestroy,0);
173: PetscObjectChangeTypeName((PetscObject)p,"Sliced");
174: p->ops->createglobalvector = SlicedCreateGlobalVector;
175: p->ops->createlocalvector = SlicedCreateLocalVector;
176: p->ops->globaltolocalbegin = SlicedGlobalToLocalBegin;
177: p->ops->globaltolocalend = SlicedGlobalToLocalEnd;
178: p->ops->getmatrix = SlicedGetMatrix;
179: p->ops->destroy = SlicedDestroy;
180: *slice = p;
181: return(0);
182: }
188: /*@C
189: SlicedDestroy - Destroys a vector slice.
191: Collective on Sliced
193: Input Parameter:
194: . slice - the slice object
196: Level: advanced
198: .seealso SlicedCreate(), SlicedCreateGlobalVector(), SlicedGetGlobalIndices()
200: @*/
201: PetscErrorCode SlicedDestroy(Sliced slice)
202: {
204: PetscTruth done;
207: DMDestroy_Private((DM)slice,&done);
208: if (!done) return(0);
210: if (slice->globalvector) {VecDestroy(slice->globalvector);}
211: PetscFree(slice->ghosts);
212: PetscHeaderDestroy(slice);
213: return(0);
214: }
219: /*@C
220: SlicedCreateGlobalVector - Creates a vector of the correct size to be gathered into
221: by the slice.
223: Collective on Sliced
225: Input Parameter:
226: . slice - the slice object
228: Output Parameters:
229: . gvec - the global vector
231: Level: advanced
233: Notes: Once this has been created you cannot add additional arrays or vectors to be packed.
235: .seealso SlicedDestroy(), SlicedCreate(), SlicedGetGlobalIndices()
237: @*/
238: PetscErrorCode SlicedCreateGlobalVector(Sliced slice,Vec *gvec)
239: {
240: PetscErrorCode ierr;
244: if (slice->globalvector) {
245: VecDuplicate(slice->globalvector,gvec);
246: } else {
247: VecCreateGhostBlock(((PetscObject)slice)->comm,slice->bs,slice->n,PETSC_DETERMINE,slice->Nghosts,slice->ghosts,&slice->globalvector);
248: *gvec = slice->globalvector;
249: PetscObjectReference((PetscObject)*gvec);
250: }
251: return(0);
252: }
256: /*@C
257: SlicedCreateLocalVector - Creates a vector of the correct size to be gatherer from
258: by the slice.
260: Collective on Sliced
262: Input Parameter:
263: . slice - the slice object
265: Output Parameters:
266: . gvec - the global vector
268: Level: advanced
270: Notes: Once this has been created you cannot add additional arrays or vectors to be packed.
272: .seealso SlicedDestroy(), SlicedCreate(), SlicedGetGlobalIndices(), SlicedCreateGlobalVector()
274: @*/
275: PetscErrorCode SlicedCreateLocalVector(Sliced slice,Vec *lvec)
276: {
278: Vec gvec;
281: if (slice->globalvector) {
282: VecDuplicate(slice->globalvector,&gvec);
283: } else {
284: VecCreateGhostBlock(((PetscObject)slice)->comm,slice->bs,slice->n,PETSC_DETERMINE,slice->Nghosts,slice->ghosts,&slice->globalvector);
285: PetscObjectReference((PetscObject)slice->globalvector);
286: gvec = slice->globalvector;
287: }
288: VecGhostGetLocalForm(gvec,lvec);
289: return(0);
290: }
294: /*@C
295: SlicedGetGlobalIndices - Gets the global indices for all the local entries
297: Collective on Sliced
299: Input Parameter:
300: . slice - the slice object
302: Output Parameters:
303: . idx - the individual indices for each packed vector/array
304:
305: Level: advanced
307: Notes:
308: The idx parameters should be freed by the calling routine with PetscFree()
310: .seealso SlicedDestroy(), SlicedCreateGlobalVector(), SlicedCreate()
312: @*/
313: PetscErrorCode SlicedGetGlobalIndices(Sliced slice,PetscInt *idx[])
314: {
315: return(0);
316: }
320: /*@C
321: SlicedGlobalToLocalBegin - Begins the communication from a global sliced vector to a local one
323: Collective on DA
325: Input Parameters:
326: + sliced - the sliced context
327: . g - the global vector
328: - mode - one of INSERT_VALUES or ADD_VALUES
330: Output Parameter:
331: . l - the local values
333: Level: beginner
335: .keywords: distributed array, global to local, begin
337: .seealso: SlicedCreate(), SlicedGlobalToLocalEnd()
338:
340: @*/
341: PetscErrorCode SlicedGlobalToLocalBegin(Sliced sliced,Vec g,InsertMode mode,Vec l)
342: {
344: Vec lform;
350: /* only works if local vector l is shared with global vector */
351: VecGhostGetLocalForm(g,&lform);
352: VecGhostRestoreLocalForm(g,&lform);
353: if (lform != l) SETERRQ(PETSC_ERR_ARG_INCOMP,"Local vector must be local form of global vector (see VecGhostUpdate())");
354: VecGhostUpdateBegin(g,mode,SCATTER_FORWARD);
355: return(0);
356: }
360: /*@C
361: SlicedGlobalToLocalEnd - Ends the communication from a global sliced vector to a local one
363: Collective on DA
365: Input Parameters:
366: + sliced - the sliced context
367: . g - the global vector
368: - mode - one of INSERT_VALUES or ADD_VALUES
370: Output Parameter:
371: . l - the local values
373: Level: beginner
375: .keywords: distributed array, global to local, begin
377: .seealso: SlicedCreate(), SlicedGlobalToLocalEnd()
378:
380: @*/
381: PetscErrorCode SlicedGlobalToLocalEnd(Sliced sliced,Vec g,InsertMode mode,Vec l)
382: {
384: Vec lform;
390: /* only works if local vector l is shared with global vector */
391: VecGhostGetLocalForm(g,&lform);
392: VecGhostRestoreLocalForm(g,&lform);
393: if (lform != l) SETERRQ(PETSC_ERR_ARG_INCOMP,"Local vector must be local form of global vector (see VecGhostUpdate())");
394: VecGhostUpdateEnd(g,mode,SCATTER_FORWARD);
395: return(0);
396: }