Actual source code: nladic.c

  1: #define PETSCMAT_DLL
  2: /*
  3:     ADIC based nonlinear operator object that can be used with FAS

  5:     This does not really belong in the matrix directories but since it 
  6:     was cloned off of Mat_DAAD I'm leaving it here until I have a better place

  8: */
 9:  #include petscsys.h
 10:  #include petscda.h

 13: #include "adic/ad_utils.h"

 16: #include "../src/dm/da/daimpl.h"
 17:  #include ../src/mat/blockinvert.h

 19: struct NLF_DAAD {
 20:   DA         da;
 21:   void       *ctx;
 22:   Vec        residual;
 23:   int        newton_its;
 24: };

 26: /*
 27:       Solves the one dimensional equation using Newton's method 
 28: */
 31: PetscErrorCode NLFNewton_DAAD(NLF A,DALocalInfo *info,MatStencil *stencil,void *ad_vu,PetscScalar *ad_vustart,int nI,int gI,PetscScalar residual)
 32: {
 34:   PetscInt       cnt = A->newton_its;
 35:   PetscScalar    ad_f[2],J,f;

 38:   ad_vustart[1+2*gI] = 1.0;
 39:   do {
 40:     /* compute the function and Jacobian */
 41:     (*A->da->adicmf_lfi)(info,stencil,ad_vu,ad_f,A->ctx);
 42:     J    = -ad_f[1];
 43:     f    = -ad_f[0] + residual;
 44:     if (f != f) SETERRQ(1,"nan");
 45:     ad_vustart[2*gI] =  ad_vustart[2*gI] - f/J;
 46:   } while (--cnt > 0 && PetscAbsScalar(f) > 1.e-14);

 48:   ad_vustart[1+2*gI] = 0.0;
 49:   return(0);
 50: }

 52: /*
 53:       Solves the four dimensional equation using Newton's method 
 54: */
 57: PetscErrorCode NLFNewton_DAAD4(NLF A,DALocalInfo *info,MatStencil *stencil,void *ad_vu,PetscScalar *ad_vustart,int nI,int gI,PetscScalar *residual)
 58: {
 60:   PetscInt       cnt = A->newton_its;
 61:   PetscScalar    ad_f[20], J[16],f[4], res, dd[5];


 65:   /* This sets the identity as the seed matrix for ADIC */
 66:     CHKMEMQ;
 67:   ad_vustart[1+5*gI   ] = 1.0;
 68:     CHKMEMQ;
 69:   ad_vustart[2+5*gI+5 ] = 1.0;
 70:     CHKMEMQ;
 71:   ad_vustart[3+5*gI+10] = 1.0;
 72:     CHKMEMQ;
 73:   ad_vustart[4+5*gI+15] = 1.0;
 74:     CHKMEMQ;

 76:   do {
 77:     /* compute the function and Jacobian */
 78:     CHKMEMQ;
 79:     (*A->da->adicmf_lfib)(info,stencil,ad_vu,ad_f,A->ctx);
 80:        CHKMEMQ;
 81:     /* copy ADIC formated Jacobian into regular C array */
 82:     J[0] = ad_f[1] ; J[1] = ad_f[2] ; J[2] = ad_f[3] ; J[3] = ad_f[4] ;
 83:     J[4] = ad_f[6] ; J[5] = ad_f[7] ; J[6] = ad_f[8] ; J[7] = ad_f[9] ;
 84:     J[8] = ad_f[11]; J[9] = ad_f[12]; J[10]= ad_f[13]; J[11]= ad_f[14];
 85:     J[12]= ad_f[16]; J[13]= ad_f[17]; J[14]= ad_f[18]; J[15]= ad_f[19];
 86:     CHKMEMQ;
 87:     f[0]    = -ad_f[0]   + residual[0];
 88:     f[1]    = -ad_f[5]   + residual[1];
 89:     f[2]    = -ad_f[10]  + residual[2];
 90:     f[3]    = -ad_f[15]  + residual[3];

 92:     /* solve Jacobian * dd = ff */

 94:     /* could use PETSc kernel code to solve system with pivoting */

 96:     /* could put code in here to compute the solution directly using ADIC data structures instead of copying first */
 97:     dd[0]=J[0]*(J[5]*(J[10]*J[15]-J[11]*J[14])-J[6]*(J[9]*J[15]-J[11]*J[13])+J[7]*(J[9]*J[14]-J[10]*J[13]))-
 98:           J[1]*(J[4]*(J[10]*J[15]-J[11]*J[14])-J[6]*(J[8]*J[15]-J[11]*J[12])+J[7]*(J[8]*J[14]-J[10]*J[12]))+
 99:           J[2]*(J[4]*(J[ 9]*J[15]-J[11]*J[13])-J[5]*(J[8]*J[15]-J[11]*J[12])+J[7]*(J[8]*J[13]-J[ 9]*J[12]))-
100:           J[3]*(J[4]*(J[ 9]*J[14]-J[10]*J[13])-J[5]*(J[8]*J[14]-J[10]*J[12])+J[6]*(J[8]*J[13]-J[ 9]*J[12]));

102:     dd[1]=(f[0]*(J[5]*(J[10]*J[15]-J[11]*J[14])-J[6]*(J[9]*J[15]-J[11]*J[13])+J[7]*(J[9]*J[14]-J[10]*J[13]))-
103:           J[1]*(f[1]*(J[10]*J[15]-J[11]*J[14])-J[6]*(f[2]*J[15]-J[11]*f[ 3])+J[7]*(f[2]*J[14]-J[10]*f[ 3]))+
104:           J[2]*(f[1]*(J[ 9]*J[15]-J[11]*J[13])-J[5]*(f[2]*J[15]-J[11]*f[ 3])+J[7]*(f[2]*J[13]-J[ 9]*f[ 3]))-
105:            J[3]*(f[1]*(J[ 9]*J[14]-J[10]*J[13])-J[5]*(f[2]*J[14]-J[10]*f[ 3])+J[6]*(f[2]*J[13]-J[ 9]*f[ 3])))/dd[0];

107:     dd[2]=(J[0]*(f[1]*(J[10]*J[15]-J[11]*J[14])-J[6]*(f[2]*J[15]-J[11]*f[ 3])+J[7]*(f[2]*J[14]-J[10]*f[ 3]))-
108:           f[0]*(J[4]*(J[10]*J[15]-J[11]*J[14])-J[6]*(J[8]*J[15]-J[11]*J[12])+J[7]*(J[8]*J[14]-J[10]*J[12]))+
109:           J[2]*(J[4]*(f[ 2]*J[15]-J[11]*f[ 3])-f[1]*(J[8]*J[15]-J[11]*J[12])+J[7]*(J[8]*f[ 3]-f[ 2]*J[12]))-
110:           J[3]*(J[4]*(f[ 2]*J[14]-J[10]*f[ 3])-f[2]*(J[8]*J[14]-J[10]*J[12])+J[6]*(J[8]*f[ 3]-f[ 2]*J[12])))/dd[0];

112:     dd[3]=(J[0]*(J[5]*(f[ 2]*J[15]-J[11]*f[ 3])-f[1]*(J[9]*J[15]-J[11]*J[13])+J[7]*(J[9]*f[ 3]-f[ 2]*J[13]))-
113:           J[1]*(J[4]*(f[ 2]*J[15]-J[11]*f[ 3])-f[1]*(J[8]*J[15]-J[11]*J[12])+J[7]*(J[8]*f[ 3]-f[ 2]*J[12]))+
114:           f[0]*(J[4]*(J[ 9]*J[15]-J[11]*J[13])-J[5]*(J[8]*J[15]-J[11]*J[12])+J[7]*(J[8]*J[13]-J[ 9]*J[12]))-
115:            J[3]*(J[4]*(J[ 9]*f[ 3]-f[ 2]*J[13])-J[5]*(J[8]*f[ 3]-f[ 2]*J[12])+f[1]*(J[8]*J[13]-J[ 9]*J[12])))/dd[0];

117:     dd[4]=(J[0]*(J[5]*(J[10]*f[ 3]-f[ 2]*J[14])-J[6]*(J[9]*f[ 3]-f[ 2]*J[13])+f[1]*(J[9]*J[14]-J[10]*J[13]))-
118:           J[1]*(J[4]*(J[10]*f[ 3]-f[ 2]*J[14])-J[6]*(J[8]*f[ 3]-f[ 2]*J[12])+f[1]*(J[8]*J[14]-J[10]*J[12]))+
119:           J[2]*(J[4]*(J[ 9]*f[ 3]-f[ 2]*J[13])-J[5]*(J[8]*f[ 3]-f[ 2]*J[12])+f[1]*(J[8]*J[13]-J[ 9]*J[12]))-
120:           f[0]*(J[4]*(J[ 9]*J[14]-J[10]*J[13])-J[5]*(J[8]*J[14]-J[10]*J[12])+J[6]*(J[8]*J[13]-J[ 9]*J[12])))/dd[0];
121:     CHKMEMQ;
122:     /* copy solution back into ADIC data structure */
123:     ad_vustart[5*(gI+0)] += dd[1];
124:     ad_vustart[5*(gI+1)] += dd[2];
125:     ad_vustart[5*(gI+2)] += dd[3];
126:     ad_vustart[5*(gI+3)] += dd[4];
127:     CHKMEMQ;
128:     res =  f[0]*f[0];
129:     res += f[1]*f[1];
130:     res += f[2]*f[2];
131:     res += f[3]*f[3];
132:     res =  sqrt(res);
133: 
134:   } while (--cnt > 0 && res > 1.e-14);

136:   /* zero out this part of the seed matrix that was set initially */
137:   ad_vustart[1+5*gI   ] = 0.0;
138:   ad_vustart[2+5*gI+5 ] = 0.0;
139:   ad_vustart[3+5*gI+10] = 0.0;
140:   ad_vustart[4+5*gI+15] = 0.0;
141:     CHKMEMQ;
142:   return(0);
143: }

