Actual source code: adelement.c
1: #ifdef ad_GRAD_MAX
2: #undef ad_GRAD_MAX
3: #endif
4: #define ad_GRAD_MAX 4
6: #include taodaapplication.h
7: #include taoapp.h
9: //#include "ad_deriv.h"
11: /* Function Gradient Evaluation over single element */
12: typedef struct {
13: int (*computeadicfunctiongradient)(int[2],DERIV_TYPE[4],DERIV_TYPE*,void*);
14: void *elementfgctx;
15: int elementfgflops;
16: DERIV_TYPE adX[4];
17: } TaoDA2D1DOFADICFGCtx;
22: static int TaoDA2dLoopADFunctionGradient(TAO_APPLICATION tao, DA da, Vec X, double *f, Vec G, void * ctx) {
24: TaoDA2D1DOFADICFGCtx *myapp = (TaoDA2D1DOFADICFGCtx*) ctx;
25: MPI_Comm comm;
26: Vec localX, localG;
27: int info, i, j, coor[2];
28: int xs, xm, gxs, gxm, xe, ys, ym, gys, gym, ye;
29: PetscScalar **x, **g;
30: PetscScalar floc = 0.0;
31: PetscScalar zero = 0.0;
32: DERIV_TYPE adF,*adX=myapp->adX;
34: info = DAGetLocalVector(da, &localX); CHKERRQ(info);
35: info = DAGetLocalVector(da, &localG); CHKERRQ(info);
36: info = VecSet(G, zero); CHKERRQ(info);
37: info = VecSet(localG, zero); CHKERRQ(info);
39: info = DAGlobalToLocalBegin(da, X, INSERT_VALUES, localX); CHKERRQ(info);
40: info = DAGlobalToLocalEnd(da, X, INSERT_VALUES, localX); CHKERRQ(info);
42: info = DAVecGetArray(da, localX, (void**)&x); CHKERRQ(info);
43: info = DAVecGetArray(da, localG, (void**)&g); CHKERRQ(info);
45: info = DAGetCorners(da, &xs, &ys, PETSC_NULL, &xm, &ym, PETSC_NULL); CHKERRQ(info);
46: info = DAGetGhostCorners(da, &gxs, &gys, PETSC_NULL, &gxm, &gym, PETSC_NULL); CHKERRQ(info);
48: xe = gxs + gxm - 1;
49: ye = gys + gym - 1;
51: for (j = ys; j < ye; j++) {
52: for (i = xs; i < xe; i++) {
54: DERIV_val(adX[0]) = x[j][i];
55: DERIV_val(adX[1]) = x[j][i+1];
56: DERIV_val(adX[2]) = x[j+1][i];
57: DERIV_val(adX[3]) = x[j+1][i+1];
58: coor[0] = i; coor[1] = j;
60: info = myapp->computeadicfunctiongradient(coor,adX,&adF,myapp->elementfgctx);
61: CHKERRQ(info);
63: floc += DERIV_val(adF);
65: g[j][i] += DERIV_grad(adF)[0];
66: g[j][i+1] += DERIV_grad(adF)[1];
67: g[j+1][i] += DERIV_grad(adF)[2];
68: g[j+1][i+1] += DERIV_grad(adF)[3];
69: }
70: }
72: PetscLogFlops((ye-ys)*(xe-xs)*(myapp->elementfgflops + 5));
74: PetscObjectGetComm((PetscObject)da,&comm); CHKERRQ(info);
75: info = MPI_Allreduce(&floc, f, 1, MPI_DOUBLE, MPI_SUM, comm); CHKERRQ(info);
77: info = DAVecRestoreArray(da, localX, (void**)&x); CHKERRQ(info);
78: info = DAVecRestoreArray(da, localG, (void**)&g); CHKERRQ(info);
80: info = DALocalToGlobalBegin(da, localG, G); CHKERRQ(info);
81: info = DALocalToGlobalEnd(da, localG, G); CHKERRQ(info);
83: info = DARestoreLocalVector(da, &localX); CHKERRQ(info);
84: info = DARestoreLocalVector(da, &localG); CHKERRQ(info);
86: return(0);
87: } /* TaoDA2dLoopADFunctionGradient */
91: static int TaoShutDownADICQQQQ(void* ctx){
93: ad_AD_Final();
94: return(0);
95: }
99: /*@
100: DAAppSetADElementFunctionGradient - Set routine that evaluates the
101: local part of a function on a 2-dimensional DA with 1 degree of freedom.
103: Collective on TAO_APPLICATION
105: Input Parameters:
106: + daapp - the TAO_DA_APPLICATION solver context
107: . funcgrad - local function gradient routine
108: . flops - the number of flops done performed in the funcgrad routine
109: - fgctx - [optional] user-defined context for private data for the evaluation.
111: Calling sequence of funcgrad:
112: $ int funcgrad(int coordinates[2], PetscScalar x[4], double *f, PetscScalar g[4], void* ctx)
114: + coord - the global coordinates [i j] in each direction of the DA
115: . x - the variables on the DA ( da[j][i], da[j][j+1], da[j+1][i], da[j+1][i+1] ) (bottom left, bottom right, top left, top right)
116: . g - the ADIC differentiated objective function with respect to each variable
117: - ctx - user defined context
120: Level: intermediate
122: .keywords: DA, gradient, ADIC
124: .seealso: DAAppSetObjectiveAndGradientRoutine(), DAAppSetElementObjectiveAndGradientRoutine()
125: @*/
126: int DAAppSetADElementFunctionGradient(TAO_APPLICATION daapplication,
127: int (*funcgrad)(int[2],DERIV_TYPE[4],DERIV_TYPE*,void*),
128: int flops, void *ctx){
129: int i,n,info;
130: int dim,dof,s;
131: DAStencilType st;
132: TaoDA2D1DOFADICFGCtx *fgctx;
133: DA da;
136: info=DAAppGetNumberOfDAGrids(daapplication,&n); CHKERRQ(info);
137: for (i=0;i<n;i++){
138: info = DAAppGetDA(daapplication, i, &da); CHKERRQ(info);
139: info = DAGetInfo(da,&dim,0,0,0,0,0,0,&dof,&s,0,&st); CHKERRQ(info);
140: if (dim!=2){
141: SETERRQ(1,"TAO DA ERROR: DA must have dimension of 2");}
142: if (dof!=1){
143: SETERRQ(1,"TAO DA ERROR: DA must have exactly 1 degree of freedom per node");}
144: if (s!=1){
145: SETERRQ(1,"TAO DA ERROR: DA stencil width must equal 1"); }
146: if (st!=DA_STENCIL_BOX){
147: SETERRQ(1,"TAO DA ERROR: DA stencil must be DA_STENCIL_BOX");}
148: }
149: PetscNew(TaoDA2D1DOFADICFGCtx,&fgctx);
150: // ad_AD_Init(4);
151: ad_AD_Init(ad_GRAD_MAX);
152: ad_AD_SetIndepArray(fgctx->adX,4);
153: ad_AD_SetIndepDone();
154: fgctx->computeadicfunctiongradient=funcgrad;
155: fgctx->elementfgctx=ctx;
156: fgctx->elementfgflops=flops;
157: info = DAAppSetObjectiveAndGradientRoutine(daapplication, TaoDA2dLoopADFunctionGradient, (void*)fgctx);
158: CHKERRQ(info);
159: info = TaoAppSetDestroyRoutine(daapplication,TaoApplicationFreeMemory, (void*)fgctx); CHKERRQ(info);
160: info = TaoAppSetDestroyRoutine(daapplication, TaoShutDownADICQQQQ, 0); CHKERRQ(info);
161: info = PetscLogInfo((daapplication,"Set objective function pointer for TAO_DA_APPLICATION object.\n")); CHKERRQ(info);
162: return(0);
163: }