Actual source code: dalocal.c
1: /*
2: Code for manipulating distributed regular arrays in parallel.
3: */
5: #include src/dm/da/daimpl.h
7: /*
8: This allows the DA vectors to properly tell Matlab their dimensions
9: */
10: #if defined(PETSC_HAVE_MATLAB) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
11: #include "engine.h" /* Matlab include file */
12: #include "mex.h" /* Matlab include file */
16: PetscErrorCode VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine)
17: {
19: PetscInt n,m;
20: Vec vec = (Vec)obj;
21: PetscScalar *array;
22: mxArray *mat;
23: DA da;
26: PetscObjectQuery((PetscObject)vec,"DA",(PetscObject*)&da);
27: if (!da) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DA");
28: DAGetGhostCorners(da,0,0,0,&m,&n,0);
30: VecGetArray(vec,&array);
31: #if !defined(PETSC_USE_COMPLEX)
32: mat = mxCreateDoubleMatrix(m,n,mxREAL);
33: #else
34: mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX);
35: #endif
36: PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));
37: PetscObjectName(obj);
38: engPutVariable((Engine *)mengine,obj->name,mat);
39:
40: VecRestoreArray(vec,&array);
41: return(0);
42: }
44: #endif
49: /*@C
50: DACreateLocalVector - Creates a Seq PETSc vector that
51: may be used with the DAXXX routines.
53: Not Collective
55: Input Parameter:
56: . da - the distributed array
58: Output Parameter:
59: . g - the local vector
61: Level: beginner
63: Note:
64: The output parameter, g, is a regular PETSc vector that should be destroyed
65: with a call to VecDestroy() when usage is finished.
67: .keywords: distributed array, create, local, vector
69: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
70: DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToLocalBegin(),
71: DAGlobalToLocalEnd(), DALocalToGlobal(), DAGetLocalVector(), DARestoreLocalVector()
72: @*/
73: PetscErrorCode DACreateLocalVector(DA da,Vec* g)
74: {
80: VecCreateSeq(PETSC_COMM_SELF,da->nlocal,g);
81: VecSetBlockSize(*g,da->w);
82: PetscObjectCompose((PetscObject)*g,"DA",(PetscObject)da);
83: #if defined(PETSC_HAVE_MATLAB) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
84: if (da->w == 1 && da->dim == 2) {
85: PetscObjectComposeFunctionDynamic((PetscObject)*g,"PetscMatlabEnginePut_C","VecMatlabEnginePut_DA2d",VecMatlabEnginePut_DA2d);
86: }
87: #endif
88: return(0);
89: }
93: /*@C
94: DAGetLocalVector - Gets a Seq PETSc vector that
95: may be used with the DAXXX routines.
97: Not Collective
99: Input Parameter:
100: . da - the distributed array
102: Output Parameter:
103: . g - the local vector
105: Level: beginner
107: Note:
108: The vector values are NOT initialized and may have garbage in them, so you may need
109: to zero them.
111: The output parameter, g, is a regular PETSc vector that should be returned with
112: DARestoreLocalVector() DO NOT call VecDestroy() on it.
114: VecStride*() operations can be useful when using DA with dof > 1
116: .keywords: distributed array, create, local, vector
118: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
119: DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToLocalBegin(),
120: DAGlobalToLocalEnd(), DALocalToGlobal(), DACreateLocalVector(), DARestoreLocalVector(),
121: VecStrideMax(), VecStrideMin(), VecStrideNorm()
122: @*/
123: PetscErrorCode DAGetLocalVector(DA da,Vec* g)
124: {
125: PetscErrorCode ierr,i;
130: for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
131: if (da->localin[i]) {
132: *g = da->localin[i];
133: da->localin[i] = PETSC_NULL;
134: goto alldone;
135: }
136: }
137: DACreateLocalVector(da,g);
139: alldone:
140: for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
141: if (!da->localout[i]) {
142: da->localout[i] = *g;
143: break;
144: }
145: }
146: return(0);
147: }
151: /*@C
152: DARestoreLocalVector - Returns a Seq PETSc vector that
153: obtained from DAGetLocalVector(). Do not use with vector obtained via
154: DACreateLocalVector().
156: Not Collective
158: Input Parameter:
159: + da - the distributed array
160: - g - the local vector
162: Level: beginner
164: .keywords: distributed array, create, local, vector
166: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
167: DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToLocalBegin(),
168: DAGlobalToLocalEnd(), DALocalToGlobal(), DACreateLocalVector(), DAGetLocalVector()
169: @*/
170: PetscErrorCode DARestoreLocalVector(DA da,Vec* g)
171: {
173: PetscInt i,j;
178: for (j=0; j<DA_MAX_WORK_VECTORS; j++) {
179: if (*g == da->localout[j]) {
180: da->localout[j] = PETSC_NULL;
181: for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
182: if (!da->localin[i]) {
183: da->localin[i] = *g;
184: goto alldone;
185: }
186: }
187: }
188: }
189: VecDestroy(*g);
190: alldone:
191: return(0);
192: }
196: /*@C
197: DAGetGlobalVector - Gets a MPI PETSc vector that
198: may be used with the DAXXX routines.
200: Collective on DA
202: Input Parameter:
203: . da - the distributed array
205: Output Parameter:
206: . g - the global vector
208: Level: beginner
210: Note:
211: The vector values are NOT initialized and may have garbage in them, so you may need
212: to zero them.
214: The output parameter, g, is a regular PETSc vector that should be returned with
215: DARestoreGlobalVector() DO NOT call VecDestroy() on it.
217: VecStride*() operations can be useful when using DA with dof > 1
219: .keywords: distributed array, create, Global, vector
221: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
222: DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToLocalBegin(),
223: DAGlobalToLocalEnd(), DALocalToGlobal(), DACreateLocalVector(), DARestoreLocalVector()
224: VecStrideMax(), VecStrideMin(), VecStrideNorm()
226: @*/
227: PetscErrorCode DAGetGlobalVector(DA da,Vec* g)
228: {
230: PetscInt i;
235: for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
236: if (da->globalin[i]) {
237: *g = da->globalin[i];
238: da->globalin[i] = PETSC_NULL;
239: goto alldone;
240: }
241: }
242: DACreateGlobalVector(da,g);
244: alldone:
245: for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
246: if (!da->globalout[i]) {
247: da->globalout[i] = *g;
248: break;
249: }
250: }
251: return(0);
252: }
256: /*@C
257: DARestoreGlobalVector - Returns a Seq PETSc vector that
258: obtained from DAGetGlobalVector(). Do not use with vector obtained via
259: DACreateGlobalVector().
261: Not Collective
263: Input Parameter:
264: + da - the distributed array
265: - g - the global vector
267: Level: beginner
269: .keywords: distributed array, create, global, vector
271: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
272: DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToGlobalBegin(),
273: DAGlobalToGlobalEnd(), DAGlobalToGlobal(), DACreateLocalVector(), DAGetGlobalVector()
274: @*/
275: PetscErrorCode DARestoreGlobalVector(DA da,Vec* g)
276: {
278: PetscInt i,j;
283: for (j=0; j<DA_MAX_WORK_VECTORS; j++) {
284: if (*g == da->globalout[j]) {
285: da->globalout[j] = PETSC_NULL;
286: for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
287: if (!da->globalin[i]) {
288: da->globalin[i] = *g;
289: goto alldone;
290: }
291: }
292: }
293: }
294: VecDestroy(*g);
295: alldone:
296: return(0);
297: }
299: /* ------------------------------------------------------------------- */
300: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
303: #include "adic/ad_utils.h"
308: /*@C
309: DAGetAdicArray - Gets an array of derivative types for a DA
310:
311: Input Parameter:
312: + da - information about my local patch
313: - ghosted - do you want arrays for the ghosted or nonghosted patch
315: Output Parameters:
316: + ptr - array data structured to be passed to ad_FormFunctionLocal()
317: . array_start - actual start of 1d array of all values that adiC can access directly (may be null)
318: - tdof - total number of degrees of freedom represented in array_start (may be null)
320: Notes:
321: The vector values are NOT initialized and may have garbage in them, so you may need
322: to zero them.
324: Returns the same type of object as the DAVecGetArray() except its elements are
325: derivative types instead of PetscScalars
327: Level: advanced
329: .seealso: DARestoreAdicArray()
331: @*/
332: PetscErrorCode DAGetAdicArray(DA da,PetscTruth ghosted,void **iptr,void **array_start,PetscInt *tdof)
333: {
335: PetscInt j,i,deriv_type_size,xs,ys,xm,ym,zs,zm,itdof;
336: char *iarray_start;
340: if (ghosted) {
341: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
342: if (da->adarrayghostedin[i]) {
343: *iptr = da->adarrayghostedin[i];
344: iarray_start = (char*)da->adstartghostedin[i];
345: itdof = da->ghostedtdof;
346: da->adarrayghostedin[i] = PETSC_NULL;
347: da->adstartghostedin[i] = PETSC_NULL;
348:
349: goto done;
350: }
351: }
352: xs = da->Xs;
353: ys = da->Ys;
354: zs = da->Zs;
355: xm = da->Xe-da->Xs;
356: ym = da->Ye-da->Ys;
357: zm = da->Ze-da->Zs;
358: } else {
359: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
360: if (da->adarrayin[i]) {
361: *iptr = da->adarrayin[i];
362: iarray_start = (char*)da->adstartin[i];
363: itdof = da->tdof;
364: da->adarrayin[i] = PETSC_NULL;
365: da->adstartin[i] = PETSC_NULL;
366:
367: goto done;
368: }
369: }
370: xs = da->xs;
371: ys = da->ys;
372: zs = da->zs;
373: xm = da->xe-da->xs;
374: ym = da->ye-da->ys;
375: zm = da->ze-da->zs;
376: }
377: deriv_type_size = PetscADGetDerivTypeSize();
379: switch (da->dim) {
380: case 1: {
381: void *ptr;
382: itdof = xm;
384: PetscMalloc(xm*deriv_type_size,&iarray_start);
386: ptr = (void*)(iarray_start - xs*deriv_type_size);
387: *iptr = (void*)ptr;
388: break;}
389: case 2: {
390: void **ptr;
391: itdof = xm*ym;
393: PetscMalloc((ym+1)*sizeof(void*)+xm*ym*deriv_type_size,&iarray_start);
395: ptr = (void**)(iarray_start + xm*ym*deriv_type_size - ys*sizeof(void*));
396: for(j=ys;j<ys+ym;j++) {
397: ptr[j] = iarray_start + deriv_type_size*(xm*(j-ys) - xs);
398: }
399: *iptr = (void*)ptr;
400: break;}
401: case 3: {
402: void ***ptr,**bptr;
403: itdof = xm*ym*zm;
405: PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*deriv_type_size,&iarray_start);
407: ptr = (void***)(iarray_start + xm*ym*zm*deriv_type_size - zs*sizeof(void*));
408: bptr = (void**)(iarray_start + xm*ym*zm*deriv_type_size + zm*sizeof(void**));
410: for(i=zs;i<zs+zm;i++) {
411: ptr[i] = bptr + ((i-zs)*ym - ys);
412: }
413: for (i=zs; i<zs+zm; i++) {
414: for (j=ys; j<ys+ym; j++) {
415: ptr[i][j] = iarray_start + deriv_type_size*(xm*ym*(i-zs) + xm*(j-ys) - xs);
416: }
417: }
419: *iptr = (void*)ptr;
420: break;}
421: default:
422: SETERRQ1(PETSC_ERR_SUP,"Dimension %D not supported",da->dim);
423: }
425: done:
426: /* add arrays to the checked out list */
427: if (ghosted) {
428: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
429: if (!da->adarrayghostedout[i]) {
430: da->adarrayghostedout[i] = *iptr ;
431: da->adstartghostedout[i] = iarray_start;
432: da->ghostedtdof = itdof;
433: break;
434: }
435: }
436: } else {
437: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
438: if (!da->adarrayout[i]) {
439: da->adarrayout[i] = *iptr ;
440: da->adstartout[i] = iarray_start;
441: da->tdof = itdof;
442: break;
443: }
444: }
445: }
446: if (i == DA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_ERR_SUP,"Too many DA ADIC arrays obtained");
447: if (tdof) *tdof = itdof;
448: if (array_start) *array_start = iarray_start;
449: return(0);
450: }
454: /*@C
455: DARestoreAdicArray - Restores an array of derivative types for a DA
456:
457: Input Parameter:
458: + da - information about my local patch
459: - ghosted - do you want arrays for the ghosted or nonghosted patch
461: Output Parameters:
462: + ptr - array data structured to be passed to ad_FormFunctionLocal()
463: . array_start - actual start of 1d array of all values that adiC can access directly
464: - tdof - total number of degrees of freedom represented in array_start
466: Level: advanced
468: .seealso: DAGetAdicArray()
470: @*/
471: PetscErrorCode DARestoreAdicArray(DA da,PetscTruth ghosted,void **iptr,void **array_start,PetscInt *tdof)
472: {
473: PetscInt i;
474: void *iarray_start = 0;
475:
478: if (ghosted) {
479: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
480: if (da->adarrayghostedout[i] == *iptr) {
481: iarray_start = da->adstartghostedout[i];
482: da->adarrayghostedout[i] = PETSC_NULL;
483: da->adstartghostedout[i] = PETSC_NULL;
484: break;
485: }
486: }
487: if (!iarray_start) SETERRQ(PETSC_ERR_ARG_WRONG,"Could not find array in checkout list");
488: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
489: if (!da->adarrayghostedin[i]){
490: da->adarrayghostedin[i] = *iptr;
491: da->adstartghostedin[i] = iarray_start;
492: break;
493: }
494: }
495: } else {
496: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
497: if (da->adarrayout[i] == *iptr) {
498: iarray_start = da->adstartout[i];
499: da->adarrayout[i] = PETSC_NULL;
500: da->adstartout[i] = PETSC_NULL;
501: break;
502: }
503: }
504: if (!iarray_start) SETERRQ(PETSC_ERR_ARG_WRONG,"Could not find array in checkout list");
505: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
506: if (!da->adarrayin[i]){
507: da->adarrayin[i] = *iptr;
508: da->adstartin[i] = iarray_start;
509: break;
510: }
511: }
512: }
513: return(0);
514: }
518: PetscErrorCode ad_DAGetArray(DA da,PetscTruth ghosted,void **iptr)
519: {
522: DAGetAdicArray(da,ghosted,iptr,0,0);
523: return(0);
524: }
528: PetscErrorCode ad_DARestoreArray(DA da,PetscTruth ghosted,void **iptr)
529: {
532: DARestoreAdicArray(da,ghosted,iptr,0,0);
533: return(0);
534: }
536: #endif
540: /*@C
541: DAGetArray - Gets a work array for a DA
542:
543: Input Parameter:
544: + da - information about my local patch
545: - ghosted - do you want arrays for the ghosted or nonghosted patch
547: Output Parameters:
548: . ptr - array data structured
550: Note: The vector values are NOT initialized and may have garbage in them, so you may need
551: to zero them.
553: Level: advanced
555: .seealso: DARestoreArray(), DAGetAdicArray()
557: @*/
558: PetscErrorCode DAGetArray(DA da,PetscTruth ghosted,void **iptr)
559: {
561: PetscInt j,i,xs,ys,xm,ym,zs,zm;
562: char *iarray_start;
566: if (ghosted) {
567: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
568: if (da->arrayghostedin[i]) {
569: *iptr = da->arrayghostedin[i];
570: iarray_start = (char*)da->startghostedin[i];
571: da->arrayghostedin[i] = PETSC_NULL;
572: da->startghostedin[i] = PETSC_NULL;
573:
574: goto done;
575: }
576: }
577: xs = da->Xs;
578: ys = da->Ys;
579: zs = da->Zs;
580: xm = da->Xe-da->Xs;
581: ym = da->Ye-da->Ys;
582: zm = da->Ze-da->Zs;
583: } else {
584: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
585: if (da->arrayin[i]) {
586: *iptr = da->arrayin[i];
587: iarray_start = (char*)da->startin[i];
588: da->arrayin[i] = PETSC_NULL;
589: da->startin[i] = PETSC_NULL;
590:
591: goto done;
592: }
593: }
594: xs = da->xs;
595: ys = da->ys;
596: zs = da->zs;
597: xm = da->xe-da->xs;
598: ym = da->ye-da->ys;
599: zm = da->ze-da->zs;
600: }
602: switch (da->dim) {
603: case 1: {
604: void *ptr;
606: PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);
608: ptr = (void*)(iarray_start - xs*sizeof(PetscScalar));
609: *iptr = (void*)ptr;
610: break;}
611: case 2: {
612: void **ptr;
614: PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);
616: ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*));
617: for(j=ys;j<ys+ym;j++) {
618: ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs);
619: }
620: *iptr = (void*)ptr;
621: break;}
622: case 3: {
623: void ***ptr,**bptr;
625: PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);
627: ptr = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*));
628: bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**));
629: for(i=zs;i<zs+zm;i++) {
630: ptr[i] = bptr + ((i-zs)*ym - ys);
631: }
632: for (i=zs; i<zs+zm; i++) {
633: for (j=ys; j<ys+ym; j++) {
634: ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
635: }
636: }
638: *iptr = (void*)ptr;
639: break;}
640: default:
641: SETERRQ1(PETSC_ERR_SUP,"Dimension %D not supported",da->dim);
642: }
644: done:
645: /* add arrays to the checked out list */
646: if (ghosted) {
647: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
648: if (!da->arrayghostedout[i]) {
649: da->arrayghostedout[i] = *iptr ;
650: da->startghostedout[i] = iarray_start;
651: break;
652: }
653: }
654: } else {
655: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
656: if (!da->arrayout[i]) {
657: da->arrayout[i] = *iptr ;
658: da->startout[i] = iarray_start;
659: break;
660: }
661: }
662: }
663: return(0);
664: }
668: /*@C
669: DARestoreArray - Restores an array of derivative types for a DA
670:
671: Input Parameter:
672: + da - information about my local patch
673: . ghosted - do you want arrays for the ghosted or nonghosted patch
674: - ptr - array data structured to be passed to ad_FormFunctionLocal()
676: Level: advanced
678: .seealso: DAGetArray(), DAGetAdicArray()
680: @*/
681: PetscErrorCode DARestoreArray(DA da,PetscTruth ghosted,void **iptr)
682: {
683: PetscInt i;
684: void *iarray_start = 0;
685:
688: if (ghosted) {
689: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
690: if (da->arrayghostedout[i] == *iptr) {
691: iarray_start = da->startghostedout[i];
692: da->arrayghostedout[i] = PETSC_NULL;
693: da->startghostedout[i] = PETSC_NULL;
694: break;
695: }
696: }
697: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
698: if (!da->arrayghostedin[i]){
699: da->arrayghostedin[i] = *iptr;
700: da->startghostedin[i] = iarray_start;
701: break;
702: }
703: }
704: } else {
705: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
706: if (da->arrayout[i] == *iptr) {
707: iarray_start = da->startout[i];
708: da->arrayout[i] = PETSC_NULL;
709: da->startout[i] = PETSC_NULL;
710: break;
711: }
712: }
713: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
714: if (!da->arrayin[i]){
715: da->arrayin[i] = *iptr;
716: da->startin[i] = iarray_start;
717: break;
718: }
719: }
720: }
721: return(0);
722: }
726: /*@C
727: DAGetAdicMFArray - Gets an array of derivative types for a DA for matrix-free ADIC.
728:
729: Input Parameter:
730: + da - information about my local patch
731: - ghosted - do you want arrays for the ghosted or nonghosted patch?
733: Output Parameters:
734: + iptr - array data structured to be passed to ad_FormFunctionLocal()
735: . array_start - actual start of 1d array of all values that adiC can access directly (may be null)
736: - tdof - total number of degrees of freedom represented in array_start (may be null)
738: Notes:
739: The vector values are NOT initialized and may have garbage in them, so you may need
740: to zero them.
742: This routine returns the same type of object as the DAVecGetArray(), except its
743: elements are derivative types instead of PetscScalars.
745: Level: advanced
747: .seealso: DARestoreAdicMFArray(), DAGetArray(), DAGetAdicArray()
749: @*/
750: PetscErrorCode DAGetAdicMFArray(DA da,PetscTruth ghosted,void **iptr,void **array_start,PetscInt *tdof)
751: {
753: PetscInt j,i,xs,ys,xm,ym,zs,zm,itdof;
754: char *iarray_start;
758: if (ghosted) {
759: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
760: if (da->admfarrayghostedin[i]) {
761: *iptr = da->admfarrayghostedin[i];
762: iarray_start = (char*)da->admfstartghostedin[i];
763: itdof = da->ghostedtdof;
764: da->admfarrayghostedin[i] = PETSC_NULL;
765: da->admfstartghostedin[i] = PETSC_NULL;
766:
767: goto done;
768: }
769: }
770: xs = da->Xs;
771: ys = da->Ys;
772: zs = da->Zs;
773: xm = da->Xe-da->Xs;
774: ym = da->Ye-da->Ys;
775: zm = da->Ze-da->Zs;
776: } else {
777: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
778: if (da->admfarrayin[i]) {
779: *iptr = da->admfarrayin[i];
780: iarray_start = (char*)da->admfstartin[i];
781: itdof = da->tdof;
782: da->admfarrayin[i] = PETSC_NULL;
783: da->admfstartin[i] = PETSC_NULL;
784:
785: goto done;
786: }
787: }
788: xs = da->xs;
789: ys = da->ys;
790: zs = da->zs;
791: xm = da->xe-da->xs;
792: ym = da->ye-da->ys;
793: zm = da->ze-da->zs;
794: }
796: switch (da->dim) {
797: case 1: {
798: void *ptr;
799: itdof = xm;
801: PetscMalloc(xm*2*sizeof(PetscScalar),&iarray_start);
803: ptr = (void*)(iarray_start - xs*2*sizeof(PetscScalar));
804: *iptr = (void*)ptr;
805: break;}
806: case 2: {
807: void **ptr;
808: itdof = xm*ym;
810: PetscMalloc((ym+1)*sizeof(void*)+xm*ym*2*sizeof(PetscScalar),&iarray_start);
812: ptr = (void**)(iarray_start + xm*ym*2*sizeof(PetscScalar) - ys*sizeof(void*));
813: for(j=ys;j<ys+ym;j++) {
814: ptr[j] = iarray_start + 2*sizeof(PetscScalar)*(xm*(j-ys) - xs);
815: }
816: *iptr = (void*)ptr;
817: break;}
818: case 3: {
819: void ***ptr,**bptr;
820: itdof = xm*ym*zm;
822: PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*2*sizeof(PetscScalar),&iarray_start);
824: ptr = (void***)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) - zs*sizeof(void*));
825: bptr = (void**)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) + zm*sizeof(void**));
826: for(i=zs;i<zs+zm;i++) {
827: ptr[i] = bptr + ((i-zs)*ym* - ys)*sizeof(void*);
828: }
829: for (i=zs; i<zs+zm; i++) {
830: for (j=ys; j<ys+ym; j++) {
831: ptr[i][j] = iarray_start + 2*sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
832: }
833: }
835: *iptr = (void*)ptr;
836: break;}
837: default:
838: SETERRQ1(PETSC_ERR_SUP,"Dimension %D not supported",da->dim);
839: }
841: done:
842: /* add arrays to the checked out list */
843: if (ghosted) {
844: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
845: if (!da->admfarrayghostedout[i]) {
846: da->admfarrayghostedout[i] = *iptr ;
847: da->admfstartghostedout[i] = iarray_start;
848: da->ghostedtdof = itdof;
849: break;
850: }
851: }
852: } else {
853: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
854: if (!da->admfarrayout[i]) {
855: da->admfarrayout[i] = *iptr ;
856: da->admfstartout[i] = iarray_start;
857: da->tdof = itdof;
858: break;
859: }
860: }
861: }
862: if (i == DA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_ERR_ARG_WRONG,"Too many DA ADIC arrays obtained");
863: if (tdof) *tdof = itdof;
864: if (array_start) *array_start = iarray_start;
865: return(0);
866: }
870: /*@C
871: DARestoreAdicMFArray - Restores an array of derivative types for a DA.
872:
873: Input Parameter:
874: + da - information about my local patch
875: - ghosted - do you want arrays for the ghosted or nonghosted patch?
877: Output Parameters:
878: + ptr - array data structure to be passed to ad_FormFunctionLocal()
879: . array_start - actual start of 1d array of all values that adiC can access directly
880: - tdof - total number of degrees of freedom represented in array_start
882: Level: advanced
884: .seealso: DAGetAdicArray()
886: @*/
887: PetscErrorCode DARestoreAdicMFArray(DA da,PetscTruth ghosted,void **iptr,void **array_start,PetscInt *tdof)
888: {
889: PetscInt i;
890: void *iarray_start = 0;
891:
894: if (ghosted) {
895: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
896: if (da->admfarrayghostedout[i] == *iptr) {
897: iarray_start = da->admfstartghostedout[i];
898: da->admfarrayghostedout[i] = PETSC_NULL;
899: da->admfstartghostedout[i] = PETSC_NULL;
900: break;
901: }
902: }
903: if (!iarray_start) SETERRQ(PETSC_ERR_ARG_WRONG,"Could not find array in checkout list");
904: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
905: if (!da->admfarrayghostedin[i]){
906: da->admfarrayghostedin[i] = *iptr;
907: da->admfstartghostedin[i] = iarray_start;
908: break;
909: }
910: }
911: } else {
912: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
913: if (da->admfarrayout[i] == *iptr) {
914: iarray_start = da->admfstartout[i];
915: da->admfarrayout[i] = PETSC_NULL;
916: da->admfstartout[i] = PETSC_NULL;
917: break;
918: }
919: }
920: if (!iarray_start) SETERRQ(PETSC_ERR_ARG_WRONG,"Could not find array in checkout list");
921: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
922: if (!da->admfarrayin[i]){
923: da->admfarrayin[i] = *iptr;
924: da->admfstartin[i] = iarray_start;
925: break;
926: }
927: }
928: }
929: return(0);
930: }
934: PetscErrorCode admf_DAGetArray(DA da,PetscTruth ghosted,void **iptr)
935: {
938: DAGetAdicMFArray(da,ghosted,iptr,0,0);
939: return(0);
940: }
944: PetscErrorCode admf_DARestoreArray(DA da,PetscTruth ghosted,void **iptr)
945: {
948: DARestoreAdicMFArray(da,ghosted,iptr,0,0);
949: return(0);
950: }