147: PetscErrorCode NLFNewton_DAAD9(NLF A,DALocalInfo *info,MatStencil *stencil,void *ad_vu,PetscScalar *ad_vustart,int nI,int gI,PetscScalar *residual)
148: {
150:   PetscInt       cnt = A->newton_its;
151:   PetscScalar    ad_f[100], J[81],f[9], res;
152:   PetscInt       i,j,ngI[9];
154: 
155:   // the order of the nodes
156:    /*
157:          (6)      (7)         (8)
158:       i-1,j+1 --- i,j+1 --- i+1,j+1
159:         |         |           |
160:         |         |           |
161:       i-1,j   --- i,j  --- i+1,j
162:         |(3)      |(4)        |(5)
163:         |         |           |
164:       i-1,j-1 --- i,j-1--- i+1,j-1
165:        (0)       (1)         (2)
166:   */
167: 
168:   // the order of the derivative for the center nodes
169:    /*
170:          (7)      (8)         (9)
171:       i-1,j+1 --- i,j+1 --- i+1,j+1
172:         |         |           |
173:         |         |           |
174:       i-1,j   --- i,j  --- i+1,j
175:         |(4)      |(5)        |(6)
176:         |         |           |
177:       i-1,j-1 --- i,j-1--- i+1,j-1
178:        (1)       (2)         (3)
179:   */
180:   if( (*stencil).i==0 || (*stencil).i==1||(*stencil).i==(*info).gxs+(*info).gxm-1 || (*stencil).i==(*info).gxs+(*info).gxm-2  || (*stencil).j==0 ||  (*stencil).j==1 ||(*stencil).j==(*info).gys+(*info).gym-1 || (*stencil).j==(*info).gys+(*info).gym -2) {

182:   ad_vustart[1+10*gI] = 1.0;
183: 
184:   do {
185:     /* compute the function and Jacobian */
186:     (*A->da->adicmf_lfi)(info,stencil,ad_vu,ad_f,A->ctx);
187:     J[0]    = -ad_f[1];
188:     f[0]    = -ad_f[0] + residual[gI];
189:     ad_vustart[10*gI] =  ad_vustart[10*gI] - f[0]/J[0];
190:   } while (--cnt > 0 && PetscAbsScalar(f[0]) > 1.e-14);

192:   ad_vustart[1+10*gI] = 0.0;
193:   return(0);


196:   }
197: 
198:   ngI[0]  =  ((*stencil).i -1 - (*info).gxs)*(*info).dof + ((*stencil).j -1 - (*info).gys)*(*info).dof*(*info).gxm + ((*stencil).k - (*info).gzs)*(*info).dof*(*info).gxm*(*info).gym;
199:   ngI[1]  =  ngI[0] + 1;
200:   ngI[2]  =  ngI[1] + 1;
201:   ngI[3]  =  gI     - 1;
202:   ngI[4]  =  gI        ;
203:   ngI[5]  =  gI     + 1;
204:   ngI[6]  =  ((*stencil).i -1 - (*info).gxs)*(*info).dof + ((*stencil).j +1 - (*info).gys)*(*info).dof*(*info).gxm + ((*stencil).k - (*info).gzs)*(*info).dof*(*info).gxm*(*info).gym;
205:   ngI[7]  =  ngI[6] + 1;
206:   ngI[8]  =  ngI[7] + 1;
207: 

209:   for(j=0 ; j<9; j++){
210:     ad_vustart[ngI[j]*10+j+1] = 1.0;
211:   }
212: 
213:   do{
214:     /* compute the function and the Jacobian */
215: 
216:     (*A->da->adicmf_lfi)(info,stencil,ad_vu,ad_f,A->ctx);
217: 
218:     for(i=0; i<9; i++){
219:       for(j=0; j<9; j++){
220:         J[i*9+j] = -ad_f[i*10+j+1];
221:       }
222:       f[i]= -ad_f[i*10] + residual[ngI[i]];
223:     }

225:     Kernel_A_gets_inverse_A_9(J,0.0);
226: 
227:     res =0 ;
228:     for(i=0; i<9; i++){
229:       for(j=0;j<9;j++){
230:         ad_vustart[10*ngI[i]]= ad_vustart[10*ngI[i]] - J[i*9 +j]*f[j];
231:       }
232:       res = res + f[i]*f[i];
233:     }
234:     res = sqrt(res);
235:   } while(--cnt>0 && res>1.e-14);

237:   for(j=0; j<9; j++){
238:     ad_vustart[10*ngI[j]+j+1]=0.0;
239:   }

241:   return(0);
242: }

