xref: /petsc/src/dm/impls/da/dalocal.c (revision dfe153156869f5931ee6b299eecdcfd8470c0dd0)
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 = VecSetType(*g,da->vectype);CHKERRQ(ierr);
61   ierr = VecSetBlockSize(*g,dd->w);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__ "DMGetLocalVector"
73 /*@
74    DMGetLocalVector - Gets a Seq PETSc vector that
75    may be used with the DMXXX routines. This vector has spaces for the ghost values.
76 
77    Not Collective
78 
79    Input Parameter:
80 .  dm - the distributed array
81 
82    Output Parameter:
83 .  g - the local vector
84 
85    Level: beginner
86 
87    Note:
88    The vector values are NOT initialized and may have garbage in them, so you may need
89    to zero them.
90 
91    The output parameter, g, is a regular PETSc vector that should be returned with
92    DMRestoreLocalVector() DO NOT call VecDestroy() on it.
93 
94    VecStride*() operations can be useful when using DM with dof > 1
95 
96 .keywords: distributed array, create, local, vector
97 
98 .seealso: DMCreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
99           DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGlobalToLocalBegin(),
100           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMCreateLocalVector(), DMRestoreLocalVector(),
101           VecStrideMax(), VecStrideMin(), VecStrideNorm()
102 @*/
103 PetscErrorCode  DMGetLocalVector(DM dm,Vec* g)
104 {
105   PetscErrorCode ierr,i;
106 
107   PetscFunctionBegin;
108   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
109   PetscValidPointer(g,2);
110   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
111     if (dm->localin[i]) {
112       *g             = dm->localin[i];
113       dm->localin[i] = PETSC_NULL;
114       goto alldone;
115     }
116   }
117   ierr = DMCreateLocalVector(dm,g);CHKERRQ(ierr);
118 
119   alldone:
120   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
121     if (!dm->localout[i]) {
122       dm->localout[i] = *g;
123       break;
124     }
125   }
126   PetscFunctionReturn(0);
127 }
128 
129 #undef __FUNCT__
130 #define __FUNCT__ "DMRestoreLocalVector"
131 /*@
132    DMRestoreLocalVector - Returns a Seq PETSc vector that
133      obtained from DMGetLocalVector(). Do not use with vector obtained via
134      DMCreateLocalVector().
135 
136    Not Collective
137 
138    Input Parameter:
139 +  dm - the distributed array
140 -  g - the local vector
141 
142    Level: beginner
143 
144 .keywords: distributed array, create, local, vector
145 
146 .seealso: DMCreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
147           DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGlobalToLocalBegin(),
148           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMCreateLocalVector(), DMGetLocalVector()
149 @*/
150 PetscErrorCode  DMRestoreLocalVector(DM dm,Vec* g)
151 {
152   PetscErrorCode ierr;
153   PetscInt       i,j;
154 
155   PetscFunctionBegin;
156   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
157   PetscValidPointer(g,2);
158   for (j=0; j<DM_MAX_WORK_VECTORS; j++) {
159     if (*g == dm->localout[j]) {
160       dm->localout[j] = PETSC_NULL;
161       for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
162         if (!dm->localin[i]) {
163           dm->localin[i] = *g;
164           goto alldone;
165         }
166       }
167     }
168   }
169   ierr = VecDestroy(g);CHKERRQ(ierr);
170   alldone:
171   PetscFunctionReturn(0);
172 }
173 
174 #undef __FUNCT__
175 #define __FUNCT__ "DMGetGlobalVector"
176 /*@
177    DMGetGlobalVector - Gets a MPI PETSc vector that
178    may be used with the DMXXX routines.
179 
180    Collective on DM
181 
182    Input Parameter:
183 .  dm - the distributed array
184 
185    Output Parameter:
186 .  g - the global vector
187 
188    Level: beginner
189 
190    Note:
191    The vector values are NOT initialized and may have garbage in them, so you may need
192    to zero them.
193 
194    The output parameter, g, is a regular PETSc vector that should be returned with
195    DMRestoreGlobalVector() DO NOT call VecDestroy() on it.
196 
197    VecStride*() operations can be useful when using DM with dof > 1
198 
199 .keywords: distributed array, create, Global, vector
200 
201 .seealso: DMCreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
202           DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGlobalToLocalBegin(),
203           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMCreateLocalVector(), DMRestoreLocalVector()
204           VecStrideMax(), VecStrideMin(), VecStrideNorm()
205 
206 @*/
207 PetscErrorCode  DMGetGlobalVector(DM dm,Vec* g)
208 {
209   PetscErrorCode ierr;
210   PetscInt       i;
211 
212   PetscFunctionBegin;
213   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
214   PetscValidPointer(g,2);
215   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
216     if (dm->globalin[i]) {
217       *g             = dm->globalin[i];
218       dm->globalin[i] = PETSC_NULL;
219       goto alldone;
220     }
221   }
222   ierr = DMCreateGlobalVector(dm,g);CHKERRQ(ierr);
223 
224   alldone:
225   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
226     if (!dm->globalout[i]) {
227       dm->globalout[i] = *g;
228       break;
229     }
230   }
231   PetscFunctionReturn(0);
232 }
233 
234 #undef __FUNCT__
235 #define __FUNCT__ "DMClearGlobalVectors"
236 /*@
237    DMClearGlobalVectors - Destroys all the global vectors that have been stashed in this DM
238 
239    Collective on DM
240 
241    Input Parameter:
242 .  dm - the distributed array
243 
244    Level: developer
245 
246 .keywords: distributed array, create, Global, vector
247 
248 .seealso: DMCreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
249           DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGlobalToLocalBegin(),
250           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMCreateLocalVector(), DMRestoreLocalVector()
251           VecStrideMax(), VecStrideMin(), VecStrideNorm()
252 
253 @*/
254 PetscErrorCode  DMClearGlobalVectors(DM dm)
255 {
256   PetscErrorCode ierr;
257   PetscInt       i;
258 
259   PetscFunctionBegin;
260   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
261   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
262     if (dm->globalout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Clearing DM of global vectors that has a global vector obtained with DMGetGlobalVector()");
263     ierr = VecDestroy(&dm->globalin[i]);CHKERRQ(ierr);
264   }
265   PetscFunctionReturn(0);
266 }
267 
268 #undef __FUNCT__
269 #define __FUNCT__ "DMRestoreGlobalVector"
270 /*@
271    DMRestoreGlobalVector - Returns a Seq PETSc vector that
272      obtained from DMGetGlobalVector(). Do not use with vector obtained via
273      DMCreateGlobalVector().
274 
275    Not Collective
276 
277    Input Parameter:
278 +  dm - the distributed array
279 -  g - the global vector
280 
281    Level: beginner
282 
283 .keywords: distributed array, create, global, vector
284 
285 .seealso: DMCreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
286           DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGlobalToGlobalBegin(),
287           DMGlobalToGlobalEnd(), DMGlobalToGlobal(), DMCreateLocalVector(), DMGetGlobalVector()
288 @*/
289 PetscErrorCode  DMRestoreGlobalVector(DM dm,Vec* g)
290 {
291   PetscErrorCode ierr;
292   PetscInt       i,j;
293 
294   PetscFunctionBegin;
295   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
296   PetscValidPointer(g,2);
297   for (j=0; j<DM_MAX_WORK_VECTORS; j++) {
298     if (*g == dm->globalout[j]) {
299       dm->globalout[j] = PETSC_NULL;
300       for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
301         if (!dm->globalin[i]) {
302           dm->globalin[i] = *g;
303           goto alldone;
304         }
305       }
306     }
307   }
308   ierr = VecDestroy(g);CHKERRQ(ierr);
309   alldone:
310   PetscFunctionReturn(0);
311 }
312 
313 #undef __FUNCT__
314 #define __FUNCT__ "DMGetNamedGlobalVector"
315 /*@C
316    DMGetNamedGlobalVector - get access to a named, persistent global vector
317 
318    Collective on DM
319 
320    Input Arguments:
321 +  dm - DM to hold named vectors
322 -  name - unique name for Vec
323 
324    Output Arguments:
325 .  X - named Vec
326 
327    Level: developer
328 
329    Note: If a Vec with the given name does not exist, it is created.
330 
331 .seealso: DMRestoreNamedGlobalVector()
332 @*/
333 PetscErrorCode DMGetNamedGlobalVector(DM dm,const char *name,Vec *X)
334 {
335   PetscErrorCode ierr;
336   DMNamedVecLink link;
337 
338   PetscFunctionBegin;
339   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
340   PetscValidCharPointer(name,2);
341   PetscValidPointer(X,3);
342   for (link=dm->namedglobal; link; link=link->next) {
343     PetscBool match;
344     ierr = PetscStrcmp(name,link->name,&match);CHKERRQ(ierr);
345     if (match) {
346       if (link->status != DMVEC_STATUS_IN) SETERRQ1(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Vec name '%s' already checked out",name);
347       goto found;
348     }
349   }
350 
351   /* Create the Vec */
352   ierr = PetscMalloc(sizeof *link,&link);CHKERRQ(ierr);
353   ierr = PetscStrallocpy(name,&link->name);CHKERRQ(ierr);
354   ierr = DMCreateGlobalVector(dm,&link->X);CHKERRQ(ierr);
355   link->next = dm->namedglobal;
356   dm->namedglobal = link;
357 
358   found:
359   *X = link->X;
360   link->status = DMVEC_STATUS_OUT;
361   PetscFunctionReturn(0);
362 }
363 
364 #undef __FUNCT__
365 #define __FUNCT__ "DMRestoreNamedGlobalVector"
366 /*@C
367    DMRestoreNamedGlobalVector - restore access to a named, persistent global vector
368 
369    Collective on DM
370 
371    Input Arguments:
372 +  dm - DM on which the vector was gotten
373 .  name - name under which the vector was gotten
374 -  X - Vec to restore
375 
376    Output Arguments:
377 
378    Level: developer
379 
380 .seealso: DMGetNamedGlobalVector()
381 @*/
382 PetscErrorCode DMRestoreNamedGlobalVector(DM dm,const char *name,Vec *X)
383 {
384   PetscErrorCode ierr;
385   DMNamedVecLink link;
386 
387   PetscFunctionBegin;
388   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
389   PetscValidCharPointer(name,2);
390   PetscValidPointer(X,3);
391   PetscValidHeaderSpecific(*X,VEC_CLASSID,3);
392   for (link=dm->namedglobal; link; link=link->next) {
393     PetscBool match;
394     ierr = PetscStrcmp(name,link->name,&match);CHKERRQ(ierr);
395     if (match) {
396       if (link->status != DMVEC_STATUS_OUT) SETERRQ1(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Vec name '%s' was not checked out",name);
397       if (link->X != *X) SETERRQ1(((PetscObject)dm)->comm,PETSC_ERR_ARG_INCOMP,"Attempt to restore Vec name '%s', but Vec does not match the cache",name);
398       link->status = DMVEC_STATUS_IN;
399       *X = PETSC_NULL;
400       PetscFunctionReturn(0);
401     }
402   }
403   SETERRQ1(((PetscObject)dm)->comm,PETSC_ERR_ARG_INCOMP,"Could not find Vec name '%s' to restore",name);
404   PetscFunctionReturn(0);
405 }
406 
407 /* ------------------------------------------------------------------- */
408 #if defined(PETSC_HAVE_ADIC)
409 
410 EXTERN_C_BEGIN
411 #include <adic/ad_utils.h>
412 EXTERN_C_END
413 
414 #undef __FUNCT__
415 #define __FUNCT__ "DMDAGetAdicArray"
416 /*@C
417      DMDAGetAdicArray - Gets an array of derivative types for a DMDA
418 
419     Input Parameter:
420 +    da - information about my local patch
421 -    ghosted - do you want arrays for the ghosted or nonghosted patch
422 
423     Output Parameters:
424 +    vptr - array data structured to be passed to ad_FormFunctionLocal()
425 .    array_start - actual start of 1d array of all values that adiC can access directly (may be null)
426 -    tdof - total number of degrees of freedom represented in array_start (may be null)
427 
428      Notes:
429        The vector values are NOT initialized and may have garbage in them, so you may need
430        to zero them.
431 
432        Returns the same type of object as the DMDAVecGetArray() except its elements are
433            derivative types instead of PetscScalars
434 
435      Level: advanced
436 
437 .seealso: DMDARestoreAdicArray()
438 
439 @*/
440 PetscErrorCode  DMDAGetAdicArray(DM da,PetscBool  ghosted,void *vptr,void *array_start,PetscInt *tdof)
441 {
442   PetscErrorCode ierr;
443   PetscInt       j,i,deriv_type_size,xs,ys,xm,ym,zs,zm,itdof;
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_AD_ARRAYS; i++) {
452       if (dd->adarrayghostedin[i]) {
453         *iptr                   = dd->adarrayghostedin[i];
454         iarray_start            = (char*)dd->adstartghostedin[i];
455         itdof                   = dd->ghostedtdof;
456         dd->adarrayghostedin[i] = PETSC_NULL;
457         dd->adstartghostedin[i] = PETSC_NULL;
458 
459         goto done;
460       }
461     }
462     xs = dd->Xs;
463     ys = dd->Ys;
464     zs = dd->Zs;
465     xm = dd->Xe-dd->Xs;
466     ym = dd->Ye-dd->Ys;
467     zm = dd->Ze-dd->Zs;
468   } else {
469     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
470       if (dd->adarrayin[i]) {
471         *iptr            = dd->adarrayin[i];
472         iarray_start     = (char*)dd->adstartin[i];
473         itdof            = dd->tdof;
474         dd->adarrayin[i] = PETSC_NULL;
475         dd->adstartin[i] = PETSC_NULL;
476 
477         goto done;
478       }
479     }
480     xs = dd->xs;
481     ys = dd->ys;
482     zs = dd->zs;
483     xm = dd->xe-dd->xs;
484     ym = dd->ye-dd->ys;
485     zm = dd->ze-dd->zs;
486   }
487   deriv_type_size = PetscADGetDerivTypeSize();
488 
489   switch (dd->dim) {
490     case 1: {
491       void *ptr;
492       itdof = xm;
493 
494       ierr  = PetscMalloc(xm*deriv_type_size,&iarray_start);CHKERRQ(ierr);
495 
496       ptr   = (void*)(iarray_start - xs*deriv_type_size);
497       *iptr = (void*)ptr;
498       break;}
499     case 2: {
500       void **ptr;
501       itdof = xm*ym;
502 
503       ierr  = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*deriv_type_size,&iarray_start);CHKERRQ(ierr);
504 
505       ptr  = (void**)(iarray_start + xm*ym*deriv_type_size - ys*sizeof(void*));
506       for(j=ys;j<ys+ym;j++) {
507         ptr[j] = iarray_start + deriv_type_size*(xm*(j-ys) - xs);
508       }
509       *iptr = (void*)ptr;
510       break;}
511     case 3: {
512       void ***ptr,**bptr;
513       itdof = xm*ym*zm;
514 
515       ierr  = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*deriv_type_size,&iarray_start);CHKERRQ(ierr);
516 
517       ptr  = (void***)(iarray_start + xm*ym*zm*deriv_type_size - zs*sizeof(void*));
518       bptr = (void**)(iarray_start + xm*ym*zm*deriv_type_size + zm*sizeof(void**));
519 
520       for(i=zs;i<zs+zm;i++) {
521         ptr[i] = bptr + ((i-zs)*ym - ys);
522       }
523       for (i=zs; i<zs+zm; i++) {
524         for (j=ys; j<ys+ym; j++) {
525           ptr[i][j] = iarray_start + deriv_type_size*(xm*ym*(i-zs) + xm*(j-ys) - xs);
526         }
527       }
528 
529       *iptr = (void*)ptr;
530       break;}
531     default:
532       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
533   }
534 
535   done:
536   /* add arrays to the checked out list */
537   if (ghosted) {
538     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
539       if (!dd->adarrayghostedout[i]) {
540         dd->adarrayghostedout[i] = *iptr ;
541         dd->adstartghostedout[i] = iarray_start;
542         dd->ghostedtdof          = itdof;
543         break;
544       }
545     }
546   } else {
547     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
548       if (!dd->adarrayout[i]) {
549         dd->adarrayout[i] = *iptr ;
550         dd->adstartout[i] = iarray_start;
551         dd->tdof          = itdof;
552         break;
553       }
554     }
555   }
556   if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Too many DMDA ADIC arrays obtained");
557   if (tdof)        *tdof = itdof;
558   if (array_start) *(void**)array_start = iarray_start;
559   PetscFunctionReturn(0);
560 }
561 
562 #undef __FUNCT__
563 #define __FUNCT__ "DMDARestoreAdicArray"
564 /*@C
565      DMDARestoreAdicArray - Restores an array of derivative types for a DMDA
566 
567     Input Parameter:
568 +    da - information about my local patch
569 -    ghosted - do you want arrays for the ghosted or nonghosted patch
570 
571     Output Parameters:
572 +    ptr - array data structured to be passed to ad_FormFunctionLocal()
573 .    array_start - actual start of 1d array of all values that adiC can access directly
574 -    tdof - total number of degrees of freedom represented in array_start
575 
576      Level: advanced
577 
578 .seealso: DMDAGetAdicArray()
579 
580 @*/
581 PetscErrorCode  DMDARestoreAdicArray(DM da,PetscBool  ghosted,void *ptr,void *array_start,PetscInt *tdof)
582 {
583   PetscInt  i;
584   void      **iptr = (void**)ptr,iarray_start = 0;
585 
586   PetscFunctionBegin;
587   PetscValidHeaderSpecific(da,DM_CLASSID,1);
588   if (ghosted) {
589     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
590       if (dd->adarrayghostedout[i] == *iptr) {
591         iarray_start             = dd->adstartghostedout[i];
592         dd->adarrayghostedout[i] = PETSC_NULL;
593         dd->adstartghostedout[i] = PETSC_NULL;
594         break;
595       }
596     }
597     if (!iarray_start) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Could not find array in checkout list");
598     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
599       if (!dd->adarrayghostedin[i]){
600         dd->adarrayghostedin[i] = *iptr;
601         dd->adstartghostedin[i] = iarray_start;
602         break;
603       }
604     }
605   } else {
606     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
607       if (dd->adarrayout[i] == *iptr) {
608         iarray_start      = dd->adstartout[i];
609         dd->adarrayout[i] = PETSC_NULL;
610         dd->adstartout[i] = PETSC_NULL;
611         break;
612       }
613     }
614     if (!iarray_start) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Could not find array in checkout list");
615     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
616       if (!dd->adarrayin[i]){
617         dd->adarrayin[i]   = *iptr;
618         dd->adstartin[i]   = iarray_start;
619         break;
620       }
621     }
622   }
623   PetscFunctionReturn(0);
624 }
625 
626 #undef __FUNCT__
627 #define __FUNCT__ "ad_DAGetArray"
628 PetscErrorCode  ad_DAGetArray(DM da,PetscBool  ghosted,void *iptr)
629 {
630   PetscErrorCode ierr;
631   PetscFunctionBegin;
632   ierr = DMDAGetAdicArray(da,ghosted,iptr,0,0);CHKERRQ(ierr);
633   PetscFunctionReturn(0);
634 }
635 
636 #undef __FUNCT__
637 #define __FUNCT__ "ad_DARestoreArray"
638 PetscErrorCode  ad_DARestoreArray(DM da,PetscBool  ghosted,void *iptr)
639 {
640   PetscErrorCode ierr;
641   PetscFunctionBegin;
642   ierr = DMDARestoreAdicArray(da,ghosted,iptr,0,0);CHKERRQ(ierr);
643   PetscFunctionReturn(0);
644 }
645 
646 #endif
647 
648 #undef __FUNCT__
649 #define __FUNCT__ "DMDAGetArray"
650 /*@C
651      DMDAGetArray - Gets a work array for a DMDA
652 
653     Input Parameter:
654 +    da - information about my local patch
655 -    ghosted - do you want arrays for the ghosted or nonghosted patch
656 
657     Output Parameters:
658 .    vptr - array data structured
659 
660     Note:  The vector values are NOT initialized and may have garbage in them, so you may need
661            to zero them.
662 
663   Level: advanced
664 
665 .seealso: DMDARestoreArray(), DMDAGetAdicArray()
666 
667 @*/
668 PetscErrorCode  DMDAGetArray(DM da,PetscBool  ghosted,void *vptr)
669 {
670   PetscErrorCode ierr;
671   PetscInt       j,i,xs,ys,xm,ym,zs,zm;
672   char           *iarray_start;
673   void           **iptr = (void**)vptr;
674   DM_DA          *dd = (DM_DA*)da->data;
675 
676   PetscFunctionBegin;
677   PetscValidHeaderSpecific(da,DM_CLASSID,1);
678   if (ghosted) {
679     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
680       if (dd->arrayghostedin[i]) {
681         *iptr                 = dd->arrayghostedin[i];
682         iarray_start          = (char*)dd->startghostedin[i];
683         dd->arrayghostedin[i] = PETSC_NULL;
684         dd->startghostedin[i] = PETSC_NULL;
685 
686         goto done;
687       }
688     }
689     xs = dd->Xs;
690     ys = dd->Ys;
691     zs = dd->Zs;
692     xm = dd->Xe-dd->Xs;
693     ym = dd->Ye-dd->Ys;
694     zm = dd->Ze-dd->Zs;
695   } else {
696     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
697       if (dd->arrayin[i]) {
698         *iptr          = dd->arrayin[i];
699         iarray_start   = (char*)dd->startin[i];
700         dd->arrayin[i] = PETSC_NULL;
701         dd->startin[i] = PETSC_NULL;
702 
703         goto done;
704       }
705     }
706     xs = dd->xs;
707     ys = dd->ys;
708     zs = dd->zs;
709     xm = dd->xe-dd->xs;
710     ym = dd->ye-dd->ys;
711     zm = dd->ze-dd->zs;
712   }
713 
714   switch (dd->dim) {
715     case 1: {
716       void *ptr;
717 
718       ierr  = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
719 
720       ptr   = (void*)(iarray_start - xs*sizeof(PetscScalar));
721       *iptr = (void*)ptr;
722       break;}
723     case 2: {
724       void **ptr;
725 
726       ierr  = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
727 
728       ptr  = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*));
729       for(j=ys;j<ys+ym;j++) {
730         ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs);
731       }
732       *iptr = (void*)ptr;
733       break;}
734     case 3: {
735       void ***ptr,**bptr;
736 
737       ierr  = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
738 
739       ptr  = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*));
740       bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**));
741       for(i=zs;i<zs+zm;i++) {
742         ptr[i] = bptr + ((i-zs)*ym - ys);
743       }
744       for (i=zs; i<zs+zm; i++) {
745         for (j=ys; j<ys+ym; j++) {
746           ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
747         }
748       }
749 
750       *iptr = (void*)ptr;
751       break;}
752     default:
753       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
754   }
755 
756   done:
757   /* add arrays to the checked out list */
758   if (ghosted) {
759     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
760       if (!dd->arrayghostedout[i]) {
761         dd->arrayghostedout[i] = *iptr ;
762         dd->startghostedout[i] = iarray_start;
763         break;
764       }
765     }
766   } else {
767     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
768       if (!dd->arrayout[i]) {
769         dd->arrayout[i] = *iptr ;
770         dd->startout[i] = iarray_start;
771         break;
772       }
773     }
774   }
775   PetscFunctionReturn(0);
776 }
777 
778 #undef __FUNCT__
779 #define __FUNCT__ "DMDARestoreArray"
780 /*@C
781      DMDARestoreArray - Restores an array of derivative types for a DMDA
782 
783     Input Parameter:
784 +    da - information about my local patch
785 .    ghosted - do you want arrays for the ghosted or nonghosted patch
786 -    vptr - array data structured to be passed to ad_FormFunctionLocal()
787 
788      Level: advanced
789 
790 .seealso: DMDAGetArray(), DMDAGetAdicArray()
791 
792 @*/
793 PetscErrorCode  DMDARestoreArray(DM da,PetscBool  ghosted,void *vptr)
794 {
795   PetscInt  i;
796   void      **iptr = (void**)vptr,*iarray_start = 0;
797   DM_DA     *dd = (DM_DA*)da->data;
798 
799   PetscFunctionBegin;
800   PetscValidHeaderSpecific(da,DM_CLASSID,1);
801   if (ghosted) {
802     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
803       if (dd->arrayghostedout[i] == *iptr) {
804         iarray_start           = dd->startghostedout[i];
805         dd->arrayghostedout[i] = PETSC_NULL;
806         dd->startghostedout[i] = PETSC_NULL;
807         break;
808       }
809     }
810     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
811       if (!dd->arrayghostedin[i]){
812         dd->arrayghostedin[i] = *iptr;
813         dd->startghostedin[i] = iarray_start;
814         break;
815       }
816     }
817   } else {
818     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
819       if (dd->arrayout[i] == *iptr) {
820         iarray_start    = dd->startout[i];
821         dd->arrayout[i] = PETSC_NULL;
822         dd->startout[i] = PETSC_NULL;
823         break;
824       }
825     }
826     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
827       if (!dd->arrayin[i]){
828         dd->arrayin[i]  = *iptr;
829         dd->startin[i]  = iarray_start;
830         break;
831       }
832     }
833   }
834   PetscFunctionReturn(0);
835 }
836 
837 #undef __FUNCT__
838 #define __FUNCT__ "DMDAGetAdicMFArray"
839 /*@C
840      DMDAGetAdicMFArray - Gets an array of derivative types for a DMDA for matrix-free ADIC.
841 
842      Input Parameter:
843 +    da - information about my local patch
844 -    ghosted - do you want arrays for the ghosted or nonghosted patch?
845 
846      Output Parameters:
847 +    vptr - array data structured to be passed to ad_FormFunctionLocal()
848 .    array_start - actual start of 1d array of all values that adiC can access directly (may be null)
849 -    tdof - total number of degrees of freedom represented in array_start (may be null)
850 
851      Notes:
852      The vector values are NOT initialized and may have garbage in them, so you may need
853      to zero them.
854 
855      This routine returns the same type of object as the DMDAVecGetArray(), except its
856      elements are derivative types instead of PetscScalars.
857 
858      Level: advanced
859 
860 .seealso: DMDARestoreAdicMFArray(), DMDAGetArray(), DMDAGetAdicArray()
861 
862 @*/
863 PetscErrorCode  DMDAGetAdicMFArray(DM da,PetscBool  ghosted,void *vptr,void *array_start,PetscInt *tdof)
864 {
865   PetscErrorCode ierr;
866   PetscInt       j,i,xs,ys,xm,ym,zs,zm,itdof = 0;
867   char           *iarray_start;
868   void           **iptr = (void**)vptr;
869   DM_DA          *dd = (DM_DA*)da->data;
870 
871   PetscFunctionBegin;
872   PetscValidHeaderSpecific(da,DM_CLASSID,1);
873   if (ghosted) {
874     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
875       if (dd->admfarrayghostedin[i]) {
876         *iptr                     = dd->admfarrayghostedin[i];
877         iarray_start              = (char*)dd->admfstartghostedin[i];
878         itdof                     = dd->ghostedtdof;
879         dd->admfarrayghostedin[i] = PETSC_NULL;
880         dd->admfstartghostedin[i] = PETSC_NULL;
881 
882         goto done;
883       }
884     }
885     xs = dd->Xs;
886     ys = dd->Ys;
887     zs = dd->Zs;
888     xm = dd->Xe-dd->Xs;
889     ym = dd->Ye-dd->Ys;
890     zm = dd->Ze-dd->Zs;
891   } else {
892     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
893       if (dd->admfarrayin[i]) {
894         *iptr              = dd->admfarrayin[i];
895         iarray_start       = (char*)dd->admfstartin[i];
896         itdof              = dd->tdof;
897         dd->admfarrayin[i] = PETSC_NULL;
898         dd->admfstartin[i] = PETSC_NULL;
899 
900         goto done;
901       }
902     }
903     xs = dd->xs;
904     ys = dd->ys;
905     zs = dd->zs;
906     xm = dd->xe-dd->xs;
907     ym = dd->ye-dd->ys;
908     zm = dd->ze-dd->zs;
909   }
910 
911   switch (dd->dim) {
912     case 1: {
913       void *ptr;
914       itdof = xm;
915 
916       ierr  = PetscMalloc(xm*2*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
917 
918       ptr   = (void*)(iarray_start - xs*2*sizeof(PetscScalar));
919       *iptr = (void*)ptr;
920       break;}
921     case 2: {
922       void **ptr;
923       itdof = xm*ym;
924 
925       ierr  = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*2*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
926 
927       ptr  = (void**)(iarray_start + xm*ym*2*sizeof(PetscScalar) - ys*sizeof(void*));
928       for(j=ys;j<ys+ym;j++) {
929         ptr[j] = iarray_start + 2*sizeof(PetscScalar)*(xm*(j-ys) - xs);
930       }
931       *iptr = (void*)ptr;
932       break;}
933     case 3: {
934       void ***ptr,**bptr;
935       itdof = xm*ym*zm;
936 
937       ierr  = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*2*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
938 
939       ptr  = (void***)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) - zs*sizeof(void*));
940       bptr = (void**)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) + zm*sizeof(void**));
941       for(i=zs;i<zs+zm;i++) {
942         ptr[i] = bptr + ((i-zs)*ym* - ys)*sizeof(void*);
943       }
944       for (i=zs; i<zs+zm; i++) {
945         for (j=ys; j<ys+ym; j++) {
946           ptr[i][j] = iarray_start + 2*sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
947         }
948       }
949 
950       *iptr = (void*)ptr;
951       break;}
952     default:
953       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
954   }
955 
956   done:
957   /* add arrays to the checked out list */
958   if (ghosted) {
959     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
960       if (!dd->admfarrayghostedout[i]) {
961         dd->admfarrayghostedout[i] = *iptr ;
962         dd->admfstartghostedout[i] = iarray_start;
963         dd->ghostedtdof            = itdof;
964         break;
965       }
966     }
967   } else {
968     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
969       if (!dd->admfarrayout[i]) {
970         dd->admfarrayout[i] = *iptr ;
971         dd->admfstartout[i] = iarray_start;
972         dd->tdof            = itdof;
973         break;
974       }
975     }
976   }
977   if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Too many DMDA ADIC arrays obtained");
978   if (tdof)        *tdof = itdof;
979   if (array_start) *(void**)array_start = iarray_start;
980   PetscFunctionReturn(0);
981 }
982 
983 #undef __FUNCT__
984 #define __FUNCT__ "DMDAGetAdicMFArray4"
985 PetscErrorCode  DMDAGetAdicMFArray4(DM da,PetscBool  ghosted,void *vptr,void *array_start,PetscInt *tdof)
986 {
987   PetscErrorCode ierr;
988   PetscInt       j,i,xs,ys,xm,ym,itdof = 0;
989   char           *iarray_start;
990   void           **iptr = (void**)vptr;
991   DM_DA          *dd = (DM_DA*)da->data;
992 
993   PetscFunctionBegin;
994   PetscValidHeaderSpecific(da,DM_CLASSID,1);
995   if (ghosted) {
996     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
997       if (dd->admfarrayghostedin[i]) {
998         *iptr                     = dd->admfarrayghostedin[i];
999         iarray_start              = (char*)dd->admfstartghostedin[i];
1000         itdof                     = dd->ghostedtdof;
1001         dd->admfarrayghostedin[i] = PETSC_NULL;
1002         dd->admfstartghostedin[i] = PETSC_NULL;
1003 
1004         goto done;
1005       }
1006     }
1007     xs = dd->Xs;
1008     ys = dd->Ys;
1009     xm = dd->Xe-dd->Xs;
1010     ym = dd->Ye-dd->Ys;
1011   } else {
1012     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1013       if (dd->admfarrayin[i]) {
1014         *iptr              = dd->admfarrayin[i];
1015         iarray_start       = (char*)dd->admfstartin[i];
1016         itdof              = dd->tdof;
1017         dd->admfarrayin[i] = PETSC_NULL;
1018         dd->admfstartin[i] = PETSC_NULL;
1019 
1020         goto done;
1021       }
1022     }
1023     xs = dd->xs;
1024     ys = dd->ys;
1025     xm = dd->xe-dd->xs;
1026     ym = dd->ye-dd->ys;
1027   }
1028 
1029   switch (dd->dim) {
1030     case 2: {
1031       void **ptr;
1032       itdof = xm*ym;
1033 
1034       ierr  = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*5*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
1035 
1036       ptr  = (void**)(iarray_start + xm*ym*5*sizeof(PetscScalar) - ys*sizeof(void*));
1037       for(j=ys;j<ys+ym;j++) {
1038         ptr[j] = iarray_start + 5*sizeof(PetscScalar)*(xm*(j-ys) - xs);
1039       }
1040       *iptr = (void*)ptr;
1041       break;}
1042     default:
1043       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
1044   }
1045 
1046   done:
1047   /* add arrays to the checked out list */
1048   if (ghosted) {
1049     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1050       if (!dd->admfarrayghostedout[i]) {
1051         dd->admfarrayghostedout[i] = *iptr ;
1052         dd->admfstartghostedout[i] = iarray_start;
1053         dd->ghostedtdof            = itdof;
1054         break;
1055       }
1056     }
1057   } else {
1058     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1059       if (!dd->admfarrayout[i]) {
1060         dd->admfarrayout[i] = *iptr ;
1061         dd->admfstartout[i] = iarray_start;
1062         dd->tdof            = itdof;
1063         break;
1064       }
1065     }
1066   }
1067   if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Too many DMDA ADIC arrays obtained");
1068   if (tdof)        *tdof = itdof;
1069   if (array_start) *(void**)array_start = iarray_start;
1070   PetscFunctionReturn(0);
1071 }
1072 
1073 #undef __FUNCT__
1074 #define __FUNCT__ "DMDAGetAdicMFArray9"
1075 PetscErrorCode  DMDAGetAdicMFArray9(DM da,PetscBool  ghosted,void *vptr,void *array_start,PetscInt *tdof)
1076 {
1077   PetscErrorCode ierr;
1078   PetscInt       j,i,xs,ys,xm,ym,itdof = 0;
1079   char           *iarray_start;
1080   void           **iptr = (void**)vptr;
1081   DM_DA          *dd = (DM_DA*)da->data;
1082 
1083   PetscFunctionBegin;
1084   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1085   if (ghosted) {
1086     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1087       if (dd->admfarrayghostedin[i]) {
1088         *iptr                     = dd->admfarrayghostedin[i];
1089         iarray_start              = (char*)dd->admfstartghostedin[i];
1090         itdof                     = dd->ghostedtdof;
1091         dd->admfarrayghostedin[i] = PETSC_NULL;
1092         dd->admfstartghostedin[i] = PETSC_NULL;
1093 
1094         goto done;
1095       }
1096     }
1097     xs = dd->Xs;
1098     ys = dd->Ys;
1099     xm = dd->Xe-dd->Xs;
1100     ym = dd->Ye-dd->Ys;
1101   } else {
1102     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1103       if (dd->admfarrayin[i]) {
1104         *iptr              = dd->admfarrayin[i];
1105         iarray_start       = (char*)dd->admfstartin[i];
1106         itdof              = dd->tdof;
1107         dd->admfarrayin[i] = PETSC_NULL;
1108         dd->admfstartin[i] = PETSC_NULL;
1109 
1110         goto done;
1111       }
1112     }
1113     xs = dd->xs;
1114     ys = dd->ys;
1115     xm = dd->xe-dd->xs;
1116     ym = dd->ye-dd->ys;
1117   }
1118 
1119   switch (dd->dim) {
1120     case 2: {
1121       void **ptr;
1122       itdof = xm*ym;
1123 
1124       ierr  = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*10*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
1125 
1126       ptr  = (void**)(iarray_start + xm*ym*10*sizeof(PetscScalar) - ys*sizeof(void*));
1127       for(j=ys;j<ys+ym;j++) {
1128         ptr[j] = iarray_start + 10*sizeof(PetscScalar)*(xm*(j-ys) - xs);
1129       }
1130       *iptr = (void*)ptr;
1131       break;}
1132     default:
1133       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
1134   }
1135 
1136   done:
1137   /* add arrays to the checked out list */
1138   if (ghosted) {
1139     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1140       if (!dd->admfarrayghostedout[i]) {
1141         dd->admfarrayghostedout[i] = *iptr ;
1142         dd->admfstartghostedout[i] = iarray_start;
1143         dd->ghostedtdof            = itdof;
1144         break;
1145       }
1146     }
1147   } else {
1148     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1149       if (!dd->admfarrayout[i]) {
1150         dd->admfarrayout[i] = *iptr ;
1151         dd->admfstartout[i] = iarray_start;
1152         dd->tdof            = itdof;
1153         break;
1154       }
1155     }
1156   }
1157   if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Too many DMDA ADIC arrays obtained");
1158   if (tdof)        *tdof = itdof;
1159   if (array_start) *(void**)array_start = iarray_start;
1160   PetscFunctionReturn(0);
1161 }
1162 
1163 #undef __FUNCT__
1164 #define __FUNCT__ "DMDAGetAdicMFArrayb"
1165 /*@C
1166      DMDAGetAdicMFArrayb - Gets an array of derivative types for a DMDA for matrix-free ADIC.
1167 
1168      Input Parameter:
1169 +    da - information about my local patch
1170 -    ghosted - do you want arrays for the ghosted or nonghosted patch?
1171 
1172      Output Parameters:
1173 +    vptr - array data structured to be passed to ad_FormFunctionLocal()
1174 .    array_start - actual start of 1d array of all values that adiC can access directly (may be null)
1175 -    tdof - total number of degrees of freedom represented in array_start (may be null)
1176 
1177      Notes:
1178      The vector values are NOT initialized and may have garbage in them, so you may need
1179      to zero them.
1180 
1181      This routine returns the same type of object as the DMDAVecGetArray(), except its
1182      elements are derivative types instead of PetscScalars.
1183 
1184      Level: advanced
1185 
1186 .seealso: DMDARestoreAdicMFArray(), DMDAGetArray(), DMDAGetAdicArray()
1187 
1188 @*/
1189 PetscErrorCode  DMDAGetAdicMFArrayb(DM da,PetscBool  ghosted,void *vptr,void *array_start,PetscInt *tdof)
1190 {
1191   PetscErrorCode ierr;
1192   PetscInt       j,i,xs,ys,xm,ym,zs,zm,itdof = 0;
1193   char           *iarray_start;
1194   void           **iptr = (void**)vptr;
1195   DM_DA          *dd = (DM_DA*)da->data;
1196   PetscInt       bs = dd->w,bs1 = bs+1;
1197 
1198   PetscFunctionBegin;
1199   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1200   if (ghosted) {
1201     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1202       if (dd->admfarrayghostedin[i]) {
1203         *iptr                     = dd->admfarrayghostedin[i];
1204         iarray_start              = (char*)dd->admfstartghostedin[i];
1205         itdof                     = dd->ghostedtdof;
1206         dd->admfarrayghostedin[i] = PETSC_NULL;
1207         dd->admfstartghostedin[i] = PETSC_NULL;
1208 
1209         goto done;
1210       }
1211     }
1212     xs = dd->Xs;
1213     ys = dd->Ys;
1214     zs = dd->Zs;
1215     xm = dd->Xe-dd->Xs;
1216     ym = dd->Ye-dd->Ys;
1217     zm = dd->Ze-dd->Zs;
1218   } else {
1219     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1220       if (dd->admfarrayin[i]) {
1221         *iptr              = dd->admfarrayin[i];
1222         iarray_start       = (char*)dd->admfstartin[i];
1223         itdof              = dd->tdof;
1224         dd->admfarrayin[i] = PETSC_NULL;
1225         dd->admfstartin[i] = PETSC_NULL;
1226 
1227         goto done;
1228       }
1229     }
1230     xs = dd->xs;
1231     ys = dd->ys;
1232     zs = dd->zs;
1233     xm = dd->xe-dd->xs;
1234     ym = dd->ye-dd->ys;
1235     zm = dd->ze-dd->zs;
1236   }
1237 
1238   switch (dd->dim) {
1239     case 1: {
1240       void *ptr;
1241       itdof = xm;
1242 
1243       ierr  = PetscMalloc(xm*bs1*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
1244 
1245       ptr   = (void*)(iarray_start - xs*bs1*sizeof(PetscScalar));
1246       *iptr = (void*)ptr;
1247       break;}
1248     case 2: {
1249       void **ptr;
1250       itdof = xm*ym;
1251 
1252       ierr  = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*bs1*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
1253 
1254       ptr  = (void**)(iarray_start + xm*ym*bs1*sizeof(PetscScalar) - ys*sizeof(void*));
1255       for(j=ys;j<ys+ym;j++) {
1256         ptr[j] = iarray_start + bs1*sizeof(PetscScalar)*(xm*(j-ys) - xs);
1257       }
1258       *iptr = (void*)ptr;
1259       break;}
1260     case 3: {
1261       void ***ptr,**bptr;
1262       itdof = xm*ym*zm;
1263 
1264       ierr  = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*bs1*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
1265 
1266       ptr  = (void***)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) - zs*sizeof(void*));
1267       bptr = (void**)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) + zm*sizeof(void**));
1268       for(i=zs;i<zs+zm;i++) {
1269         ptr[i] = bptr + ((i-zs)*ym* - ys)*sizeof(void*);
1270       }
1271       for (i=zs; i<zs+zm; i++) {
1272         for (j=ys; j<ys+ym; j++) {
1273           ptr[i][j] = iarray_start + bs1*sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
1274         }
1275       }
1276 
1277       *iptr = (void*)ptr;
1278       break;}
1279     default:
1280       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
1281   }
1282 
1283   done:
1284   /* add arrays to the checked out list */
1285   if (ghosted) {
1286     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1287       if (!dd->admfarrayghostedout[i]) {
1288         dd->admfarrayghostedout[i] = *iptr ;
1289         dd->admfstartghostedout[i] = iarray_start;
1290         dd->ghostedtdof            = itdof;
1291         break;
1292       }
1293     }
1294   } else {
1295     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1296       if (!dd->admfarrayout[i]) {
1297         dd->admfarrayout[i] = *iptr ;
1298         dd->admfstartout[i] = iarray_start;
1299         dd->tdof            = itdof;
1300         break;
1301       }
1302     }
1303   }
1304   if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Too many DMDA ADIC arrays obtained");
1305   if (tdof)        *tdof = itdof;
1306   if (array_start) *(void**)array_start = iarray_start;
1307   PetscFunctionReturn(0);
1308 }
1309 
1310 #undef __FUNCT__
1311 #define __FUNCT__ "DMDARestoreAdicMFArray"
1312 /*@C
1313      DMDARestoreAdicMFArray - Restores an array of derivative types for a DMDA.
1314 
1315      Input Parameter:
1316 +    da - information about my local patch
1317 -    ghosted - do you want arrays for the ghosted or nonghosted patch?
1318 
1319      Output Parameters:
1320 +    ptr - array data structure to be passed to ad_FormFunctionLocal()
1321 .    array_start - actual start of 1d array of all values that adiC can access directly
1322 -    tdof - total number of degrees of freedom represented in array_start
1323 
1324      Level: advanced
1325 
1326 .seealso: DMDAGetAdicArray()
1327 
1328 @*/
1329 PetscErrorCode  DMDARestoreAdicMFArray(DM da,PetscBool  ghosted,void *vptr,void *array_start,PetscInt *tdof)
1330 {
1331   PetscInt  i;
1332   void      **iptr = (void**)vptr,*iarray_start = 0;
1333   DM_DA     *dd = (DM_DA*)da->data;
1334 
1335   PetscFunctionBegin;
1336   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1337   if (ghosted) {
1338     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1339       if (dd->admfarrayghostedout[i] == *iptr) {
1340         iarray_start               = dd->admfstartghostedout[i];
1341         dd->admfarrayghostedout[i] = PETSC_NULL;
1342         dd->admfstartghostedout[i] = PETSC_NULL;
1343         break;
1344       }
1345     }
1346     if (!iarray_start) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Could not find array in checkout list");
1347     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1348       if (!dd->admfarrayghostedin[i]){
1349         dd->admfarrayghostedin[i] = *iptr;
1350         dd->admfstartghostedin[i] = iarray_start;
1351         break;
1352       }
1353     }
1354   } else {
1355     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1356       if (dd->admfarrayout[i] == *iptr) {
1357         iarray_start        = dd->admfstartout[i];
1358         dd->admfarrayout[i] = PETSC_NULL;
1359         dd->admfstartout[i] = PETSC_NULL;
1360         break;
1361       }
1362     }
1363     if (!iarray_start) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Could not find array in checkout list");
1364     for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) {
1365       if (!dd->admfarrayin[i]){
1366         dd->admfarrayin[i] = *iptr;
1367         dd->admfstartin[i] = iarray_start;
1368         break;
1369       }
1370     }
1371   }
1372   PetscFunctionReturn(0);
1373 }
1374 
1375 #undef __FUNCT__
1376 #define __FUNCT__ "admf_DAGetArray"
1377 PetscErrorCode  admf_DAGetArray(DM da,PetscBool  ghosted,void *iptr)
1378 {
1379   PetscErrorCode ierr;
1380   PetscFunctionBegin;
1381   ierr = DMDAGetAdicMFArray(da,ghosted,iptr,0,0);CHKERRQ(ierr);
1382   PetscFunctionReturn(0);
1383 }
1384 
1385 #undef __FUNCT__
1386 #define __FUNCT__ "admf_DARestoreArray"
1387 PetscErrorCode  admf_DARestoreArray(DM da,PetscBool  ghosted,void *iptr)
1388 {
1389   PetscErrorCode ierr;
1390   PetscFunctionBegin;
1391   ierr = DMDARestoreAdicMFArray(da,ghosted,iptr,0,0);CHKERRQ(ierr);
1392   PetscFunctionReturn(0);
1393 }
1394 
1395