xref: /petsc/src/dm/impls/da/dalocal.c (revision 2981ebdba252c1ab1bebf9afd5693bb0df49865f)
1 
2 /*
3   Code for manipulating distributed regular arrays in parallel.
4 */
5 
6 #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/
7 
8 /*
9    This allows the DMDA vectors to properly tell MATLAB their dimensions
10 */
11 #if defined(PETSC_HAVE_MATLAB_ENGINE)
12 #include <engine.h>   /* MATLAB include file */
13 #include <mex.h>      /* MATLAB include file */
14 EXTERN_C_BEGIN
15 #undef __FUNCT__
16 #define __FUNCT__ "VecMatlabEnginePut_DA2d"
17 PetscErrorCode  VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine)
18 {
19   PetscErrorCode ierr;
20   PetscInt       n,m;
21   Vec            vec = (Vec)obj;
22   PetscScalar    *array;
23   mxArray        *mat;
24   DM             da;
25 
26   PetscFunctionBegin;
27   ierr = PetscObjectQuery((PetscObject)vec,"DM",(PetscObject*)&da);CHKERRQ(ierr);
28   if (!da) SETERRQ(((PetscObject)vec)->comm,PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA");
29   ierr = DMDAGetGhostCorners(da,0,0,0,&m,&n,0);CHKERRQ(ierr);
30 
31   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
32 #if !defined(PETSC_USE_COMPLEX)
33   mat  = mxCreateDoubleMatrix(m,n,mxREAL);
34 #else
35   mat  = mxCreateDoubleMatrix(m,n,mxCOMPLEX);
36 #endif
37   ierr = PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));CHKERRQ(ierr);
38   ierr = PetscObjectName(obj);CHKERRQ(ierr);
39   engPutVariable((Engine *)mengine,obj->name,mat);
40 
41   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
42   PetscFunctionReturn(0);
43 }
44 EXTERN_C_END
45 #endif
46 
47 
48 #undef __FUNCT__
49 #define __FUNCT__ "DMCreateLocalVector_DA"
50 PetscErrorCode  DMCreateLocalVector_DA(DM da,Vec* g)
51 {
52   PetscErrorCode ierr;
53   DM_DA          *dd = (DM_DA*)da->data;
54 
55   PetscFunctionBegin;
56   PetscValidHeaderSpecific(da,DM_CLASSID,1);
57   PetscValidPointer(g,2);
58   ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr);
59   ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr);
60   ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr);
61   ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr);
62   ierr = PetscObjectCompose((PetscObject)*g,"DM",(PetscObject)da);CHKERRQ(ierr);
63 #if defined(PETSC_HAVE_MATLAB_ENGINE)
64   if (dd->w == 1  && dd->dim == 2) {
65     ierr = PetscObjectComposeFunctionDynamic((PetscObject)*g,"PetscMatlabEnginePut_C","VecMatlabEnginePut_DA2d",VecMatlabEnginePut_DA2d);CHKERRQ(ierr);
66   }
67 #endif
68   PetscFunctionReturn(0);
69 }
70 
71 #undef __FUNCT__
72 #define __FUNCT__ "DMDACreateSection"
73 /*@C
74   DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with
75   different numbers of dofs on vertices, cells, and faces in each direction.
76 
77   Input Parameters:
78 + dm- The DMDA
79 . numFields - The number of fields
80 . numComp - The number of components in each field, or PETSC_NULL for 1
81 . numVertexDof - The number of dofs per vertex for each field, or PETSC_NULL
82 . numFaceDof - The number of dofs per face for each field and direction, or PETSC_NULL
83 - numCellDof - The number of dofs per cell for each field, or PETSC_NULL
84 
85   Level: developer
86 
87   Note:
88   The default DMDA numbering is as follows:
89 
90     - Cells:    [0,             nC)
91     - Vertices: [nC,            nC+nV)
92     - X-Faces:  [nC+nV,         nC+nV+nXF)
93     - Y-Faces:  [nC+nV+nXF,     nC+nV+nXF+nYF)
94     - Z-Faces:  [nC+nV+nXF+nYF, nC+nV+nXF+nYF+nZF)
95 
96   We interpret the default DMDA partition as a cell partition, and the data assignment as a cell assignment.
97 @*/
98 PetscErrorCode DMDACreateSection(DM dm, PetscInt numFields, PetscInt numComp[], PetscInt numVertexDof[], PetscInt numFaceDof[], PetscInt numCellDof[])
99 {
100   DM_DA         *da  = (DM_DA *) dm->data;
101   PetscInt       mx  = da->Xe - da->Xs, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
102   PetscInt       nC  = (mx  )*(da->dim > 1 ? (my  )*(da->dim > 2 ? (mz  ) : 1) : 1);
103   PetscInt       nV  = (mx+1)*(da->dim > 1 ? (my+1)*(da->dim > 2 ? (mz+1) : 1) : 1);
104   PetscInt       nXF = (mx+1)*(da->dim > 1 ? (my  )*(da->dim > 2 ? (mz  ) : 1) : 1);
105   PetscInt       nYF = (mx  )*(da->dim > 1 ? (my+1)*(da->dim > 2 ? (mz  ) : 0) : 1);
106   PetscInt       nZF = (mx  )*(da->dim > 1 ? (my  )*(da->dim > 2 ? (mz+1) : 0) : 0);
107   PetscInt       cStart  = 0,     cEnd  = cStart+nC;
108   PetscInt       vStart  = cEnd,  vEnd  = vStart+nV;
109   PetscInt       xfStart = vEnd,  xfEnd = xfStart+nXF;
110   PetscInt       yfStart = xfEnd, yfEnd = yfStart+nYF;
111   PetscInt       zfStart = yfEnd, zfEnd = zfStart+nZF;
112   PetscInt       pStart  = 0,     pEnd  = zfEnd;
113   PetscInt       numVertexTotDof = 0, numCellTotDof = 0, numFaceTotDof[3] = {0, 0, 0};
114   PetscInt       f, v, c, xf, yf, zf;
115   PetscErrorCode ierr;
116 
117   PetscFunctionBegin;
118   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
119   for(f = 0; f < numFields; ++f) {
120     if (numVertexDof) {numVertexTotDof  += numVertexDof[f];}
121     if (numCellDof)   {numCellTotDof    += numCellDof[f];}
122     if (numFaceDof)   {numFaceTotDof[0] += numFaceDof[f*da->dim+0];
123                        numFaceTotDof[1] += da->dim > 1 ? numFaceDof[f*da->dim+1] : 0;
124                        numFaceTotDof[2] += da->dim > 2 ? numFaceDof[f*da->dim+2] : 0;}
125   }
126   ierr = PetscSectionCreate(((PetscObject) dm)->comm, &da->defaultSection);CHKERRQ(ierr);
127   if (numFields > 1) {
128     ierr = PetscSectionSetNumFields(da->defaultSection, numFields);CHKERRQ(ierr);
129     if (numComp) {
130       for(f = 0; f < numFields; ++f) {
131         ierr = PetscSectionSetFieldComponents(da->defaultSection, f, numComp[f]);CHKERRQ(ierr);
132       }
133     }
134   } else {
135     numFields = 0;
136   }
137   ierr = PetscSectionSetChart(da->defaultSection, pStart, pEnd);CHKERRQ(ierr);
138   if (numVertexDof) {
139     for(v = vStart; v < vEnd; ++v) {
140       for(f = 0; f < numFields; ++f) {
141         ierr = PetscSectionSetFieldDof(da->defaultSection, v, f, numVertexDof[f]);CHKERRQ(ierr);
142       }
143       ierr = PetscSectionSetDof(da->defaultSection, v, numVertexTotDof);CHKERRQ(ierr);
144     }
145   }
146   if (numFaceDof) {
147     for(xf = xfStart; xf < xfEnd; ++xf) {
148       for(f = 0; f < numFields; ++f) {
149         ierr = PetscSectionSetFieldDof(da->defaultSection, xf, f, numFaceDof[f*da->dim+0]);CHKERRQ(ierr);
150       }
151       ierr = PetscSectionSetDof(da->defaultSection, xf, numFaceTotDof[0]);CHKERRQ(ierr);
152     }
153     for(yf = yfStart; yf < yfEnd; ++yf) {
154       for(f = 0; f < numFields; ++f) {
155         ierr = PetscSectionSetFieldDof(da->defaultSection, yf, f, numFaceDof[f*da->dim+1]);CHKERRQ(ierr);
156       }
157       ierr = PetscSectionSetDof(da->defaultSection, yf, numFaceTotDof[1]);CHKERRQ(ierr);
158     }
159     for(zf = zfStart; zf < zfEnd; ++zf) {
160       for(f = 0; f < numFields; ++f) {
161         ierr = PetscSectionSetFieldDof(da->defaultSection, zf, f, numFaceDof[f*da->dim+2]);CHKERRQ(ierr);
162       }
163       ierr = PetscSectionSetDof(da->defaultSection, zf, numFaceTotDof[2]);CHKERRQ(ierr);
164     }
165   }
166   if (numCellDof) {
167     for(c = cStart; c < cEnd; ++c) {
168       for(f = 0; f < numFields; ++f) {
169         ierr = PetscSectionSetFieldDof(da->defaultSection, c, f, numCellDof[f]);CHKERRQ(ierr);
170       }
171       ierr = PetscSectionSetDof(da->defaultSection, c, numCellTotDof);CHKERRQ(ierr);
172     }
173   }
174   ierr = PetscSectionSetUp(da->defaultSection);CHKERRQ(ierr);
175   da->defaultGlobalSection = da->defaultSection;
176   PetscFunctionReturn(0);
177 }
178 
179 /* ------------------------------------------------------------------- */
180 #if defined(PETSC_HAVE_ADIC)
181 
182 EXTERN_C_BEGIN
183 #include <adic/ad_utils.h>
184 EXTERN_C_END
185 
186 #undef __FUNCT__
187 #define __FUNCT__ "DMDAGetAdicArray"
188 /*@C
189      DMDAGetAdicArray - Gets an array of derivative types for a DMDA
190 
191     Input Parameter:
192 +    da - information about my local patch
193 -    ghosted - do you want arrays for the ghosted or nonghosted patch
194 
195     Output Parameters:
196 +    vptr - array data structured to be passed to ad_FormFunctionLocal()
197 .    array_start - actual start of 1d array of all values that adiC can access directly (may be null)
198 -    tdof - total number of degrees of freedom represented in array_start (may be null)
199 
200      Notes:
201        The vector values are NOT initialized and may have garbage in them, so you may need
202        to zero them.
203 
204        Returns the same type of object as the DMDAVecGetArray() except its elements are
205            derivative types instead of PetscScalars
206 
207      Level: advanced
208 
209 .seealso: DMDARestoreAdicArray()
210 
211 @*/
212 PetscErrorCode  DMDAGetAdicArray(DM da,PetscBool  ghosted,void *vptr,void *array_start,PetscInt *tdof)
213 {
214   PetscErrorCode ierr;
215   PetscInt       j,i,deriv_type_size,xs,ys,xm,ym,zs,zm,itdof;
216   char           *iarray_start;
217   void           **iptr = (void**)vptr;
218   DM_DA          *dd = (DM_DA*)da->data;
219 
220   PetscFunctionBegin;
221   PetscValidHeaderSpecific(da,DM_CLASSID,1);
222   if (ghosted) {
223     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
224       if (dd->adarrayghostedin[i]) {
225         *iptr                   = dd->adarrayghostedin[i];
226         iarray_start            = (char*)dd->adstartghostedin[i];
227         itdof                   = dd->ghostedtdof;
228         dd->adarrayghostedin[i] = PETSC_NULL;
229         dd->adstartghostedin[i] = PETSC_NULL;
230 
231         goto done;
232       }
233     }
234     xs = dd->Xs;
235     ys = dd->Ys;
236     zs = dd->Zs;
237     xm = dd->Xe-dd->Xs;
238     ym = dd->Ye-dd->Ys;
239     zm = dd->Ze-dd->Zs;
240   } else {
241     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
242       if (dd->adarrayin[i]) {
243         *iptr            = dd->adarrayin[i];
244         iarray_start     = (char*)dd->adstartin[i];
245         itdof            = dd->tdof;
246         dd->adarrayin[i] = PETSC_NULL;
247         dd->adstartin[i] = PETSC_NULL;
248 
249         goto done;
250       }
251     }
252     xs = dd->xs;
253     ys = dd->ys;
254     zs = dd->zs;
255     xm = dd->xe-dd->xs;
256     ym = dd->ye-dd->ys;
257     zm = dd->ze-dd->zs;
258   }
259   deriv_type_size = PetscADGetDerivTypeSize();
260 
261   switch (dd->dim) {
262     case 1: {
263       void *ptr;
264       itdof = xm;
265 
266       ierr  = PetscMalloc(xm*deriv_type_size,&iarray_start);CHKERRQ(ierr);
267 
268       ptr   = (void*)(iarray_start - xs*deriv_type_size);
269       *iptr = (void*)ptr;
270       break;}
271     case 2: {
272       void **ptr;
273       itdof = xm*ym;
274 
275       ierr  = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*deriv_type_size,&iarray_start);CHKERRQ(ierr);
276 
277       ptr  = (void**)(iarray_start + xm*ym*deriv_type_size - ys*sizeof(void*));
278       for(j=ys;j<ys+ym;j++) {
279         ptr[j] = iarray_start + deriv_type_size*(xm*(j-ys) - xs);
280       }
281       *iptr = (void*)ptr;
282       break;}
283     case 3: {
284       void ***ptr,**bptr;
285       itdof = xm*ym*zm;
286 
287       ierr  = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*deriv_type_size,&iarray_start);CHKERRQ(ierr);
288 
289       ptr  = (void***)(iarray_start + xm*ym*zm*deriv_type_size - zs*sizeof(void*));
290       bptr = (void**)(iarray_start + xm*ym*zm*deriv_type_size + zm*sizeof(void**));
291 
292       for(i=zs;i<zs+zm;i++) {
293         ptr[i] = bptr + ((i-zs)*ym - ys);
294       }
295       for (i=zs; i<zs+zm; i++) {
296         for (j=ys; j<ys+ym; j++) {
297           ptr[i][j] = iarray_start + deriv_type_size*(xm*ym*(i-zs) + xm*(j-ys) - xs);
298         }
299       }
300 
301       *iptr = (void*)ptr;
302       break;}
303     default:
304       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
305   }
306 
307   done:
308   /* add arrays to the checked out list */
309   if (ghosted) {
310     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
311       if (!dd->adarrayghostedout[i]) {
312         dd->adarrayghostedout[i] = *iptr ;
313         dd->adstartghostedout[i] = iarray_start;
314         dd->ghostedtdof          = itdof;
315         break;
316       }
317     }
318   } else {
319     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
320       if (!dd->adarrayout[i]) {
321         dd->adarrayout[i] = *iptr ;
322         dd->adstartout[i] = iarray_start;
323         dd->tdof          = itdof;
324         break;
325       }
326     }
327   }
328   if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Too many DMDA ADIC arrays obtained");
329   if (tdof)        *tdof = itdof;
330   if (array_start) *(void**)array_start = iarray_start;
331   PetscFunctionReturn(0);
332 }
333 
334 #undef __FUNCT__
335 #define __FUNCT__ "DMDARestoreAdicArray"
336 /*@C
337      DMDARestoreAdicArray - Restores an array of derivative types for a DMDA
338 
339     Input Parameter:
340 +    da - information about my local patch
341 -    ghosted - do you want arrays for the ghosted or nonghosted patch
342 
343     Output Parameters:
344 +    ptr - array data structured to be passed to ad_FormFunctionLocal()
345 .    array_start - actual start of 1d array of all values that adiC can access directly
346 -    tdof - total number of degrees of freedom represented in array_start
347 
348      Level: advanced
349 
350 .seealso: DMDAGetAdicArray()
351 
352 @*/
353 PetscErrorCode  DMDARestoreAdicArray(DM da,PetscBool  ghosted,void *ptr,void *array_start,PetscInt *tdof)
354 {
355   PetscInt  i;
356   void      **iptr = (void**)ptr,iarray_start = 0;
357 
358   PetscFunctionBegin;
359   PetscValidHeaderSpecific(da,DM_CLASSID,1);
360   if (ghosted) {
361     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
362       if (dd->adarrayghostedout[i] == *iptr) {
363         iarray_start             = dd->adstartghostedout[i];
364         dd->adarrayghostedout[i] = PETSC_NULL;
365         dd->adstartghostedout[i] = PETSC_NULL;
366         break;
367       }
368     }
369     if (!iarray_start) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Could not find array in checkout list");
370     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
371       if (!dd->adarrayghostedin[i]){
372         dd->adarrayghostedin[i] = *iptr;
373         dd->adstartghostedin[i] = iarray_start;
374         break;
375       }
376     }
377   } else {
378     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
379       if (dd->adarrayout[i] == *iptr) {
380         iarray_start      = dd->adstartout[i];
381         dd->adarrayout[i] = PETSC_NULL;
382         dd->adstartout[i] = PETSC_NULL;
383         break;
384       }
385     }
386     if (!iarray_start) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Could not find array in checkout list");
387     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
388       if (!dd->adarrayin[i]){
389         dd->adarrayin[i]   = *iptr;
390         dd->adstartin[i]   = iarray_start;
391         break;
392       }
393     }
394   }
395   PetscFunctionReturn(0);
396 }
397 
398 #undef __FUNCT__
399 #define __FUNCT__ "ad_DAGetArray"
400 PetscErrorCode  ad_DAGetArray(DM da,PetscBool  ghosted,void *iptr)
401 {
402   PetscErrorCode ierr;
403   PetscFunctionBegin;
404   ierr = DMDAGetAdicArray(da,ghosted,iptr,0,0);CHKERRQ(ierr);
405   PetscFunctionReturn(0);
406 }
407 
408 #undef __FUNCT__
409 #define __FUNCT__ "ad_DARestoreArray"
410 PetscErrorCode  ad_DARestoreArray(DM da,PetscBool  ghosted,void *iptr)
411 {
412   PetscErrorCode ierr;
413   PetscFunctionBegin;
414   ierr = DMDARestoreAdicArray(da,ghosted,iptr,0,0);CHKERRQ(ierr);
415   PetscFunctionReturn(0);
416 }
417 
418 #endif
419 
420 #undef __FUNCT__
421 #define __FUNCT__ "DMDAGetArray"
422 /*@C
423      DMDAGetArray - Gets a work array for a DMDA
424 
425     Input Parameter:
426 +    da - information about my local patch
427 -    ghosted - do you want arrays for the ghosted or nonghosted patch
428 
429     Output Parameters:
430 .    vptr - array data structured
431 
432     Note:  The vector values are NOT initialized and may have garbage in them, so you may need
433            to zero them.
434 
435   Level: advanced
436 
437 .seealso: DMDARestoreArray(), DMDAGetAdicArray()
438 
439 @*/
440 PetscErrorCode  DMDAGetArray(DM da,PetscBool  ghosted,void *vptr)
441 {
442   PetscErrorCode ierr;
443   PetscInt       j,i,xs,ys,xm,ym,zs,zm;
444   char           *iarray_start;
445   void           **iptr = (void**)vptr;
446   DM_DA          *dd = (DM_DA*)da->data;
447 
448   PetscFunctionBegin;
449   PetscValidHeaderSpecific(da,DM_CLASSID,1);
450   if (ghosted) {
451     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
452       if (dd->arrayghostedin[i]) {
453         *iptr                 = dd->arrayghostedin[i];
454         iarray_start          = (char*)dd->startghostedin[i];
455         dd->arrayghostedin[i] = PETSC_NULL;
456         dd->startghostedin[i] = PETSC_NULL;
457 
458         goto done;
459       }
460     }
461     xs = dd->Xs;
462     ys = dd->Ys;
463     zs = dd->Zs;
464     xm = dd->Xe-dd->Xs;
465     ym = dd->Ye-dd->Ys;
466     zm = dd->Ze-dd->Zs;
467   } else {
468     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
469       if (dd->arrayin[i]) {
470         *iptr          = dd->arrayin[i];
471         iarray_start   = (char*)dd->startin[i];
472         dd->arrayin[i] = PETSC_NULL;
473         dd->startin[i] = PETSC_NULL;
474 
475         goto done;
476       }
477     }
478     xs = dd->xs;
479     ys = dd->ys;
480     zs = dd->zs;
481     xm = dd->xe-dd->xs;
482     ym = dd->ye-dd->ys;
483     zm = dd->ze-dd->zs;
484   }
485 
486   switch (dd->dim) {
487     case 1: {
488       void *ptr;
489 
490       ierr  = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
491 
492       ptr   = (void*)(iarray_start - xs*sizeof(PetscScalar));
493       *iptr = (void*)ptr;
494       break;}
495     case 2: {
496       void **ptr;
497 
498       ierr  = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
499 
500       ptr  = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*));
501       for(j=ys;j<ys+ym;j++) {
502         ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs);
503       }
504       *iptr = (void*)ptr;
505       break;}
506     case 3: {
507       void ***ptr,**bptr;
508 
509       ierr  = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
510 
511       ptr  = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*));
512       bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**));
513       for(i=zs;i<zs+zm;i++) {
514         ptr[i] = bptr + ((i-zs)*ym - ys);
515       }
516       for (i=zs; i<zs+zm; i++) {
517         for (j=ys; j<ys+ym; j++) {
518           ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
519         }
520       }
521 
522       *iptr = (void*)ptr;
523       break;}
524     default:
525       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
526   }
527 
528   done:
529   /* add arrays to the checked out list */
530   if (ghosted) {
531     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
532       if (!dd->arrayghostedout[i]) {
533         dd->arrayghostedout[i] = *iptr ;
534         dd->startghostedout[i] = iarray_start;
535         break;
536       }
537     }
538   } else {
539     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
540       if (!dd->arrayout[i]) {
541         dd->arrayout[i] = *iptr ;
542         dd->startout[i] = iarray_start;
543         break;
544       }
545     }
546   }
547   PetscFunctionReturn(0);
548 }
549 
550 #undef __FUNCT__
551 #define __FUNCT__ "DMDARestoreArray"
552 /*@C
553      DMDARestoreArray - Restores an array of derivative types for a DMDA
554 
555     Input Parameter:
556 +    da - information about my local patch
557 .    ghosted - do you want arrays for the ghosted or nonghosted patch
558 -    vptr - array data structured to be passed to ad_FormFunctionLocal()
559 
560      Level: advanced
561 
562 .seealso: DMDAGetArray(), DMDAGetAdicArray()
563 
564 @*/
565 PetscErrorCode  DMDARestoreArray(DM da,PetscBool  ghosted,void *vptr)
566 {
567   PetscInt  i;
568   void      **iptr = (void**)vptr,*iarray_start = 0;
569   DM_DA     *dd = (DM_DA*)da->data;
570 
571   PetscFunctionBegin;
572   PetscValidHeaderSpecific(da,DM_CLASSID,1);
573   if (ghosted) {
574     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
575       if (dd->arrayghostedout[i] == *iptr) {
576         iarray_start           = dd->startghostedout[i];
577         dd->arrayghostedout[i] = PETSC_NULL;
578         dd->startghostedout[i] = PETSC_NULL;
579         break;
580       }
581     }
582     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
583       if (!dd->arrayghostedin[i]){
584         dd->arrayghostedin[i] = *iptr;
585         dd->startghostedin[i] = iarray_start;
586         break;
587       }
588     }
589   } else {
590     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
591       if (dd->arrayout[i] == *iptr) {
592         iarray_start    = dd->startout[i];
593         dd->arrayout[i] = PETSC_NULL;
594         dd->startout[i] = PETSC_NULL;
595         break;
596       }
597     }
598     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
599       if (!dd->arrayin[i]){
600         dd->arrayin[i]  = *iptr;
601         dd->startin[i]  = iarray_start;
602         break;
603       }
604     }
605   }
606   PetscFunctionReturn(0);
607 }
608 
609 #undef __FUNCT__
610 #define __FUNCT__ "DMDAGetAdicMFArray"
611 /*@C
612      DMDAGetAdicMFArray - Gets an array of derivative types for a DMDA for matrix-free ADIC.
613 
614      Input Parameter:
615 +    da - information about my local patch
616 -    ghosted - do you want arrays for the ghosted or nonghosted patch?
617 
618      Output Parameters:
619 +    vptr - array data structured to be passed to ad_FormFunctionLocal()
620 .    array_start - actual start of 1d array of all values that adiC can access directly (may be null)
621 -    tdof - total number of degrees of freedom represented in array_start (may be null)
622 
623      Notes:
624      The vector values are NOT initialized and may have garbage in them, so you may need
625      to zero them.
626 
627      This routine returns the same type of object as the DMDAVecGetArray(), except its
628      elements are derivative types instead of PetscScalars.
629 
630      Level: advanced
631 
632 .seealso: DMDARestoreAdicMFArray(), DMDAGetArray(), DMDAGetAdicArray()
633 
634 @*/
635 PetscErrorCode  DMDAGetAdicMFArray(DM da,PetscBool  ghosted,void *vptr,void *array_start,PetscInt *tdof)
636 {
637   PetscErrorCode ierr;
638   PetscInt       j,i,xs,ys,xm,ym,zs,zm,itdof = 0;
639   char           *iarray_start;
640   void           **iptr = (void**)vptr;
641   DM_DA          *dd = (DM_DA*)da->data;
642 
643   PetscFunctionBegin;
644   PetscValidHeaderSpecific(da,DM_CLASSID,1);
645   if (ghosted) {
646     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
647       if (dd->admfarrayghostedin[i]) {
648         *iptr                     = dd->admfarrayghostedin[i];
649         iarray_start              = (char*)dd->admfstartghostedin[i];
650         itdof                     = dd->ghostedtdof;
651         dd->admfarrayghostedin[i] = PETSC_NULL;
652         dd->admfstartghostedin[i] = PETSC_NULL;
653 
654         goto done;
655       }
656     }
657     xs = dd->Xs;
658     ys = dd->Ys;
659     zs = dd->Zs;
660     xm = dd->Xe-dd->Xs;
661     ym = dd->Ye-dd->Ys;
662     zm = dd->Ze-dd->Zs;
663   } else {
664     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
665       if (dd->admfarrayin[i]) {
666         *iptr              = dd->admfarrayin[i];
667         iarray_start       = (char*)dd->admfstartin[i];
668         itdof              = dd->tdof;
669         dd->admfarrayin[i] = PETSC_NULL;
670         dd->admfstartin[i] = PETSC_NULL;
671 
672         goto done;
673       }
674     }
675     xs = dd->xs;
676     ys = dd->ys;
677     zs = dd->zs;
678     xm = dd->xe-dd->xs;
679     ym = dd->ye-dd->ys;
680     zm = dd->ze-dd->zs;
681   }
682 
683   switch (dd->dim) {
684     case 1: {
685       void *ptr;
686       itdof = xm;
687 
688       ierr  = PetscMalloc(xm*2*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
689 
690       ptr   = (void*)(iarray_start - xs*2*sizeof(PetscScalar));
691       *iptr = (void*)ptr;
692       break;}
693     case 2: {
694       void **ptr;
695       itdof = xm*ym;
696 
697       ierr  = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*2*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
698 
699       ptr  = (void**)(iarray_start + xm*ym*2*sizeof(PetscScalar) - ys*sizeof(void*));
700       for(j=ys;j<ys+ym;j++) {
701         ptr[j] = iarray_start + 2*sizeof(PetscScalar)*(xm*(j-ys) - xs);
702       }
703       *iptr = (void*)ptr;
704       break;}
705     case 3: {
706       void ***ptr,**bptr;
707       itdof = xm*ym*zm;
708 
709       ierr  = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*2*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
710 
711       ptr  = (void***)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) - zs*sizeof(void*));
712       bptr = (void**)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) + zm*sizeof(void**));
713       for(i=zs;i<zs+zm;i++) {
714         ptr[i] = bptr + ((i-zs)*ym* - ys)*sizeof(void*);
715       }
716       for (i=zs; i<zs+zm; i++) {
717         for (j=ys; j<ys+ym; j++) {
718           ptr[i][j] = iarray_start + 2*sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
719         }
720       }
721 
722       *iptr = (void*)ptr;
723       break;}
724     default:
725       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
726   }
727 
728   done:
729   /* add arrays to the checked out list */
730   if (ghosted) {
731     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
732       if (!dd->admfarrayghostedout[i]) {
733         dd->admfarrayghostedout[i] = *iptr ;
734         dd->admfstartghostedout[i] = iarray_start;
735         dd->ghostedtdof            = itdof;
736         break;
737       }
738     }
739   } else {
740     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
741       if (!dd->admfarrayout[i]) {
742         dd->admfarrayout[i] = *iptr ;
743         dd->admfstartout[i] = iarray_start;
744         dd->tdof            = itdof;
745         break;
746       }
747     }
748   }
749   if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Too many DMDA ADIC arrays obtained");
750   if (tdof)        *tdof = itdof;
751   if (array_start) *(void**)array_start = iarray_start;
752   PetscFunctionReturn(0);
753 }
754 
755 #undef __FUNCT__
756 #define __FUNCT__ "DMDAGetAdicMFArray4"
757 PetscErrorCode  DMDAGetAdicMFArray4(DM da,PetscBool  ghosted,void *vptr,void *array_start,PetscInt *tdof)
758 {
759   PetscErrorCode ierr;
760   PetscInt       j,i,xs,ys,xm,ym,itdof = 0;
761   char           *iarray_start;
762   void           **iptr = (void**)vptr;
763   DM_DA          *dd = (DM_DA*)da->data;
764 
765   PetscFunctionBegin;
766   PetscValidHeaderSpecific(da,DM_CLASSID,1);
767   if (ghosted) {
768     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
769       if (dd->admfarrayghostedin[i]) {
770         *iptr                     = dd->admfarrayghostedin[i];
771         iarray_start              = (char*)dd->admfstartghostedin[i];
772         itdof                     = dd->ghostedtdof;
773         dd->admfarrayghostedin[i] = PETSC_NULL;
774         dd->admfstartghostedin[i] = PETSC_NULL;
775 
776         goto done;
777       }
778     }
779     xs = dd->Xs;
780     ys = dd->Ys;
781     xm = dd->Xe-dd->Xs;
782     ym = dd->Ye-dd->Ys;
783   } else {
784     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
785       if (dd->admfarrayin[i]) {
786         *iptr              = dd->admfarrayin[i];
787         iarray_start       = (char*)dd->admfstartin[i];
788         itdof              = dd->tdof;
789         dd->admfarrayin[i] = PETSC_NULL;
790         dd->admfstartin[i] = PETSC_NULL;
791 
792         goto done;
793       }
794     }
795     xs = dd->xs;
796     ys = dd->ys;
797     xm = dd->xe-dd->xs;
798     ym = dd->ye-dd->ys;
799   }
800 
801   switch (dd->dim) {
802     case 2: {
803       void **ptr;
804       itdof = xm*ym;
805 
806       ierr  = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*5*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
807 
808       ptr  = (void**)(iarray_start + xm*ym*5*sizeof(PetscScalar) - ys*sizeof(void*));
809       for(j=ys;j<ys+ym;j++) {
810         ptr[j] = iarray_start + 5*sizeof(PetscScalar)*(xm*(j-ys) - xs);
811       }
812       *iptr = (void*)ptr;
813       break;}
814     default:
815       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
816   }
817 
818   done:
819   /* add arrays to the checked out list */
820   if (ghosted) {
821     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
822       if (!dd->admfarrayghostedout[i]) {
823         dd->admfarrayghostedout[i] = *iptr ;
824         dd->admfstartghostedout[i] = iarray_start;
825         dd->ghostedtdof            = itdof;
826         break;
827       }
828     }
829   } else {
830     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
831       if (!dd->admfarrayout[i]) {
832         dd->admfarrayout[i] = *iptr ;
833         dd->admfstartout[i] = iarray_start;
834         dd->tdof            = itdof;
835         break;
836       }
837     }
838   }
839   if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Too many DMDA ADIC arrays obtained");
840   if (tdof)        *tdof = itdof;
841   if (array_start) *(void**)array_start = iarray_start;
842   PetscFunctionReturn(0);
843 }
844 
845 #undef __FUNCT__
846 #define __FUNCT__ "DMDAGetAdicMFArray9"
847 PetscErrorCode  DMDAGetAdicMFArray9(DM da,PetscBool  ghosted,void *vptr,void *array_start,PetscInt *tdof)
848 {
849   PetscErrorCode ierr;
850   PetscInt       j,i,xs,ys,xm,ym,itdof = 0;
851   char           *iarray_start;
852   void           **iptr = (void**)vptr;
853   DM_DA          *dd = (DM_DA*)da->data;
854 
855   PetscFunctionBegin;
856   PetscValidHeaderSpecific(da,DM_CLASSID,1);
857   if (ghosted) {
858     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
859       if (dd->admfarrayghostedin[i]) {
860         *iptr                     = dd->admfarrayghostedin[i];
861         iarray_start              = (char*)dd->admfstartghostedin[i];
862         itdof                     = dd->ghostedtdof;
863         dd->admfarrayghostedin[i] = PETSC_NULL;
864         dd->admfstartghostedin[i] = PETSC_NULL;
865 
866         goto done;
867       }
868     }
869     xs = dd->Xs;
870     ys = dd->Ys;
871     xm = dd->Xe-dd->Xs;
872     ym = dd->Ye-dd->Ys;
873   } else {
874     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
875       if (dd->admfarrayin[i]) {
876         *iptr              = dd->admfarrayin[i];
877         iarray_start       = (char*)dd->admfstartin[i];
878         itdof              = dd->tdof;
879         dd->admfarrayin[i] = PETSC_NULL;
880         dd->admfstartin[i] = PETSC_NULL;
881 
882         goto done;
883       }
884     }
885     xs = dd->xs;
886     ys = dd->ys;
887     xm = dd->xe-dd->xs;
888     ym = dd->ye-dd->ys;
889   }
890 
891   switch (dd->dim) {
892     case 2: {
893       void **ptr;
894       itdof = xm*ym;
895 
896       ierr  = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*10*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
897 
898       ptr  = (void**)(iarray_start + xm*ym*10*sizeof(PetscScalar) - ys*sizeof(void*));
899       for(j=ys;j<ys+ym;j++) {
900         ptr[j] = iarray_start + 10*sizeof(PetscScalar)*(xm*(j-ys) - xs);
901       }
902       *iptr = (void*)ptr;
903       break;}
904     default:
905       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
906   }
907 
908   done:
909   /* add arrays to the checked out list */
910   if (ghosted) {
911     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
912       if (!dd->admfarrayghostedout[i]) {
913         dd->admfarrayghostedout[i] = *iptr ;
914         dd->admfstartghostedout[i] = iarray_start;
915         dd->ghostedtdof            = itdof;
916         break;
917       }
918     }
919   } else {
920     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
921       if (!dd->admfarrayout[i]) {
922         dd->admfarrayout[i] = *iptr ;
923         dd->admfstartout[i] = iarray_start;
924         dd->tdof            = itdof;
925         break;
926       }
927     }
928   }
929   if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Too many DMDA ADIC arrays obtained");
930   if (tdof)        *tdof = itdof;
931   if (array_start) *(void**)array_start = iarray_start;
932   PetscFunctionReturn(0);
933 }
934 
935 #undef __FUNCT__
936 #define __FUNCT__ "DMDAGetAdicMFArrayb"
937 /*@C
938      DMDAGetAdicMFArrayb - Gets an array of derivative types for a DMDA for matrix-free ADIC.
939 
940      Input Parameter:
941 +    da - information about my local patch
942 -    ghosted - do you want arrays for the ghosted or nonghosted patch?
943 
944      Output Parameters:
945 +    vptr - array data structured to be passed to ad_FormFunctionLocal()
946 .    array_start - actual start of 1d array of all values that adiC can access directly (may be null)
947 -    tdof - total number of degrees of freedom represented in array_start (may be null)
948 
949      Notes:
950      The vector values are NOT initialized and may have garbage in them, so you may need
951      to zero them.
952 
953      This routine returns the same type of object as the DMDAVecGetArray(), except its
954      elements are derivative types instead of PetscScalars.
955 
956      Level: advanced
957 
958 .seealso: DMDARestoreAdicMFArray(), DMDAGetArray(), DMDAGetAdicArray()
959 
960 @*/
961 PetscErrorCode  DMDAGetAdicMFArrayb(DM da,PetscBool  ghosted,void *vptr,void *array_start,PetscInt *tdof)
962 {
963   PetscErrorCode ierr;
964   PetscInt       j,i,xs,ys,xm,ym,zs,zm,itdof = 0;
965   char           *iarray_start;
966   void           **iptr = (void**)vptr;
967   DM_DA          *dd = (DM_DA*)da->data;
968   PetscInt       bs = dd->w,bs1 = bs+1;
969 
970   PetscFunctionBegin;
971   PetscValidHeaderSpecific(da,DM_CLASSID,1);
972   if (ghosted) {
973     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
974       if (dd->admfarrayghostedin[i]) {
975         *iptr                     = dd->admfarrayghostedin[i];
976         iarray_start              = (char*)dd->admfstartghostedin[i];
977         itdof                     = dd->ghostedtdof;
978         dd->admfarrayghostedin[i] = PETSC_NULL;
979         dd->admfstartghostedin[i] = PETSC_NULL;
980 
981         goto done;
982       }
983     }
984     xs = dd->Xs;
985     ys = dd->Ys;
986     zs = dd->Zs;
987     xm = dd->Xe-dd->Xs;
988     ym = dd->Ye-dd->Ys;
989     zm = dd->Ze-dd->Zs;
990   } else {
991     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
992       if (dd->admfarrayin[i]) {
993         *iptr              = dd->admfarrayin[i];
994         iarray_start       = (char*)dd->admfstartin[i];
995         itdof              = dd->tdof;
996         dd->admfarrayin[i] = PETSC_NULL;
997         dd->admfstartin[i] = PETSC_NULL;
998 
999         goto done;
1000       }
1001     }
1002     xs = dd->xs;
1003     ys = dd->ys;
1004     zs = dd->zs;
1005     xm = dd->xe-dd->xs;
1006     ym = dd->ye-dd->ys;
1007     zm = dd->ze-dd->zs;
1008   }
1009 
1010   switch (dd->dim) {
1011     case 1: {
1012       void *ptr;
1013       itdof = xm;
1014 
1015       ierr  = PetscMalloc(xm*bs1*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
1016 
1017       ptr   = (void*)(iarray_start - xs*bs1*sizeof(PetscScalar));
1018       *iptr = (void*)ptr;
1019       break;}
1020     case 2: {
1021       void **ptr;
1022       itdof = xm*ym;
1023 
1024       ierr  = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*bs1*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
1025 
1026       ptr  = (void**)(iarray_start + xm*ym*bs1*sizeof(PetscScalar) - ys*sizeof(void*));
1027       for(j=ys;j<ys+ym;j++) {
1028         ptr[j] = iarray_start + bs1*sizeof(PetscScalar)*(xm*(j-ys) - xs);
1029       }
1030       *iptr = (void*)ptr;
1031       break;}
1032     case 3: {
1033       void ***ptr,**bptr;
1034       itdof = xm*ym*zm;
1035 
1036       ierr  = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*bs1*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
1037 
1038       ptr  = (void***)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) - zs*sizeof(void*));
1039       bptr = (void**)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) + zm*sizeof(void**));
1040       for(i=zs;i<zs+zm;i++) {
1041         ptr[i] = bptr + ((i-zs)*ym* - ys)*sizeof(void*);
1042       }
1043       for (i=zs; i<zs+zm; i++) {
1044         for (j=ys; j<ys+ym; j++) {
1045           ptr[i][j] = iarray_start + bs1*sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
1046         }
1047       }
1048 
1049       *iptr = (void*)ptr;
1050       break;}
1051     default:
1052       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
1053   }
1054 
1055   done:
1056   /* add arrays to the checked out list */
1057   if (ghosted) {
1058     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1059       if (!dd->admfarrayghostedout[i]) {
1060         dd->admfarrayghostedout[i] = *iptr ;
1061         dd->admfstartghostedout[i] = iarray_start;
1062         dd->ghostedtdof            = itdof;
1063         break;
1064       }
1065     }
1066   } else {
1067     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1068       if (!dd->admfarrayout[i]) {
1069         dd->admfarrayout[i] = *iptr ;
1070         dd->admfstartout[i] = iarray_start;
1071         dd->tdof            = itdof;
1072         break;
1073       }
1074     }
1075   }
1076   if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Too many DMDA ADIC arrays obtained");
1077   if (tdof)        *tdof = itdof;
1078   if (array_start) *(void**)array_start = iarray_start;
1079   PetscFunctionReturn(0);
1080 }
1081 
1082 #undef __FUNCT__
1083 #define __FUNCT__ "DMDARestoreAdicMFArray"
1084 /*@C
1085      DMDARestoreAdicMFArray - Restores an array of derivative types for a DMDA.
1086 
1087      Input Parameter:
1088 +    da - information about my local patch
1089 -    ghosted - do you want arrays for the ghosted or nonghosted patch?
1090 
1091      Output Parameters:
1092 +    ptr - array data structure to be passed to ad_FormFunctionLocal()
1093 .    array_start - actual start of 1d array of all values that adiC can access directly
1094 -    tdof - total number of degrees of freedom represented in array_start
1095 
1096      Level: advanced
1097 
1098 .seealso: DMDAGetAdicArray()
1099 
1100 @*/
1101 PetscErrorCode  DMDARestoreAdicMFArray(DM da,PetscBool  ghosted,void *vptr,void *array_start,PetscInt *tdof)
1102 {
1103   PetscInt  i;
1104   void      **iptr = (void**)vptr,*iarray_start = 0;
1105   DM_DA     *dd = (DM_DA*)da->data;
1106 
1107   PetscFunctionBegin;
1108   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1109   if (ghosted) {
1110     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1111       if (dd->admfarrayghostedout[i] == *iptr) {
1112         iarray_start               = dd->admfstartghostedout[i];
1113         dd->admfarrayghostedout[i] = PETSC_NULL;
1114         dd->admfstartghostedout[i] = PETSC_NULL;
1115         break;
1116       }
1117     }
1118     if (!iarray_start) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Could not find array in checkout list");
1119     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1120       if (!dd->admfarrayghostedin[i]){
1121         dd->admfarrayghostedin[i] = *iptr;
1122         dd->admfstartghostedin[i] = iarray_start;
1123         break;
1124       }
1125     }
1126   } else {
1127     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1128       if (dd->admfarrayout[i] == *iptr) {
1129         iarray_start        = dd->admfstartout[i];
1130         dd->admfarrayout[i] = PETSC_NULL;
1131         dd->admfstartout[i] = PETSC_NULL;
1132         break;
1133       }
1134     }
1135     if (!iarray_start) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Could not find array in checkout list");
1136     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1137       if (!dd->admfarrayin[i]){
1138         dd->admfarrayin[i] = *iptr;
1139         dd->admfstartin[i] = iarray_start;
1140         break;
1141       }
1142     }
1143   }
1144   PetscFunctionReturn(0);
1145 }
1146 
1147 #undef __FUNCT__
1148 #define __FUNCT__ "admf_DAGetArray"
1149 PetscErrorCode  admf_DAGetArray(DM da,PetscBool  ghosted,void *iptr)
1150 {
1151   PetscErrorCode ierr;
1152   PetscFunctionBegin;
1153   ierr = DMDAGetAdicMFArray(da,ghosted,iptr,0,0);CHKERRQ(ierr);
1154   PetscFunctionReturn(0);
1155 }
1156 
1157 #undef __FUNCT__
1158 #define __FUNCT__ "admf_DARestoreArray"
1159 PetscErrorCode  admf_DARestoreArray(DM da,PetscBool  ghosted,void *iptr)
1160 {
1161   PetscErrorCode ierr;
1162   PetscFunctionBegin;
1163   ierr = DMDARestoreAdicMFArray(da,ghosted,iptr,0,0);CHKERRQ(ierr);
1164   PetscFunctionReturn(0);
1165 }
1166 
1167