244: /*
245:         Nonlinear relax on all the equations with an initial guess in x
246: */
250: PetscErrorCode  NLFRelax_DAAD(NLF A,MatSORType flag,int its,Vec xx)
251: {
253:   PetscInt       j,gtdof,nI,gI;
254:   PetscScalar    *avu,*av,*ad_vustart,*residual;
255:   Vec            localxx;
256:   DALocalInfo    info;
257:   MatStencil     stencil;
258:   void*          *ad_vu;

261:   if (its <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D positive",its);

263:   DAGetLocalVector(A->da,&localxx);
264:   /* get space for derivative object.  */
265:   DAGetAdicMFArray(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,&gtdof);
266:   VecGetArray(A->residual,&residual);


269:   /* tell ADIC we will be computing one dimensional Jacobians */
270:   PetscADResetIndep();
271:   PetscADIncrementTotalGradSize(1);
272:   PetscADSetIndepDone();

274:   DAGetLocalInfo(A->da,&info);
275:   while (its--) {

277:     /* get initial solution properly ghosted */
278:     DAGlobalToLocalBegin(A->da,xx,INSERT_VALUES,localxx);
279:     DAGlobalToLocalEnd(A->da,xx,INSERT_VALUES,localxx);

281:     /* copy input vector into derivative object */
282:     VecGetArray(localxx,&avu);
283:     for (j=0; j<gtdof; j++) {
284:       ad_vustart[2*j]   = avu[j];
285:       ad_vustart[2*j+1] = 0.0;
286:     }
287: 
288:     VecRestoreArray(localxx,&avu);

290:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
291:       nI = 0;
292:       for (stencil.k = info.zs; stencil.k<info.zs+info.zm; stencil.k++) {
293:         for (stencil.j = info.ys; stencil.j<info.ys+info.ym; stencil.j++) {
294:           for (stencil.i = info.xs; stencil.i<info.xs+info.xm; stencil.i++) {
295:             for (stencil.c = 0; stencil.c<info.dof; stencil.c++) {
296:               gI   = stencil.c + (stencil.i - info.gxs)*info.dof + (stencil.j - info.gys)*info.dof*info.gxm + (stencil.k - info.gzs)*info.dof*info.gxm*info.gym;
297:               NLFNewton_DAAD(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual[nI]);
298:               nI++;
299:             }
300:           }
301:         }
302:       }
303:     }
304:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
305:       nI = info.dof*info.xm*info.ym*info.zm - 1;
306:       for (stencil.k = info.zs+info.zm-1; stencil.k>=info.zs; stencil.k--) {
307:         for (stencil.j = info.ys+info.ym-1; stencil.j>=info.ys; stencil.j--) {
308:           for (stencil.i = info.xs+info.xm-1; stencil.i>=info.xs; stencil.i--) {
309:             for (stencil.c = info.dof-1; stencil.c>=0; stencil.c--) {
310:               gI   = stencil.c + (stencil.i - info.gxs)*info.dof + (stencil.j - info.gys)*info.dof*info.gxm + (stencil.k - info.gzs)*info.dof*info.gxm*info.gym;
311:               NLFNewton_DAAD(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual[nI]);
312:               nI--;
313:             }
314:           }
315:         }
316:       }
317:     }

319:     /* copy solution back into ghosted vector from derivative object */
320:     VecGetArray(localxx,&av);
321:     for (j=0; j<gtdof; j++) {
322:       av[j] = ad_vustart[2*j];
323:     }
324:     VecRestoreArray(localxx,&av);
325:     /* stick relaxed solution back into global solution */
326:     DALocalToGlobal(A->da,localxx,INSERT_VALUES,xx);
327:   }


