Actual source code: dalocal.c
1: /*$Id: dalocal.c,v 1.35 2001/08/07 03:04:39 balay Exp $*/
2:
3: /*
4: Code for manipulating distributed regular arrays in parallel.
5: */
7: #include src/dm/da/daimpl.h
9: /*
10: This allows the DA vectors to properly tell Matlab their dimensions
11: */
12: #if defined(PETSC_HAVE_MATLAB) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
13: #include "engine.h" /* Matlab include file */
14: #include "mex.h" /* Matlab include file */
15: EXTERN_C_BEGIN
18: int VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine)
19: {
20: int ierr,n,m;
21: Vec vec = (Vec)obj;
22: PetscScalar *array;
23: mxArray *mat;
24: DA da;
27: PetscObjectQuery((PetscObject)vec,"DA",(PetscObject*)&da);
28: if (!da) SETERRQ(1,"Vector not associated with a DA");
29: DAGetGhostCorners(da,0,0,0,&m,&n,0);
31: VecGetArray(vec,&array);
32: #if !defined(PETSC_USE_COMPLEX)
33: mat = mxCreateDoubleMatrix(m,n,mxREAL);
34: #else
35: mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX);
36: #endif
37: PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));
38: PetscObjectName(obj);
39: engPutVariable((Engine *)mengine,obj->name,mat);
40:
41: VecRestoreArray(vec,&array);
42: return(0);
43: }
44: EXTERN_C_END
45: #endif
50: /*@C
51: DACreateLocalVector - Creates a Seq PETSc vector that
52: may be used with the DAXXX routines.
54: Not Collective
56: Input Parameter:
57: . da - the distributed array
59: Output Parameter:
60: . g - the local vector
62: Level: beginner
64: Note:
65: The output parameter, g, is a regular PETSc vector that should be destroyed
66: with a call to VecDestroy() when usage is finished.
68: .keywords: distributed array, create, local, vector
70: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
71: DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToLocalBegin(),
72: DAGlobalToLocalEnd(), DALocalToGlobal(), DAGetLocalVector(), DARestoreLocalVector()
73: @*/
74: int DACreateLocalVector(DA da,Vec* g)
75: {
81: VecCreateSeq(PETSC_COMM_SELF,da->nlocal,g);
82: VecSetBlockSize(*g,da->w);
83: PetscObjectCompose((PetscObject)*g,"DA",(PetscObject)da);
84: #if defined(PETSC_HAVE_MATLAB) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
85: if (da->w == 1 && da->dim == 2) {
86: PetscObjectComposeFunctionDynamic((PetscObject)*g,"PetscMatlabEnginePut_C","VecMatlabEnginePut_DA2d",VecMatlabEnginePut_DA2d);
87: }
88: #endif
89: return(0);
90: }
94: /*@C
95: DAGetLocalVector - Gets a Seq PETSc vector that
96: may be used with the DAXXX routines.
98: Not Collective
100: Input Parameter:
101: . da - the distributed array
103: Output Parameter:
104: . g - the local vector
106: Level: beginner
108: Note:
109: The vector values are NOT initialized and may have garbage in them, so you may need
110: to zero them.
112: The output parameter, g, is a regular PETSc vector that should be returned with
113: DARestoreLocalVector() DO NOT call VecDestroy() on it.
115: VecStride*() operations can be useful when using DA with dof > 1
117: .keywords: distributed array, create, local, vector
119: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
120: DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToLocalBegin(),
121: DAGlobalToLocalEnd(), DALocalToGlobal(), DACreateLocalVector(), DARestoreLocalVector(),
122: VecStrideMax(), VecStrideMin(), VecStrideNorm()
123: @*/
124: int DAGetLocalVector(DA da,Vec* g)
125: {
126: int ierr,i;
131: for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
132: if (da->localin[i]) {
133: *g = da->localin[i];
134: da->localin[i] = PETSC_NULL;
135: goto alldone;
136: }
137: }
138: DACreateLocalVector(da,g);
140: alldone:
141: for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
142: if (!da->localout[i]) {
143: da->localout[i] = *g;
144: break;
145: }
146: }
147: return(0);
148: }
152: /*@C
153: DARestoreLocalVector - Returns a Seq PETSc vector that
154: obtained from DAGetLocalVector(). Do not use with vector obtained via
155: DACreateLocalVector().
157: Not Collective
159: Input Parameter:
160: + da - the distributed array
161: - g - the local vector
163: Level: beginner
165: .keywords: distributed array, create, local, vector
167: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
168: DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToLocalBegin(),
169: DAGlobalToLocalEnd(), DALocalToGlobal(), DACreateLocalVector(), DAGetLocalVector()
170: @*/
171: int DARestoreLocalVector(DA da,Vec* g)
172: {
173: int ierr,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: int DAGetGlobalVector(DA da,Vec* g)
228: {
229: int ierr,i;
234: for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
235: if (da->globalin[i]) {
236: *g = da->globalin[i];
237: da->globalin[i] = PETSC_NULL;
238: goto alldone;
239: }
240: }
241: DACreateGlobalVector(da,g);
243: alldone:
244: for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
245: if (!da->globalout[i]) {
246: da->globalout[i] = *g;
247: break;
248: }
249: }
250: return(0);
251: }
255: /*@C
256: DARestoreGlobalVector - Returns a Seq PETSc vector that
257: obtained from DAGetGlobalVector(). Do not use with vector obtained via
258: DACreateGlobalVector().
260: Not Collective
262: Input Parameter:
263: + da - the distributed array
264: - g - the global vector
266: Level: beginner
268: .keywords: distributed array, create, global, vector
270: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
271: DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToGlobalBegin(),
272: DAGlobalToGlobalEnd(), DAGlobalToGlobal(), DACreateLocalVector(), DAGetGlobalVector()
273: @*/
274: int DARestoreGlobalVector(DA da,Vec* g)
275: {
276: int ierr,i,j;
281: for (j=0; j<DA_MAX_WORK_VECTORS; j++) {
282: if (*g == da->globalout[j]) {
283: da->globalout[j] = PETSC_NULL;
284: for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
285: if (!da->globalin[i]) {
286: da->globalin[i] = *g;
287: goto alldone;
288: }
289: }
290: }
291: }
292: VecDestroy(*g);
293: alldone:
294: return(0);
295: }
297: /* ------------------------------------------------------------------- */
298: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
300: EXTERN_C_BEGIN
301: #include "adic/ad_utils.h"
302: EXTERN_C_END
306: /*@C
307: DAGetAdicArray - Gets an array of derivative types for a DA
308:
309: Input Parameter:
310: + da - information about my local patch
311: - ghosted - do you want arrays for the ghosted or nonghosted patch
313: Output Parameters:
314: + ptr - array data structured to be passed to ad_FormFunctionLocal()
315: . array_start - actual start of 1d array of all values that adiC can access directly (may be null)
316: - tdof - total number of degrees of freedom represented in array_start (may be null)
318: Notes:
319: The vector values are NOT initialized and may have garbage in them, so you may need
320: to zero them.
322: Returns the same type of object as the DAVecGetArray() except its elements are
323: derivative types instead of PetscScalars
325: Level: advanced
327: .seealso: DARestoreAdicArray()
329: @*/
330: int DAGetAdicArray(DA da,PetscTruth ghosted,void **iptr,void **array_start,int *tdof)
331: {
332: int ierr,j,i,deriv_type_size,xs,ys,xm,ym,zs,zm,itdof;
333: char *iarray_start;
337: if (ghosted) {
338: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
339: if (da->adarrayghostedin[i]) {
340: *iptr = da->adarrayghostedin[i];
341: iarray_start = (char*)da->adstartghostedin[i];
342: itdof = da->ghostedtdof;
343: da->adarrayghostedin[i] = PETSC_NULL;
344: da->adstartghostedin[i] = PETSC_NULL;
345:
346: goto done;
347: }
348: }
349: xs = da->Xs;
350: ys = da->Ys;
351: zs = da->Zs;
352: xm = da->Xe-da->Xs;
353: ym = da->Ye-da->Ys;
354: zm = da->Ze-da->Zs;
355: } else {
356: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
357: if (da->adarrayin[i]) {
358: *iptr = da->adarrayin[i];
359: iarray_start = (char*)da->adstartin[i];
360: itdof = da->tdof;
361: da->adarrayin[i] = PETSC_NULL;
362: da->adstartin[i] = PETSC_NULL;
363:
364: goto done;
365: }
366: }
367: xs = da->xs;
368: ys = da->ys;
369: zs = da->zs;
370: xm = da->xe-da->xs;
371: ym = da->ye-da->ys;
372: zm = da->ze-da->zs;
373: }
374: deriv_type_size = PetscADGetDerivTypeSize();
376: switch (da->dim) {
377: case 1: {
378: void *ptr;
379: itdof = xm;
381: PetscMalloc(xm*deriv_type_size,&iarray_start);
383: ptr = (void*)(iarray_start - xs*deriv_type_size);
384: *iptr = (void*)ptr;
385: break;}
386: case 2: {
387: void **ptr;
388: itdof = xm*ym;
390: PetscMalloc((ym+1)*sizeof(void *)+xm*ym*deriv_type_size,&iarray_start);
392: ptr = (void**)(iarray_start + xm*ym*deriv_type_size - ys*sizeof(void*));
393: for(j=ys;j<ys+ym;j++) {
394: ptr[j] = iarray_start + deriv_type_size*(xm*(j-ys) - xs);
395: }
396: *iptr = (void*)ptr;
397: break;}
398: case 3: {
399: void ***ptr,**bptr;
400: itdof = xm*ym*zm;
402: PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*deriv_type_size,&iarray_start);
404: ptr = (void***)(iarray_start + xm*ym*zm*deriv_type_size - zs*sizeof(void*));
405: bptr = (void**)(iarray_start + xm*ym*zm*deriv_type_size + zm*sizeof(void**));
407: for(i=zs;i<zs+zm;i++) {
408: ptr[i] = bptr + ((i-zs)*ym - ys);
409: }
410: for (i=zs; i<zs+zm; i++) {
411: for (j=ys; j<ys+ym; j++) {
412: ptr[i][j] = iarray_start + deriv_type_size*(xm*ym*(i-zs) + xm*(j-ys) - xs);
413: }
414: }
416: *iptr = (void*)ptr;
417: break;}
418: default:
419: SETERRQ1(1,"Dimension %d not supported",da->dim);
420: }
422: done:
423: /* add arrays to the checked out list */
424: if (ghosted) {
425: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
426: if (!da->adarrayghostedout[i]) {
427: da->adarrayghostedout[i] = *iptr ;
428: da->adstartghostedout[i] = iarray_start;
429: da->ghostedtdof = itdof;
430: break;
431: }
432: }
433: } else {
434: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
435: if (!da->adarrayout[i]) {
436: da->adarrayout[i] = *iptr ;
437: da->adstartout[i] = iarray_start;
438: da->tdof = itdof;
439: break;
440: }
441: }
442: }
443: if (i == DA_MAX_AD_ARRAYS+1) SETERRQ(1,"Too many DA ADIC arrays obtained");
444: if (tdof) *tdof = itdof;
445: if (array_start) *array_start = iarray_start;
446: return(0);
447: }
451: /*@C
452: DARestoreAdicArray - Restores an array of derivative types for a DA
453:
454: Input Parameter:
455: + da - information about my local patch
456: - ghosted - do you want arrays for the ghosted or nonghosted patch
458: Output Parameters:
459: + ptr - array data structured to be passed to ad_FormFunctionLocal()
460: . array_start - actual start of 1d array of all values that adiC can access directly
461: - tdof - total number of degrees of freedom represented in array_start
463: Level: advanced
465: .seealso: DAGetAdicArray()
467: @*/
468: int DARestoreAdicArray(DA da,PetscTruth ghosted,void **iptr,void **array_start,int *tdof)
469: {
470: int i;
471: void *iarray_start = 0;
472:
475: if (ghosted) {
476: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
477: if (da->adarrayghostedout[i] == *iptr) {
478: iarray_start = da->adstartghostedout[i];
479: da->adarrayghostedout[i] = PETSC_NULL;
480: da->adstartghostedout[i] = PETSC_NULL;
481: break;
482: }
483: }
484: if (!iarray_start) SETERRQ(1,"Could not find array in checkout list");
485: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
486: if (!da->adarrayghostedin[i]){
487: da->adarrayghostedin[i] = *iptr;
488: da->adstartghostedin[i] = iarray_start;
489: break;
490: }
491: }
492: } else {
493: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
494: if (da->adarrayout[i] == *iptr) {
495: iarray_start = da->adstartout[i];
496: da->adarrayout[i] = PETSC_NULL;
497: da->adstartout[i] = PETSC_NULL;
498: break;
499: }
500: }
501: if (!iarray_start) SETERRQ(1,"Could not find array in checkout list");
502: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
503: if (!da->adarrayin[i]){
504: da->adarrayin[i] = *iptr;
505: da->adstartin[i] = iarray_start;
506: break;
507: }
508: }
509: }
510: return(0);
511: }
515: int ad_DAGetArray(DA da,PetscTruth ghosted,void **iptr)
516: {
519: DAGetAdicArray(da,ghosted,iptr,0,0);
520: return(0);
521: }
525: int ad_DARestoreArray(DA da,PetscTruth ghosted,void **iptr)
526: {
529: DARestoreAdicArray(da,ghosted,iptr,0,0);
530: return(0);
531: }
533: #endif
537: /*@C
538: DAGetArray - Gets a work array for a DA
539:
540: Input Parameter:
541: + da - information about my local patch
542: - ghosted - do you want arrays for the ghosted or nonghosted patch
544: Output Parameters:
545: . ptr - array data structured
547: Note: The vector values are NOT initialized and may have garbage in them, so you may need
548: to zero them.
550: Level: advanced
552: .seealso: DARestoreArray(), DAGetAdicArray()
554: @*/
555: int DAGetArray(DA da,PetscTruth ghosted,void **iptr)
556: {
557: int ierr,j,i,xs,ys,xm,ym,zs,zm;
558: char *iarray_start;
562: if (ghosted) {
563: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
564: if (da->arrayghostedin[i]) {
565: *iptr = da->arrayghostedin[i];
566: iarray_start = (char*)da->startghostedin[i];
567: da->arrayghostedin[i] = PETSC_NULL;
568: da->startghostedin[i] = PETSC_NULL;
569:
570: goto done;
571: }
572: }
573: xs = da->Xs;
574: ys = da->Ys;
575: zs = da->Zs;
576: xm = da->Xe-da->Xs;
577: ym = da->Ye-da->Ys;
578: zm = da->Ze-da->Zs;
579: } else {
580: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
581: if (da->arrayin[i]) {
582: *iptr = da->arrayin[i];
583: iarray_start = (char*)da->startin[i];
584: da->arrayin[i] = PETSC_NULL;
585: da->startin[i] = PETSC_NULL;
586:
587: goto done;
588: }
589: }
590: xs = da->xs;
591: ys = da->ys;
592: zs = da->zs;
593: xm = da->xe-da->xs;
594: ym = da->ye-da->ys;
595: zm = da->ze-da->zs;
596: }
598: switch (da->dim) {
599: case 1: {
600: void *ptr;
602: PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);
604: ptr = (void*)(iarray_start - xs*sizeof(PetscScalar));
605: *iptr = (void*)ptr;
606: break;}
607: case 2: {
608: void **ptr;
610: PetscMalloc((ym+1)*sizeof(void *)+xm*ym*sizeof(PetscScalar),&iarray_start);
612: ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*));
613: for(j=ys;j<ys+ym;j++) {
614: ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs);
615: }
616: *iptr = (void*)ptr;
617: break;}
618: case 3: {
619: void ***ptr,**bptr;
621: PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);
623: ptr = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*));
624: bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**));
625: for(i=zs;i<zs+zm;i++) {
626: ptr[i] = bptr + ((i-zs)*ym - ys);
627: }
628: for (i=zs; i<zs+zm; i++) {
629: for (j=ys; j<ys+ym; j++) {
630: ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
631: }
632: }
634: *iptr = (void*)ptr;
635: break;}
636: default:
637: SETERRQ1(1,"Dimension %d not supported",da->dim);
638: }
640: done:
641: /* add arrays to the checked out list */
642: if (ghosted) {
643: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
644: if (!da->arrayghostedout[i]) {
645: da->arrayghostedout[i] = *iptr ;
646: da->startghostedout[i] = iarray_start;
647: break;
648: }
649: }
650: } else {
651: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
652: if (!da->arrayout[i]) {
653: da->arrayout[i] = *iptr ;
654: da->startout[i] = iarray_start;
655: break;
656: }
657: }
658: }
659: return(0);
660: }
664: /*@C
665: DARestoreArray - Restores an array of derivative types for a DA
666:
667: Input Parameter:
668: + da - information about my local patch
669: . ghosted - do you want arrays for the ghosted or nonghosted patch
670: - ptr - array data structured to be passed to ad_FormFunctionLocal()
672: Level: advanced
674: .seealso: DAGetArray(), DAGetAdicArray()
676: @*/
677: int DARestoreArray(DA da,PetscTruth ghosted,void **iptr)
678: {
679: int i;
680: void *iarray_start = 0;
681:
684: if (ghosted) {
685: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
686: if (da->arrayghostedout[i] == *iptr) {
687: iarray_start = da->startghostedout[i];
688: da->arrayghostedout[i] = PETSC_NULL;
689: da->startghostedout[i] = PETSC_NULL;
690: break;
691: }
692: }
693: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
694: if (!da->arrayghostedin[i]){
695: da->arrayghostedin[i] = *iptr;
696: da->startghostedin[i] = iarray_start;
697: break;
698: }
699: }
700: } else {
701: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
702: if (da->arrayout[i] == *iptr) {
703: iarray_start = da->startout[i];
704: da->arrayout[i] = PETSC_NULL;
705: da->startout[i] = PETSC_NULL;
706: break;
707: }
708: }
709: for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
710: if (!da->arrayin[i]){
711: da->arrayin[i] = *iptr;
712: da->startin[i] = iarray_start;
713: break;
714: }
715: }
716: }
717: return(0);
718: }
722: /*@C
723: DAGetAdicMFArray - Gets an array of derivative types for a DA for matrix-free ADIC.
724:
725: Input Parameter:
726: + da - information about my local patch
727: - ghosted - do you want arrays for the ghosted or nonghosted patch?
729: Output Parameters:
730: + iptr - array data structured to be passed to ad_FormFunctionLocal()
731: . array_start - actual start of 1d array of all values that adiC can access directly (may be null)
732: - tdof - total number of degrees of freedom represented in array_start (may be null)
734: Notes:
735: The vector values are NOT initialized and may have garbage in them, so you may need
736: to zero them.
738: This routine returns the same type of object as the DAVecGetArray(), except its
739: elements are derivative types instead of PetscScalars.
741: Level: advanced
743: .seealso: DARestoreAdicMFArray(), DAGetArray(), DAGetAdicArray()
745: @*/
746: int DAGetAdicMFArray(DA da,PetscTruth ghosted,void **iptr,void **array_start,int *tdof)
747: {
748: int ierr,j,i,xs,ys,xm,ym,zs,zm,itdof;
749: char *iarray_start;
753: if (ghosted) {
754: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
755: if (da->admfarrayghostedin[i]) {
756: *iptr = da->admfarrayghostedin[i];
757: iarray_start = (char*)da->admfstartghostedin[i];
758: itdof = da->ghostedtdof;
759: da->admfarrayghostedin[i] = PETSC_NULL;
760: da->admfstartghostedin[i] = PETSC_NULL;
761:
762: goto done;
763: }
764: }
765: xs = da->Xs;
766: ys = da->Ys;
767: zs = da->Zs;
768: xm = da->Xe-da->Xs;
769: ym = da->Ye-da->Ys;
770: zm = da->Ze-da->Zs;
771: } else {
772: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
773: if (da->admfarrayin[i]) {
774: *iptr = da->admfarrayin[i];
775: iarray_start = (char*)da->admfstartin[i];
776: itdof = da->tdof;
777: da->admfarrayin[i] = PETSC_NULL;
778: da->admfstartin[i] = PETSC_NULL;
779:
780: goto done;
781: }
782: }
783: xs = da->xs;
784: ys = da->ys;
785: zs = da->zs;
786: xm = da->xe-da->xs;
787: ym = da->ye-da->ys;
788: zm = da->ze-da->zs;
789: }
791: switch (da->dim) {
792: case 1: {
793: void *ptr;
794: itdof = xm;
796: PetscMalloc(xm*2*sizeof(PetscScalar),&iarray_start);
798: ptr = (void*)(iarray_start - xs*2*sizeof(PetscScalar));
799: *iptr = (void*)ptr;
800: break;}
801: case 2: {
802: void **ptr;
803: itdof = xm*ym;
805: PetscMalloc((ym+1)*sizeof(void *)+xm*ym*2*sizeof(PetscScalar),&iarray_start);
807: ptr = (void**)(iarray_start + xm*ym*2*sizeof(PetscScalar) - ys*sizeof(void*));
808: for(j=ys;j<ys+ym;j++) {
809: ptr[j] = iarray_start + 2*sizeof(PetscScalar)*(xm*(j-ys) - xs);
810: }
811: *iptr = (void*)ptr;
812: break;}
813: case 3: {
814: void ***ptr,**bptr;
815: itdof = xm*ym*zm;
817: PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*2*sizeof(PetscScalar),&iarray_start);
819: ptr = (void***)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) - zs*sizeof(void*));
820: bptr = (void**)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) + zm*sizeof(void**));
821: for(i=zs;i<zs+zm;i++) {
822: ptr[i] = bptr + ((i-zs)*ym* - ys)*sizeof(void*);
823: }
824: for (i=zs; i<zs+zm; i++) {
825: for (j=ys; j<ys+ym; j++) {
826: ptr[i][j] = iarray_start + 2*sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
827: }
828: }
830: *iptr = (void*)ptr;
831: break;}
832: default:
833: SETERRQ1(1,"Dimension %d not supported",da->dim);
834: }
836: done:
837: /* add arrays to the checked out list */
838: if (ghosted) {
839: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
840: if (!da->admfarrayghostedout[i]) {
841: da->admfarrayghostedout[i] = *iptr ;
842: da->admfstartghostedout[i] = iarray_start;
843: da->ghostedtdof = itdof;
844: break;
845: }
846: }
847: } else {
848: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
849: if (!da->admfarrayout[i]) {
850: da->admfarrayout[i] = *iptr ;
851: da->admfstartout[i] = iarray_start;
852: da->tdof = itdof;
853: break;
854: }
855: }
856: }
857: if (i == DA_MAX_AD_ARRAYS+1) SETERRQ(1,"Too many DA ADIC arrays obtained");
858: if (tdof) *tdof = itdof;
859: if (array_start) *array_start = iarray_start;
860: return(0);
861: }
865: /*@C
866: DARestoreAdicMFArray - Restores an array of derivative types for a DA.
867:
868: Input Parameter:
869: + da - information about my local patch
870: - ghosted - do you want arrays for the ghosted or nonghosted patch?
872: Output Parameters:
873: + ptr - array data structure to be passed to ad_FormFunctionLocal()
874: . array_start - actual start of 1d array of all values that adiC can access directly
875: - tdof - total number of degrees of freedom represented in array_start
877: Level: advanced
879: .seealso: DAGetAdicArray()
881: @*/
882: int DARestoreAdicMFArray(DA da,PetscTruth ghosted,void **iptr,void **array_start,int *tdof)
883: {
884: int i;
885: void *iarray_start = 0;
886:
889: if (ghosted) {
890: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
891: if (da->admfarrayghostedout[i] == *iptr) {
892: iarray_start = da->admfstartghostedout[i];
893: da->admfarrayghostedout[i] = PETSC_NULL;
894: da->admfstartghostedout[i] = PETSC_NULL;
895: break;
896: }
897: }
898: if (!iarray_start) SETERRQ(1,"Could not find array in checkout list");
899: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
900: if (!da->admfarrayghostedin[i]){
901: da->admfarrayghostedin[i] = *iptr;
902: da->admfstartghostedin[i] = iarray_start;
903: break;
904: }
905: }
906: } else {
907: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
908: if (da->admfarrayout[i] == *iptr) {
909: iarray_start = da->admfstartout[i];
910: da->admfarrayout[i] = PETSC_NULL;
911: da->admfstartout[i] = PETSC_NULL;
912: break;
913: }
914: }
915: if (!iarray_start) SETERRQ(1,"Could not find array in checkout list");
916: for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
917: if (!da->admfarrayin[i]){
918: da->admfarrayin[i] = *iptr;
919: da->admfstartin[i] = iarray_start;
920: break;
921: }
922: }
923: }
924: return(0);
925: }
929: int admf_DAGetArray(DA da,PetscTruth ghosted,void **iptr)
930: {
933: DAGetAdicMFArray(da,ghosted,iptr,0,0);
934: return(0);
935: }
939: int admf_DARestoreArray(DA da,PetscTruth ghosted,void **iptr)
940: {
943: DARestoreAdicMFArray(da,ghosted,iptr,0,0);
944: return(0);
945: }