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 petsc.h
10: #include petscda.h
11: #include petscsys.h
14: #include "adic/ad_utils.h"
17: #include src/dm/da/daimpl.h
18: #include src/inline/ilu.h
20: struct NLF_DAAD {
21: DA da;
22: void *ctx;
23: Vec residual;
24: int newton_its;
25: };
27: /*
28: Solves the one dimensional equation using Newton's method
29: */
32: PetscErrorCode NLFNewton_DAAD(NLF A,DALocalInfo *info,MatStencil *stencil,void *ad_vu,PetscScalar *ad_vustart,int nI,int gI,PetscScalar residual)
33: {
35: PetscInt cnt = A->newton_its;
36: PetscScalar ad_f[2],J,f;
39: ad_vustart[1+2*gI] = 1.0;
40: do {
41: /* compute the function and Jacobian */
42: (*A->da->adicmf_lfi)(info,stencil,ad_vu,ad_f,A->ctx);
43: J = -ad_f[1];
44: f = -ad_f[0] + residual;
45: if (f != f) SETERRQ(1,"nan");
46: ad_vustart[2*gI] = ad_vustart[2*gI] - f/J;
47: } while (--cnt > 0 && PetscAbsScalar(f) > 1.e-14);
49: ad_vustart[1+2*gI] = 0.0;
50: return(0);
51: }
53: /*
54: Solves the four dimensional equation using Newton's method
55: */
58: PetscErrorCode NLFNewton_DAAD4(NLF A,DALocalInfo *info,MatStencil *stencil,void *ad_vu,PetscScalar *ad_vustart,int nI,int gI,PetscScalar *residual)
59: {
61: PetscInt cnt = A->newton_its;
62: PetscScalar ad_f[20], J[16],f[4], res, dd[5];
66: /* This sets the identity as the seed matrix for ADIC */
67: CHKMEMQ;
68: ad_vustart[1+5*gI ] = 1.0;
69: CHKMEMQ;
70: ad_vustart[2+5*gI+5 ] = 1.0;
71: CHKMEMQ;
72: ad_vustart[3+5*gI+10] = 1.0;
73: CHKMEMQ;
74: ad_vustart[4+5*gI+15] = 1.0;
75: CHKMEMQ;
77: do {
78: /* compute the function and Jacobian */
79: CHKMEMQ;
80: (*A->da->adicmf_lfib)(info,stencil,ad_vu,ad_f,A->ctx);
81: CHKMEMQ;
82: /* copy ADIC formated Jacobian into regular C array */
83: J[0] = ad_f[1] ; J[1] = ad_f[2] ; J[2] = ad_f[3] ; J[3] = ad_f[4] ;
84: J[4] = ad_f[6] ; J[5] = ad_f[7] ; J[6] = ad_f[8] ; J[7] = ad_f[9] ;
85: J[8] = ad_f[11]; J[9] = ad_f[12]; J[10]= ad_f[13]; J[11]= ad_f[14];
86: J[12]= ad_f[16]; J[13]= ad_f[17]; J[14]= ad_f[18]; J[15]= ad_f[19];
87: CHKMEMQ;
88: f[0] = -ad_f[0] + residual[0];
89: f[1] = -ad_f[5] + residual[1];
90: f[2] = -ad_f[10] + residual[2];
91: f[3] = -ad_f[15] + residual[3];
93: /* solve Jacobian * dd = ff */
95: /* could use PETSc kernel code to solve system with pivoting */
97: /* could put code in here to compute the solution directly using ADIC data structures instead of copying first */
98: 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]))-
99: 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]))+
100: 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]))-
101: 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]));
103: 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]))-
104: 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]))+
105: 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]))-
106: 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];
108: 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]))-
109: 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]))+
110: 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]))-
111: 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];
113: 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]))-
114: 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]))+
115: 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]))-
116: 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];
118: 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]))-
119: 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]))+
120: 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]))-
121: 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];
122: CHKMEMQ;
123: /* copy solution back into ADIC data structure */
124: ad_vustart[5*(gI+0)] += dd[1];
125: ad_vustart[5*(gI+1)] += dd[2];
126: ad_vustart[5*(gI+2)] += dd[3];
127: ad_vustart[5*(gI+3)] += dd[4];
128: CHKMEMQ;
129: res = f[0]*f[0];
130: res += f[1]*f[1];
131: res += f[2]*f[2];
132: res += f[3]*f[3];
133: res = sqrt(res);
134:
135: } while (--cnt > 0 && res > 1.e-14);
137: /* zero out this part of the seed matrix that was set initially */
138: ad_vustart[1+5*gI ] = 0.0;
139: ad_vustart[2+5*gI+5 ] = 0.0;
140: ad_vustart[3+5*gI+10] = 0.0;
141: ad_vustart[4+5*gI+15] = 0.0;
142: CHKMEMQ;
143: return(0);
144: }
148: PetscErrorCode NLFNewton_DAAD9(NLF A,DALocalInfo *info,MatStencil *stencil,void *ad_vu,PetscScalar *ad_vustart,int nI,int gI,PetscScalar *residual)
149: {
151: PetscInt cnt = A->newton_its;
152: PetscScalar ad_f[100], J[81],f[9], res;
153: PetscInt i,j,ngI[9];
155:
156: // the order of the nodes
157: /*
158: (6) (7) (8)
159: i-1,j+1 --- i,j+1 --- i+1,j+1
160: | | |
161: | | |
162: i-1,j --- i,j --- i+1,j
163: |(3) |(4) |(5)
164: | | |
165: i-1,j-1 --- i,j-1--- i+1,j-1
166: (0) (1) (2)
167: */
168:
169: // the order of the derivative for the center nodes
170: /*
171: (7) (8) (9)
172: i-1,j+1 --- i,j+1 --- i+1,j+1
173: | | |
174: | | |
175: i-1,j --- i,j --- i+1,j
176: |(4) |(5) |(6)
177: | | |
178: i-1,j-1 --- i,j-1--- i+1,j-1
179: (1) (2) (3)
180: */
181: 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) {
183: ad_vustart[1+10*gI] = 1.0;
184:
185: do {
186: /* compute the function and Jacobian */
187: (*A->da->adicmf_lfi)(info,stencil,ad_vu,ad_f,A->ctx);
188: J[0] = -ad_f[1];
189: f[0] = -ad_f[0] + residual[gI];
190: ad_vustart[10*gI] = ad_vustart[10*gI] - f[0]/J[0];
191: } while (--cnt > 0 && PetscAbsScalar(f[0]) > 1.e-14);
193: ad_vustart[1+10*gI] = 0.0;
194: return(0);
197: }
198:
199: 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;
200: ngI[1] = ngI[0] + 1;
201: ngI[2] = ngI[1] + 1;
202: ngI[3] = gI - 1;
203: ngI[4] = gI ;
204: ngI[5] = gI + 1;
205: 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;
206: ngI[7] = ngI[6] + 1;
207: ngI[8] = ngI[7] + 1;
208:
210: for(j=0 ; j<9; j++){
211: ad_vustart[ngI[j]*10+j+1] = 1.0;
212: }
213:
214: do{
215: /* compute the function and the Jacobian */
216:
217: (*A->da->adicmf_lfi)(info,stencil,ad_vu,ad_f,A->ctx);
218:
219: for(i=0; i<9; i++){
220: for(j=0; j<9; j++){
221: J[i*9+j] = -ad_f[i*10+j+1];
222: }
223: f[i]= -ad_f[i*10] + residual[ngI[i]];
224: }
226: Kernel_A_gets_inverse_A_9(J);
227:
228: res =0 ;
229: for(i=0; i<9; i++){
230: for(j=0;j<9;j++){
231: ad_vustart[10*ngI[i]]= ad_vustart[10*ngI[i]] - J[i*9 +j]*f[j];
232: }
233: res = res + f[i]*f[i];
234: }
235: res = sqrt(res);
236: } while(--cnt>0 && res>1.e-14);
238: for(j=0; j<9; j++){
239: ad_vustart[10*ngI[j]+j+1]=0.0;
240: }
242: return(0);
243: }
245: /*
246: Nonlinear relax on all the equations with an initial guess in x
247: */
251: PetscErrorCode NLFRelax_DAAD(NLF A,MatSORType flag,int its,Vec xx)
252: {
254: PetscInt j,gtdof,nI,gI;
255: PetscScalar *avu,*av,*ad_vustart,*residual;
256: Vec localxx;
257: DALocalInfo info;
258: MatStencil stencil;
259: void* *ad_vu;
262: if (its <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D positive",its);
264: DAGetLocalVector(A->da,&localxx);
265: /* get space for derivative object. */
266: DAGetAdicMFArray(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,>dof);
267: VecGetArray(A->residual,&residual);
270: /* tell ADIC we will be computing one dimensional Jacobians */
271: PetscADResetIndep();
272: PetscADIncrementTotalGradSize(1);
273: PetscADSetIndepDone();
275: DAGetLocalInfo(A->da,&info);
276: while (its--) {
278: /* get initial solution properly ghosted */
279: DAGlobalToLocalBegin(A->da,xx,INSERT_VALUES,localxx);
280: DAGlobalToLocalEnd(A->da,xx,INSERT_VALUES,localxx);
282: /* copy input vector into derivative object */
283: VecGetArray(localxx,&avu);
284: for (j=0; j<gtdof; j++) {
285: ad_vustart[2*j] = avu[j];
286: ad_vustart[2*j+1] = 0.0;
287: }
288:
289: VecRestoreArray(localxx,&avu);
291: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
292: nI = 0;
293: for (stencil.k = info.zs; stencil.k<info.zs+info.zm; stencil.k++) {
294: for (stencil.j = info.ys; stencil.j<info.ys+info.ym; stencil.j++) {
295: for (stencil.i = info.xs; stencil.i<info.xs+info.xm; stencil.i++) {
296: for (stencil.c = 0; stencil.c<info.dof; stencil.c++) {
297: 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;
298: NLFNewton_DAAD(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual[nI]);
299: nI++;
300: }
301: }
302: }
303: }
304: }
305: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
306: nI = info.dof*info.xm*info.ym*info.zm - 1;
307: for (stencil.k = info.zs+info.zm-1; stencil.k>=info.zs; stencil.k--) {
308: for (stencil.j = info.ys+info.ym-1; stencil.j>=info.ys; stencil.j--) {
309: for (stencil.i = info.xs+info.xm-1; stencil.i>=info.xs; stencil.i--) {
310: for (stencil.c = info.dof-1; stencil.c>=0; stencil.c--) {
311: 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;
312: NLFNewton_DAAD(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual[nI]);
313: nI--;
314: }
315: }
316: }
317: }
318: }
320: /* copy solution back into ghosted vector from derivative object */
321: VecGetArray(localxx,&av);
322: for (j=0; j<gtdof; j++) {
323: av[j] = ad_vustart[2*j];
324: }
325: VecRestoreArray(localxx,&av);
326: /* stick relaxed solution back into global solution */
327: DALocalToGlobal(A->da,localxx,INSERT_VALUES,xx);
328: }
331: VecRestoreArray(A->residual,&residual);
332: DARestoreLocalVector(A->da,&localxx);
333: DARestoreAdicMFArray(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,>dof);
334: return(0);
335: }
341: PetscErrorCode NLFRelax_DAAD4(NLF A,MatSORType flag,int its,Vec xx)
342: {
344: PetscInt j,gtdof,nI,gI;
345: PetscScalar *avu,*av,*ad_vustart,*residual;
346: Vec localxx;
347: DALocalInfo info;
348: MatStencil stencil;
349: void* *ad_vu;
352: if (its <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D positive",its);
353:
354: DAGetLocalVector(A->da,&localxx);
355: /* get space for derivative object. */
356: DAGetAdicMFArray4(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,>dof);
357: VecGetArray(A->residual,&residual);
360: /* tell ADIC we will be computing four dimensional Jacobians */
361: PetscADResetIndep();
362: PetscADIncrementTotalGradSize(4);
363: PetscADSetIndepDone();
365: DAGetLocalInfo(A->da,&info);
366: while (its--) {
368: /* get initial solution properly ghosted */
369: DAGlobalToLocalBegin(A->da,xx,INSERT_VALUES,localxx);
370: DAGlobalToLocalEnd(A->da,xx,INSERT_VALUES,localxx);
372: /* copy input vector into derivative object */
373: VecGetArray(localxx,&avu);
374: for (j=0; j<gtdof; j++) {
375: ad_vustart[5*j ] = avu[j];
376: ad_vustart[5*j+1] = 0.0;
377: ad_vustart[5*j+2] = 0.0;
378: ad_vustart[5*j+3] = 0.0;
379: ad_vustart[5*j+4] = 0.0;
380: }
381:
382: VecRestoreArray(localxx,&avu);
384: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
385: nI = 0;
386: for (stencil.k = info.zs; stencil.k<info.zs+info.zm; stencil.k++) {
387: for (stencil.j = info.ys; stencil.j<info.ys+info.ym; stencil.j++) {
388: for (stencil.i = info.xs; stencil.i<info.xs+info.xm; stencil.i++) {
389: CHKMEMQ;
390: 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;
392: NLFNewton_DAAD4(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual+nI);
393: nI=nI+4;
394: CHKMEMQ;
395: }
396: }
397: }
398: }
399: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
400: nI = info.dof*info.xm*info.ym*info.zm - 4;
401:
402: for (stencil.k = info.zs+info.zm-1; stencil.k>=info.zs; stencil.k--) {
403: for (stencil.j = info.ys+info.ym-1; stencil.j>=info.ys; stencil.j--) {
404: for (stencil.i = info.xs+info.xm-1; stencil.i>=info.xs; stencil.i--) {
405: 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;
406: NLFNewton_DAAD4(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual+nI);
407: nI=nI-4;
408: }
409: }
410: }
411: }
412:
413: /* copy solution back into ghosted vector from derivative object */
414: VecGetArray(localxx,&av);
415: for (j=0; j<gtdof; j++) {
416: av[j] = ad_vustart[5*j];
417: }
418: VecRestoreArray(localxx,&av);
419: /* stick relaxed solution back into global solution */
420: DALocalToGlobal(A->da,localxx,INSERT_VALUES,xx);
421: }
424: VecRestoreArray(A->residual,&residual);
425: DARestoreLocalVector(A->da,&localxx);
426: DARestoreAdicMFArray(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,>dof);
427: return(0);
428: }
433: PetscErrorCode NLFRelax_DAAD9(NLF A,MatSORType flag,int its,Vec xx)
434: {
436: PetscInt j,gtdof,nI,gI;
437: PetscScalar *avu,*av,*ad_vustart,*residual;
438: Vec localxx;
439: DALocalInfo info;
440: MatStencil stencil;
441: void* *ad_vu;
444: if (its <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D positive",its);
445:
446: DAGetLocalVector(A->da,&localxx);
447: /* get space for derivative object. */
448: DAGetAdicMFArray9(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,>dof);
449: VecGetArray(A->residual,&residual);
452: /* tell ADIC we will be computing nine dimensional Jacobians */
453: PetscADResetIndep();
454: PetscADIncrementTotalGradSize(9);
455: PetscADSetIndepDone();
457: DAGetLocalInfo(A->da,&info);
458: while (its--) {
460: /* get initial solution properly ghosted */
461: DAGlobalToLocalBegin(A->da,xx,INSERT_VALUES,localxx);
462: DAGlobalToLocalEnd(A->da,xx,INSERT_VALUES,localxx);
464: /* copy input vector into derivative object */
465: VecGetArray(localxx,&avu);
466: for (j=0; j<gtdof; j++) {
467: ad_vustart[10*j ] = avu[j];
468: ad_vustart[10*j+1] = 0.0;
469: ad_vustart[10*j+2] = 0.0;
470: ad_vustart[10*j+3] = 0.0;
471: ad_vustart[10*j+4] = 0.0;
472: ad_vustart[10*j+5] = 0.0;
473: ad_vustart[10*j+6] = 0.0;
474: ad_vustart[10*j+7] = 0.0;
475: ad_vustart[10*j+8] = 0.0;
476: ad_vustart[10*j+9] = 0.0;
477: }
478:
479: VecRestoreArray(localxx,&avu);
481: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
482: nI = 0;
483: for (stencil.k = info.zs; stencil.k<info.zs+info.zm; stencil.k++) {
484: for (stencil.j = info.ys; stencil.j<info.ys+info.ym; stencil.j++) {
485: for (stencil.i = info.xs; stencil.i<info.xs+info.xm; stencil.i++) {
486: CHKMEMQ;
487: 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;
488: NLFNewton_DAAD9(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual);
489: nI=nI+1;
490: CHKMEMQ;
491: }
492: }
493: }
494: }
495: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
496: nI = info.dof*info.xm*info.ym*info.zm - 4;
497:
498: for (stencil.k = info.zs+info.zm-1; stencil.k>=info.zs; stencil.k--) {
499: for (stencil.j = info.ys+info.ym-1; stencil.j>=info.ys; stencil.j--) {
500: for (stencil.i = info.xs+info.xm-1; stencil.i>=info.xs; stencil.i--) {
501: 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;
502: NLFNewton_DAAD9(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual);
503: nI=nI-1;
504: }
505: }
506: }
507: }
508:
509: /* copy solution back into ghosted vector from derivative object */
510: VecGetArray(localxx,&av);
511: for (j=0; j<gtdof; j++) {
512: av[j] = ad_vustart[10*j];
513: }
514: VecRestoreArray(localxx,&av);
515: /* stick relaxed solution back into global solution */
516: DALocalToGlobal(A->da,localxx,INSERT_VALUES,xx);
517: }
520: VecRestoreArray(A->residual,&residual);
521: DARestoreLocalVector(A->da,&localxx);
522: DARestoreAdicMFArray(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,>dof);
523: return(0);
524: }
526: /*
527: Point-block nonlinear relax on all the equations with an initial guess in xx using
528: */
532: PetscErrorCode NLFRelax_DAADb(NLF A,MatSORType flag,int its,Vec xx)
533: {
535: PetscInt i,j,gtdof,nI,gI, bs = A->da->w, bs1 = bs + 1;
536: PetscScalar *avu,*av,*ad_vustart,*residual;
537: Vec localxx;
538: DALocalInfo info;
539: MatStencil stencil;
540: void* *ad_vu;
541: PetscErrorCode (*NLFNewton_DAADb)(NLF,DALocalInfo*,MatStencil*,void*,PetscScalar*,int,int,PetscScalar*);
544: if (its <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D positive",its);
545: if (bs == 4) {
546: NLFNewton_DAADb = NLFNewton_DAAD4;
547: } else {
548: SETERRQ1(PETSC_ERR_SUP,"Point block nonlinear relaxation currently not for this block size",bs);
549: }
551: DAGetLocalVector(A->da,&localxx);
552: /* get space for derivative object. */
553: DAGetAdicMFArrayb(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,>dof);
554: VecGetArray(A->residual,&residual);
557: /* tell ADIC we will be computing bs dimensional Jacobians */
558: PetscADResetIndep();
559: PetscADIncrementTotalGradSize(bs);
560: PetscADSetIndepDone();
562: DAGetLocalInfo(A->da,&info);
563: while (its--) {
565: /* get initial solution properly ghosted */
566: DAGlobalToLocalBegin(A->da,xx,INSERT_VALUES,localxx);
567: DAGlobalToLocalEnd(A->da,xx,INSERT_VALUES,localxx);
569: /* copy input vector into derivative object */
570: VecGetArray(localxx,&avu);
571: for (j=0; j<gtdof; j++) {
572: ad_vustart[bs1*j] = avu[j];
573: for (i=0; i<bs; i++) {
574: ad_vustart[bs1*j+1+i] = 0.0;
575: }
576: }
577: VecRestoreArray(localxx,&avu);
579: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
580: nI = 0;
581: for (stencil.k = info.zs; stencil.k<info.zs+info.zm; stencil.k++) {
582: for (stencil.j = info.ys; stencil.j<info.ys+info.ym; stencil.j++) {
583: for (stencil.i = info.xs; stencil.i<info.xs+info.xm; stencil.i++) {
584: 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;
585: (*NLFNewton_DAADb)(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual+nI);
586: nI += bs;
587: }
588: }
589: }
590: }
591: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
592: nI = info.dof*info.xm*info.ym*info.zm - bs;
593: for (stencil.k = info.zs+info.zm-1; stencil.k>=info.zs; stencil.k--) {
594: for (stencil.j = info.ys+info.ym-1; stencil.j>=info.ys; stencil.j--) {
595: for (stencil.i = info.xs+info.xm-1; stencil.i>=info.xs; stencil.i--) {
596: 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;
597: (*NLFNewton_DAADb)(A,&info,&stencil,ad_vu,ad_vustart,nI,gI,residual+nI);
598: nI -= bs;
599: }
600: }
601: }
602: }
604: /* copy solution back into ghosted vector from derivative object */
605: VecGetArray(localxx,&av);
606: for (j=0; j<gtdof; j++) {
607: av[j] = ad_vustart[bs1*j];
608: }
609: VecRestoreArray(localxx,&av);
610: /* stick relaxed solution back into global solution */
611: DALocalToGlobal(A->da,localxx,INSERT_VALUES,xx);
612: }
614: VecRestoreArray(A->residual,&residual);
615: DARestoreLocalVector(A->da,&localxx);
616: DARestoreAdicMFArray(A->da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,>dof);
617: return(0);
618: }
624: PetscErrorCode NLFDestroy_DAAD(NLF A)
625: {
629: DADestroy(A->da);
630: PetscFree(A);
631: return(0);
632: }
637: PetscErrorCode NLFDAADSetDA_DAAD(NLF A,DA da)
638: {
642: A->da = da;
643: PetscObjectReference((PetscObject)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: }