330:   VecRestoreArray(A->residual,&residual);
331:   DARestoreLocalVector(A->da,&localxx);
332:   DARestoreAdicMFArray(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,&gtdof);
333:   return(0);
334: }

340: PetscErrorCode  NLFRelax_DAAD4(NLF A,MatSORType flag,int its,Vec xx)
341: {
343:   PetscInt       j,gtdof,nI,gI;
344:   PetscScalar    *avu,*av,*ad_vustart,*residual;
345:   Vec            localxx;
346:   DALocalInfo    info;
347:   MatStencil     stencil;
348:   void*          *ad_vu;

351:   if (its <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D positive",its);
352: 
353:   DAGetLocalVector(A->da,&localxx);
354:   /* get space for derivative object.  */
355:   DAGetAdicMFArray4(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,&gtdof);
356:   VecGetArray(A->residual,&residual);


359:   /* tell ADIC we will be computing four dimensional Jacobians */
360:   PetscADResetIndep();
361:   PetscADIncrementTotalGradSize(4);
362:   PetscADSetIndepDone();

364:   DAGetLocalInfo(A->da,&info);
365:   while (its--) {

367:     /* get initial solution properly ghosted */
368:     DAGlobalToLocalBegin(A->da,xx,INSERT_VALUES,localxx);
369:     DAGlobalToLocalEnd(A->da,xx,INSERT_VALUES,localxx);

371:     /* copy input vector into derivative object */
372:     VecGetArray(localxx,&avu);
373:     for (j=0; j<gtdof; j++) {
374:       ad_vustart[5*j  ]   = avu[j];
375:       ad_vustart[5*j+1]   = 0.0;
376:       ad_vustart[5*j+2]   = 0.0;
377:       ad_vustart[5*j+3]   = 0.0;
378:       ad_vustart[5*j+4]   = 0.0;
379:     }
380: 
381:     VecRestoreArray(localxx,&avu);

383:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
384:       nI = 0;
385:       for (stencil.k = info.zs; stencil.k<info.zs+info.zm; stencil.k++) {
386:         for (stencil.j = info.ys; stencil.j<info.ys+info.ym; stencil.j++) {
387:           for (stencil.i = info.xs; stencil.i<info.xs+info.xm; stencil.i++) {
388:             CHKMEMQ;
389:             gI   = (stencil.i - info.gxs)*info.dof + (stencil.j - info.gys)*info.dof*info.gxm + (stencil.k - info.gzs)*info.dof*info.gxm*info.gym;

391:             NLFNewton_DAAD4(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual+nI);
392:             nI=nI+4;
393:             CHKMEMQ;
394:           }
395:         }
396:       }
397:     }
398:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
399:       nI = info.dof*info.xm*info.ym*info.zm - 4;
400: 
401:       for (stencil.k = info.zs+info.zm-1; stencil.k>=info.zs; stencil.k--) {
402:         for (stencil.j = info.ys+info.ym-1; stencil.j>=info.ys; stencil.j--) {
403:           for (stencil.i = info.xs+info.xm-1; stencil.i>=info.xs; stencil.i--) {
404:             gI   = (stencil.i - info.gxs)*info.dof + (stencil.j - info.gys)*info.dof*info.gxm + (stencil.k - info.gzs)*info.dof*info.gxm*info.gym;
405:             NLFNewton_DAAD4(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual+nI);
406:             nI=nI-4;
407:           }
408:         }
409:       }
410:     }
411: 
412:     /* copy solution back into ghosted vector from derivative object */
413:     VecGetArray(localxx,&av);
414:     for (j=0; j<gtdof; j++) {
415:       av[j] = ad_vustart[5*j];
416:     }
417:     VecRestoreArray(localxx,&av);
418:     /* stick relaxed solution back into global solution */
419:     DALocalToGlobal(A->da,localxx,INSERT_VALUES,xx);
420:   }


