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