423:   VecRestoreArray(A->residual,&residual);
424:   DARestoreLocalVector(A->da,&localxx);
425:   DARestoreAdicMFArray(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,&gtdof);
426:   return(0);
427: }

432: PetscErrorCode  NLFRelax_DAAD9(NLF A,MatSORType flag,int its,Vec xx)
433: {
435:   PetscInt       j,gtdof,nI,gI;
436:   PetscScalar    *avu,*av,*ad_vustart,*residual;
437:   Vec            localxx;
438:   DALocalInfo    info;
439:   MatStencil     stencil;
440:   void*          *ad_vu;

443:   if (its <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D positive",its);
444: 
445:   DAGetLocalVector(A->da,&localxx);
446:   /* get space for derivative object.  */
447:   DAGetAdicMFArray9(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,&gtdof);
448:   VecGetArray(A->residual,&residual);


451:   /* tell ADIC we will be computing nine dimensional Jacobians */
452:   PetscADResetIndep();
453:   PetscADIncrementTotalGradSize(9);
454:   PetscADSetIndepDone();

456:   DAGetLocalInfo(A->da,&info);
457:   while (its--) {

459:     /* get initial solution properly ghosted */
460:     DAGlobalToLocalBegin(A->da,xx,INSERT_VALUES,localxx);
461:     DAGlobalToLocalEnd(A->da,xx,INSERT_VALUES,localxx);

463:     /* copy input vector into derivative object */
464:     VecGetArray(localxx,&avu);
465:      for (j=0; j<gtdof; j++) {
466:       ad_vustart[10*j  ]   = avu[j];
467:       ad_vustart[10*j+1]   = 0.0;
468:       ad_vustart[10*j+2]   = 0.0;
469:       ad_vustart[10*j+3]   = 0.0;
470:       ad_vustart[10*j+4]   = 0.0;
471:       ad_vustart[10*j+5]   = 0.0;
472:       ad_vustart[10*j+6]   = 0.0;
473:       ad_vustart[10*j+7]   = 0.0;
474:       ad_vustart[10*j+8]   = 0.0;
475:       ad_vustart[10*j+9]   = 0.0;
476:     }
477: 
478:     VecRestoreArray(localxx,&avu);

480:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
481:       nI = 0;
482:       for (stencil.k = info.zs; stencil.k<info.zs+info.zm; stencil.k++) {
483:         for (stencil.j = info.ys; stencil.j<info.ys+info.ym; stencil.j++) {
484:           for (stencil.i = info.xs; stencil.i<info.xs+info.xm; stencil.i++) {
485:             CHKMEMQ;
486:               gI   = (stencil.i - info.gxs)*info.dof + (stencil.j - info.gys)*info.dof*info.gxm + (stencil.k - info.gzs)*info.dof*info.gxm*info.gym;
487:  NLFNewton_DAAD9(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual);
488:  nI=nI+1;
489:             CHKMEMQ;
490:           }
491:         }
492:       }
493:     }
494:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
495:       nI = info.dof*info.xm*info.ym*info.zm - 4;
496: 
497:       for (stencil.k = info.zs+info.zm-1; stencil.k>=info.zs; stencil.k--) {
498:         for (stencil.j = info.ys+info.ym-1; stencil.j>=info.ys; stencil.j--) {
499:           for (stencil.i = info.xs+info.xm-1; stencil.i>=info.xs; stencil.i--) {
500:               gI   = (stencil.i - info.gxs)*info.dof + (stencil.j - info.gys)*info.dof*info.gxm + (stencil.k - info.gzs)*info.dof*info.gxm*info.gym;
501:               NLFNewton_DAAD9(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual);
502:     nI=nI-1;
503:           }
504:         }
505:       }
506:     }
507: 
508:     /* copy solution back into ghosted vector from derivative object */
509:     VecGetArray(localxx,&av);
510:     for (j=0; j<gtdof; j++) {
511:       av[j] = ad_vustart[10*j];
512:     }
513:     VecRestoreArray(localxx,&av);
514:     /* stick relaxed solution back into global solution */
515:     DALocalToGlobal(A->da,localxx,INSERT_VALUES,xx);
516:   }


519:   VecRestoreArray(A->residual,&residual);
520:   DARestoreLocalVector(A->da,&localxx);
521:   DARestoreAdicMFArray(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,&gtdof);
522:   return(0);
523: }

525: /*
526:         Point-block nonlinear relax on all the equations with an initial guess in xx using 
527: */
531: PetscErrorCode  NLFRelax_DAADb(NLF A,MatSORType flag,int its,Vec xx)
532: {
534:   PetscInt       i,j,gtdof,nI,gI, bs = A->da->w, bs1 = bs + 1;
535:   PetscScalar    *avu,*av,*ad_vustart,*residual;
536:   Vec            localxx;
537:   DALocalInfo    info;
538:   MatStencil     stencil;
539:   void*          *ad_vu;
540:   PetscErrorCode (*NLFNewton_DAADb)(NLF,DALocalInfo*,MatStencil*,void*,PetscScalar*,int,int,PetscScalar*);

543:   if (its <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D positive",its);
544:   if (bs == 4) {
545:     NLFNewton_DAADb       = NLFNewton_DAAD4;
546:   } else {
547:     SETERRQ1(PETSC_ERR_SUP,"Point block nonlinear relaxation currently not for this block size",bs);
548:   }

550:   DAGetLocalVector(A->da,&localxx);
551:   /* get space for derivative object.  */
552:   DAGetAdicMFArrayb(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,&gtdof);
553:   VecGetArray(A->residual,&residual);


556:   /* tell ADIC we will be computing bs dimensional Jacobians */
557:   PetscADResetIndep();
558:   PetscADIncrementTotalGradSize(bs);
559:   PetscADSetIndepDone();

561:   DAGetLocalInfo(A->da,&info);
562:   while (its--) {

564:     /* get initial solution properly ghosted */
565:     DAGlobalToLocalBegin(A->da,xx,INSERT_VALUES,localxx);
566:     DAGlobalToLocalEnd(A->da,xx,INSERT_VALUES,localxx);

568:     /* copy input vector into derivative object */
569:     VecGetArray(localxx,&avu);
570:     for (j=0; j<gtdof; j++) {
571:       ad_vustart[bs1*j]   = avu[j];
572:       for (i=0; i<bs; i++) {
573:          ad_vustart[bs1*j+1+i] = 0.0;
574:       }
575:     }
576:     VecRestoreArray(localxx,&avu);

578:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
579:       nI = 0;
580:       for (stencil.k = info.zs; stencil.k<info.zs+info.zm; stencil.k++) {
581:         for (stencil.j = info.ys; stencil.j<info.ys+info.ym; stencil.j++) {
582:           for (stencil.i = info.xs; stencil.i<info.xs+info.xm; stencil.i++) {
583:             gI   = (stencil.i - info.gxs)*info.dof + (stencil.j - info.gys)*info.dof*info.gxm + (stencil.k - info.gzs)*info.dof*info.gxm*info.gym;
584:             (*NLFNewton_DAADb)(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual+nI);
585:             nI += bs;
586:           }
587:         }
588:       }
589:     }
590:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
591:       nI = info.dof*info.xm*info.ym*info.zm - bs;
592:       for (stencil.k = info.zs+info.zm-1; stencil.k>=info.zs; stencil.k--) {
593:         for (stencil.j = info.ys+info.ym-1; stencil.j>=info.ys; stencil.j--) {
594:           for (stencil.i = info.xs+info.xm-1; stencil.i>=info.xs; stencil.i--) {
595:             gI   = (stencil.i - info.gxs)*info.dof + (stencil.j - info.gys)*info.dof*info.gxm + (stencil.k - info.gzs)*info.dof*info.gxm*info.gym;
596:             (*NLFNewton_DAADb)(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual+nI);
597:             nI -= bs;
598:           }
599:         }
600:       }
601:     }

603:     /* copy solution back into ghosted vector from derivative object */
604:     VecGetArray(localxx,&av);
605:     for (j=0; j<gtdof; j++) {
606:       av[j] = ad_vustart[bs1*j];
607:     }
608:     VecRestoreArray(localxx,&av);
609:     /* stick relaxed solution back into global solution */
610:     DALocalToGlobal(A->da,localxx,INSERT_VALUES,xx);
611:   }

613:   VecRestoreArray(A->residual,&residual);
614:   DARestoreLocalVector(A->da,&localxx);
615:   DARestoreAdicMFArray(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,&gtdof);
616:   return(0);
617: }


623: PetscErrorCode NLFDestroy_DAAD(NLF A)
624: {

628:   DADestroy(A->da);
629:   PetscFree(A);
630:   return(0);
631: }

636: PetscErrorCode  NLFDAADSetDA_DAAD(NLF A,DA da)
637: {

641:   PetscObjectReference((PetscObject)da);
642:   if (A->da) {DADestroy(A->da);}
643:   A->da = da;
644:   return(0);
645: }

651: PetscErrorCode  NLFDAADSetNewtonIterations_DAAD(NLF A,int its)
652: {
654:   A->newton_its = its;
655:   return(0);
656: }

662: PetscErrorCode  NLFDAADSetResidual_DAAD(NLF A,Vec residual)
663: {
665:   A->residual = residual;
666:   return(0);
667: }


674: PetscErrorCode  NLFDAADSetCtx_DAAD(NLF A,void *ctx)
675: {
677:   A->ctx = ctx;
678:   return(0);
679: }

685: PetscErrorCode  NLFCreate_DAAD(NLF *A)
686: {

690:   PetscNew(struct NLF_DAAD,A);
691:   (*A)->da         = 0;
692:   (*A)->ctx        = 0;
693:   (*A)->newton_its = 2;
694:   return(0);
695: }