Actual source code: damgsnesad.c
1:
2: #include <petscdmda.h> /*I "petscdmda.h" I*/
3: #include <petscpcmg.h> /*I "petscpcmg.h" I*/
4: #include <petscdmmg.h> /*I "petscdmmg.h" I*/
5: #include <../src/mat/blockinvert.h>
6: #include <../src/snes/impls/ls/lsimpl.h>
19: /*
20: DMMGComputeJacobianWithAdic - Evaluates the Jacobian via Adic when the user has provided
21: a local function evaluation routine.
22: */
23: PetscErrorCode DMMGComputeJacobianWithAdic(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
24: {
25: DMMG dmmg = (DMMG) ptr;
27: Vec localX;
28: DM da = dmmg->dm;
31: DMGetLocalVector(da,&localX);
32: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
33: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
34: DMDAComputeJacobian1WithAdic(da,localX,*B,dmmg->user);
35: DMRestoreLocalVector(da,&localX);
36: /* Assemble true Jacobian; if it is different */
37: if (*J != *B) {
38: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
39: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
40: }
41: MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
42: *flag = SAME_NONZERO_PATTERN;
43: return(0);
44: }
48: /*@
49: SNESDAComputeJacobianWithAdic - This is a universal Jacobian evaluation routine
50: that may be used with SNESSetJacobian() as long as the user context has a DMDA as
51: its first record and DMDASetLocalAdicFunction() has been called.
53: Collective on SNES
55: Input Parameters:
56: + snes - the SNES context
57: . X - input vector
58: . J - Jacobian
59: . B - Jacobian used in preconditioner (usally same as J)
60: . flag - indicates if the matrix changed its structure
61: - ptr - optional user-defined context, as set by SNESSetFunction()
63: Level: intermediate
65: .seealso: DMDASetLocalFunction(), DMDASetLocalAdicFunction(), SNESSetFunction(), SNESSetJacobian()
67: @*/
68: PetscErrorCode SNESDAComputeJacobianWithAdic(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
69: {
70: DM da = *(DM*) ptr;
72: Vec localX;
75: DMGetLocalVector(da,&localX);
76: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
77: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
78: DMDAComputeJacobian1WithAdic(da,localX,*B,ptr);
79: DMRestoreLocalVector(da,&localX);
80: /* Assemble true Jacobian; if it is different */
81: if (*J != *B) {
82: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
83: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
84: }
85: MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
86: *flag = SAME_NONZERO_PATTERN;
87: return(0);
88: }
90: #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscpcmg.h" I*/
91: /*
92: This is pre-beta FAS code. It's design should not be taken seriously!
94: R is the usual multigrid restriction (e.g. the tranpose of piecewise linear interpolation)
95: Q is either a scaled injection or the usual R
96: */
99: PetscErrorCode DMMGSolveFAS(DMMG *dmmg,PetscInt level)
100: {
102: PetscInt i,j,k;
103: PetscReal norm;
104: PC_MG **mg;
105: PC pc;
108: VecSet(dmmg[level]->r,0.0);
109: for (j=1; j<=level; j++) {
110: if (!dmmg[j]->inject) {
111: DMGetInjection(dmmg[j-1]->dm,dmmg[j]->dm,&dmmg[j]->inject);
112: }
113: }
115: KSPGetPC(dmmg[level]->ksp,&pc);
116: mg = ((PC_MG**)pc->data);
118: for(i = 0; i < 100; i++) {
120: for(j = level; j > 0; j--) {
122: /* Relax residual_fine --> F(x_fine) = 0 */
123: for(k = 0; k < dmmg[j]->presmooth; k++) {
124: NLFRelax_DAAD(dmmg[j]->nlf, SOR_SYMMETRIC_SWEEP, 1, dmmg[j]->x);
125: }
127: /* R*(residual_fine - F(x_fine)) */
128: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
129: VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);
131: if (j == level || dmmg[j]->monitorall) {
132: /* norm( residual_fine - f(x_fine) ) */
133: VecNorm(dmmg[j]->w,NORM_2,&norm);
134: if (j == level) {
135: if (norm < dmmg[level]->abstol) goto theend;
136: if (i == 0) {
137: dmmg[level]->rrtol = norm*dmmg[level]->rtol;
138: } else {
139: if (norm < dmmg[level]->rrtol) goto theend;
140: }
141: }
142: if (dmmg[j]->monitorall) {
143: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
144: PetscPrintf(dmmg[j]->comm,"FAS function norm %G\n",norm);
145: }
146: }
148: MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);
149:
150: /* F(Q*x_fine) */
151: VecScatterBegin(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
152: VecScatterEnd(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
153: DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);
155: /* residual_coarse = F(Q*x_fine) + R*(residual_fine - F(x_fine)) */
156: VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);
158: /* save Q*x_fine into b (needed when interpolating compute x back up */
159: VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);
160: }
162: for (j=0; j<dmmg[0]->coarsesmooth; j++) {
163: NLFRelax_DAAD(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);
164: }
165: if (dmmg[0]->monitorall){
166: DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);
167: VecAXPY(dmmg[0]->w,-1.0,dmmg[0]->r);
168: VecNorm(dmmg[0]->w,NORM_2,&norm);
169: for (k=0; k<level+1; k++) {PetscPrintf(dmmg[0]->comm," ");}
170: PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %G\n",norm);
171: }
173: for (j=1; j<=level; j++) {
174: /* x_fine = x_fine + R'*(x_coarse - Q*x_fine) */
175: VecAXPY(dmmg[j-1]->x,-1.0,dmmg[j-1]->b);
176: MatInterpolateAdd(mg[j]->interpolate,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);
178: if (dmmg[j]->monitorall) {
179: /* norm( F(x_fine) - residual_fine ) */
180: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
181: VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
182: VecNorm(dmmg[j]->w,NORM_2,&norm);
183: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
184: PetscPrintf(dmmg[j]->comm,"FAS function norm %G\n",norm);
185: }
187: /* Relax residual_fine - F(x_fine) = 0 */
188: for (k=0; k<dmmg[j]->postsmooth; k++) {
189: NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
190: }
192: if (dmmg[j]->monitorall) {
193: /* norm( F(x_fine) - residual_fine ) */
194: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
195: VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
196: VecNorm(dmmg[j]->w,NORM_2,&norm);
197: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
198: PetscPrintf(dmmg[j]->comm,"FAS function norm %G\n",norm);
199: }
200: }
202: if (dmmg[level]->monitor){
203: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
204: VecNorm(dmmg[level]->w,NORM_2,&norm);
205: PetscPrintf(dmmg[level]->comm,"%D FAS function norm %G\n",i+1,norm);
206: }
207: }
208: theend:
209: return(0);
210: }
212: /*
213: This is the point-block version of FAS
214: */
217: PetscErrorCode DMMGSolveFASb(DMMG *dmmg,PetscInt level)
218: {
220: PetscInt i,j,k;
221: PetscReal norm;
222: PC_MG **mg;
223: PC pc;
226: VecSet(dmmg[level]->r,0.0);
227: for (j=1; j<=level; j++) {
228: if (!dmmg[j]->inject) {
229: DMGetInjection(dmmg[j-1]->dm,dmmg[j]->dm,&dmmg[j]->inject);
230: }
231: }
233: KSPGetPC(dmmg[level]->ksp,&pc);
234: mg = ((PC_MG**)pc->data);
236: for (i=0; i<100; i++) {
238: for (j=level; j>0; j--) {
240: /* Relax residual_fine - F(x_fine) = 0 */
241: for (k=0; k<dmmg[j]->presmooth; k++) {
242: NLFRelax_DAADb(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
243: }
245: /* R*(residual_fine - F(x_fine)) */
246: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
247: VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);
249: if (j == level || dmmg[j]->monitorall) {
250: /* norm( residual_fine - f(x_fine) ) */
251: VecNorm(dmmg[j]->w,NORM_2,&norm);
252: if (j == level) {
253: if (norm < dmmg[level]->abstol) goto theend;
254: if (i == 0) {
255: dmmg[level]->rrtol = norm*dmmg[level]->rtol;
256: } else {
257: if (norm < dmmg[level]->rrtol) goto theend;
258: }
259: }
260: }
262: if (dmmg[j]->monitorall) {
263: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
264: PetscPrintf(dmmg[j]->comm,"FAS function norm %G\n",norm);
265: }
266: MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);
267:
268: /* F(Q*x_fine) */
269: VecScatterBegin(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
270: VecScatterEnd(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
271: DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);
273: /* residual_coarse = F(Q*x_fine) + R*(residual_fine - F(x_fine)) */
274: VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);
276: /* save Q*x_fine into b (needed when interpolating compute x back up */
277: VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);
278: }
280: for (j=0; j<dmmg[0]->coarsesmooth; j++) {
281: NLFRelax_DAADb(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);
282: }
283: if (dmmg[0]->monitorall){
284: DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);
285: VecAXPY(dmmg[0]->w,-1.0,dmmg[0]->r);
286: VecNorm(dmmg[0]->w,NORM_2,&norm);
287: for (k=0; k<level+1; k++) {PetscPrintf(dmmg[0]->comm," ");}
288: PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %G\n",norm);
289: }
291: for (j=1; j<=level; j++) {
292: /* x_fine = x_fine + R'*(x_coarse - Q*x_fine) */
293: VecAXPY(dmmg[j-1]->x,-1.0,dmmg[j-1]->b);
294: MatInterpolateAdd(mg[j]->interpolate,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);
296: if (dmmg[j]->monitorall) {
297: /* norm( F(x_fine) - residual_fine ) */
298: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
299: VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
300: VecNorm(dmmg[j]->w,NORM_2,&norm);
301: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
302: PetscPrintf(dmmg[j]->comm,"FAS function norm %G\n",norm);
303: }
305: /* Relax residual_fine - F(x_fine) = 0 */
306: for (k=0; k<dmmg[j]->postsmooth; k++) {
307: NLFRelax_DAADb(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
308: }
310: if (dmmg[j]->monitorall) {
311: /* norm( F(x_fine) - residual_fine ) */
312: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
313: VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
314: VecNorm(dmmg[j]->w,NORM_2,&norm);
315: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
316: PetscPrintf(dmmg[j]->comm,"FAS function norm %G\n",norm);
317: }
318: }
320: if (dmmg[level]->monitor){
321: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
322: VecNorm(dmmg[level]->w,NORM_2,&norm);
323: PetscPrintf(dmmg[level]->comm,"%D FAS function norm %G\n",i+1,norm);
324: }
325: }
326: theend:
327: return(0);
328: }
331: #include <adic/ad_utils.h>
336: PetscErrorCode PetscADView(PetscInt N,PetscInt nc,double *ptr,PetscViewer viewer)
337: {
338: PetscInt i,j,nlen = PetscADGetDerivTypeSize();
339: char *cptr = (char*)ptr;
340: double *values;
344: for (i=0; i<N; i++) {
345: PetscPrintf(PETSC_COMM_SELF,"Element %D value %G derivatives: ",i,*(double*)cptr);
346: values = PetscADGetGradArray(cptr);
347: for (j=0; j<nc; j++) {
348: PetscPrintf(PETSC_COMM_SELF,"%G ",*values++);
349: }
350: PetscPrintf(PETSC_COMM_SELF,"\n");
351: cptr += nlen;
352: }
354: return(0);
355: }
359: PetscErrorCode DMMGSolveFAS4(DMMG *dmmg,PetscInt level)
360: {
362: PetscInt i,j,k;
363: PetscReal norm;
364: PetscScalar zero = 0.0,mone = -1.0,one = 1.0;
365: PC_MG **mg;
366: PC pc;
369: VecSet(dmmg[level]->r,zero);
370: for (j=1; j<=level; j++) {
371: if (!dmmg[j]->inject) {
372: DMGetInjection(dmmg[j-1]->dm,dmmg[j]->dm,&dmmg[j]->inject);
373: }
374: }
376: KSPGetPC(dmmg[level]->ksp,&pc);
377: mg = ((PC_MG**)pc->data);
378: for (i=0; i<100; i++) {
380: for (j=level; j>0; j--) {
381: PetscPrintf(PETSC_COMM_WORLD,"I am here");
382: /* Relax residual_fine - F(x_fine) = 0 */
383: for (k=0; k<dmmg[j]->presmooth; k++) {
384: NLFRelax_DAAD4(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
385: }
387: /* R*(residual_fine - F(x_fine)) */
388: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
389: VecAYPX(dmmg[j]->w,mone,dmmg[j]->r);
391: if (j == level || dmmg[j]->monitorall) {
392: /* norm( residual_fine - f(x_fine) ) */
393: VecNorm(dmmg[j]->w,NORM_2,&norm);
394: if (j == level) {
395: if (norm < dmmg[level]->abstol) goto theend;
396: if (i == 0) {
397: dmmg[level]->rrtol = norm*dmmg[level]->rtol;
398: } else {
399: if (norm < dmmg[level]->rrtol) goto theend;
400: }
401: }
402: } PetscPrintf(PETSC_COMM_WORLD,"I am here");
404: if (dmmg[j]->monitorall) {
405: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
406: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
407: }
408: MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);
409:
410: /* F(R*x_fine) */
411: VecScatterBegin(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
412: VecScatterEnd(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
413: DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);
415: /* residual_coarse = F(R*x_fine) + R*(residual_fine - F(x_fine)) */
416: VecAYPX(dmmg[j-1]->r,one,dmmg[j-1]->w);
418: /* save R*x_fine into b (needed when interpolating compute x back up */
419: VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);
420: }
422: for (j=0; j<dmmg[0]->presmooth; j++) {
423: NLFRelax_DAAD4(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);
424: }
425: if (dmmg[0]->monitorall){
426: DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);
427: VecAXPY(dmmg[0]->w,mone,dmmg[0]->r);
428: VecNorm(dmmg[0]->w,NORM_2,&norm);
429: for (k=0; k<level+1; k++) {PetscPrintf(dmmg[0]->comm," ");}
430: PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %g\n",norm);
431: }
433: for (j=1; j<=level; j++) {
434: /* x_fine = x_fine + R'*(x_coarse - R*x_fine) */
435: VecAXPY(dmmg[j-1]->x,mone,dmmg[j-1]->b);
436: MatInterpolateAdd(mg[j]->interpolate,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);
438: if (dmmg[j]->monitorall) {
439: /* norm( F(x_fine) - residual_fine ) */
440: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
441: VecAXPY(dmmg[j]->w,mone,dmmg[j]->r);
442: VecNorm(dmmg[j]->w,NORM_2,&norm);
443: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
444: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
445: }
447: /* Relax residual_fine - F(x_fine) = 0 */
448: for (k=0; k<dmmg[j]->postsmooth; k++) {
449: NLFRelax_DAAD4(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
450: }
452: if (dmmg[j]->monitorall) {
453: /* norm( F(x_fine) - residual_fine ) */
454: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
455: VecAXPY(dmmg[j]->w,mone,dmmg[j]->r);
456: VecNorm(dmmg[j]->w,NORM_2,&norm);
457: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
458: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
459: }
460: }
462: if (dmmg[level]->monitor){
463: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
464: VecNorm(dmmg[level]->w,NORM_2,&norm);
465: PetscPrintf(dmmg[level]->comm,"%D FAS function norm %g\n",i,norm);
466: }
467: }
468: theend:
469: return(0);
470: }
472: /*
473: This function provide several FAS v_cycle iteration
475: iter: the number of FAS it run
477: */
480: PetscErrorCode DMMGSolveFASn(DMMG *dmmg,PetscInt level,PetscInt iter)
481: {
483: PetscInt i,j,k;
484: PetscReal norm;
485: PC_MG **mg;
486: PC pc;
489: VecSet(dmmg[level]->r,0.0);
490: for (j=1; j<=level; j++) {
491: if (!dmmg[j]->inject) {
492: DMGetInjection(dmmg[j-1]->dm,dmmg[j]->dm,&dmmg[j]->inject);
493: }
494: }
496: KSPGetPC(dmmg[level]->ksp,&pc);
497: mg = ((PC_MG**)pc->data);
499: for (i=0; i<iter; i++) {
501: for (j=level; j>0; j--) {
503: /* Relax residual_fine - F(x_fine) = 0 */
504: for (k=0; k<dmmg[j]->presmooth; k++) {
505: NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
506: }
508: /* R*(residual_fine - F(x_fine)) */
509: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
510: VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);
512: if (j == level || dmmg[j]->monitorall) {
513: /* norm( residual_fine - f(x_fine) ) */
514: VecNorm(dmmg[j]->w,NORM_2,&norm);
515: if (j == level) {
516: if (norm < dmmg[level]->abstol) goto theend;
517: if (i == 0) {
518: dmmg[level]->rrtol = norm*dmmg[level]->rtol;
519: } else {
520: if (norm < dmmg[level]->rrtol) goto theend;
521: }
522: }
523: }
525: if (dmmg[j]->monitorall) {
526: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
527: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
528: }
529: MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);
530:
531: /* F(RI*x_fine) */
532: VecScatterBegin(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
533: VecScatterEnd(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
534: DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);
536: /* residual_coarse = F(RI*x_fine) + R*(residual_fine - F(x_fine)) */
537: VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);
539: /* save RI*x_fine into b (needed when interpolating compute x back up */
540: VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);
541: }
543: for (j=0; j<dmmg[0]->coarsesmooth; j++) {
544: NLFRelax_DAAD(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);
545: }
546: if (dmmg[0]->monitorall){
547: DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);
548: VecAXPY(dmmg[0]->w,-1.0,dmmg[0]->r);
549: VecNorm(dmmg[0]->w,NORM_2,&norm);
550: for (k=0; k<level+1; k++) {PetscPrintf(dmmg[0]->comm," ");}
551: PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %g\n",norm);
552: }
554: for (j=1; j<=level; j++) {
555: /* x_fine = x_fine + R'*(x_coarse - RI*x_fine) */
556: VecAXPY(dmmg[j-1]->x,-1.0,dmmg[j-1]->b);
557: MatInterpolateAdd(mg[j]->interpolate,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);
559: if (dmmg[j]->monitorall) {
560: /* norm( F(x_fine) - residual_fine ) */
561: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
562: VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
563: VecNorm(dmmg[j]->w,NORM_2,&norm);
564: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
565: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
566: }
568: /* Relax residual_fine - F(x_fine) = 0 */
569: for (k=0; k<dmmg[j]->postsmooth; k++) {
570: NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
571: }
573: if (dmmg[j]->monitorall) {
574: /* norm( F(x_fine) - residual_fine ) */
575: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
576: VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
577: VecNorm(dmmg[j]->w,NORM_2,&norm);
578: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
579: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
580: }
581: }
583: if (dmmg[level]->monitor){
584: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
585: VecNorm(dmmg[level]->w,NORM_2,&norm);
586: PetscPrintf(dmmg[level]->comm,"%D FAS function norm %g\n",i+1,norm);
587: }
588: }
589: theend:
590: return(0);
591: }
592: /*
593: This is a simple FAS setup function
594: */
599: PetscErrorCode DMMGSolveFASSetUp(DMMG *dmmg,PetscInt level)
600: {
602: PetscInt j;//,nlevels=dmmg[0]->nlevels-1;
603: //PC pc;
606: VecSet(dmmg[level]->r,0.0);
607: for (j=1; j<=level; j++) {
608: if (!dmmg[j]->inject) {
609: DMGetInjection(dmmg[j-1]->dm,dmmg[j]->dm,&dmmg[j]->inject);
610: }
611: }
612: VecSet(dmmg[level]->r,0.0);
613: dmmg[level]->rrtol = 0.0001*dmmg[level]->rtol;//I want to get more precise solution with FAS
614: return(0);
615: }
618: /*
619: This is function to implement multiplicative FAS
622: Options:
623:
624: -dmmg_fas_cycles 1 : FAS v-cycle
625: 2 : FAS w-cycle
628: */
632: PetscErrorCode DMMGSolveFASMCycle(DMMG *dmmg,PetscInt level,PetscBool * converged)
633: {
635: PetscInt j,k,cycles=1,nlevels=level;//nlevels=dmmg[0]->nlevels-1;
636: // I need to put nlevels=level in order to get grid sequence correctly
637: PetscReal norm;
638: PC_MG **mg;
639: PC pc;
640:
642:
644: PetscOptionsGetInt(PETSC_NULL,"-dmmg_fas_cycles",&cycles,PETSC_NULL);
645:
646: KSPGetPC(dmmg[level]->ksp,&pc);
647: mg = ((PC_MG**)pc->data);
649: j=level;
651: if(j) {/* not the coarsest level */
652: /* Relax residual_fine - F(x_fine) = 0 */
653: for (k=0; k<dmmg[j]->presmooth; k++) {
654: NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
655: }
656:
657:
659: /* R*(residual_fine - F(x_fine)) */
660: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
661: VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);
663: if (j == nlevels || dmmg[j]->monitorall) {
664: /* norm( residual_fine - f(x_fine) ) */
665: VecNorm(dmmg[j]->w,NORM_2,&norm);
666:
667: if (j == nlevels) {
668: if (norm < dmmg[level]->abstol) {
669: *converged = PETSC_TRUE;
670: goto theend;
671: }
673: if (norm < dmmg[level]->rrtol){
674: *converged = PETSC_TRUE;
675: goto theend;
676:
677: }
678: }
679: }
681: if (dmmg[j]->monitorall) {
682: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
683: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
684: }
685: MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);
686:
687: /* F(RI*x_fine) */
688: VecScatterBegin(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
689: VecScatterEnd(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
690: DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);
692: /* residual_coarse = F(RI*x_fine) + R*(residual_fine - F(x_fine)) */
693: VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);
695: /* save RI*x_fine into b (needed when interpolating compute x back up */
696: VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);
697:
698:
699: while (cycles--) {
700: DMMGSolveFASMCycle(dmmg,level-1,converged);
701: }
702: }
703: else { /* for the coarsest level */
704: for (k=0; k<dmmg[0]->coarsesmooth; k++) {
705: NLFRelax_DAAD(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);
706: }
707: if (dmmg[0]->monitorall){
708: DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);
709: VecAXPY(dmmg[0]->w,-1.0,dmmg[0]->r);
710: VecNorm(dmmg[0]->w,NORM_2,&norm);
711: for (k=0; k<level+1; k++) {PetscPrintf(dmmg[0]->comm," ");}
712: PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %g\n",norm);
713: }
714: if (j == nlevels || dmmg[j]->monitorall) {
715: /* norm( residual_fine - f(x_fine) ) */
716: VecNorm(dmmg[j]->w,NORM_2,&norm);
717:
718: if (j == nlevels) {
719: if (norm < dmmg[level]->abstol) {
720: *converged = PETSC_TRUE;
721: goto theend;
722: }
723:
724: if (norm < dmmg[level]->rrtol){
725: *converged = PETSC_TRUE;
726: goto theend;
727:
728: }
729: }
730: }
732:
733: }
734: j=level;
735: if(j) { /* not for the coarsest level */
736: /* x_fine = x_fine + R'*(x_coarse - RI*x_fine) */
737: VecAXPY(dmmg[j-1]->x,-1.0,dmmg[j-1]->b);
738: MatInterpolateAdd(mg[j]->interpolate,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);
740: if (dmmg[j]->monitorall) {
741: /* norm( F(x_fine) - residual_fine ) */
742: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
743: VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
744: VecNorm(dmmg[j]->w,NORM_2,&norm);
745: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
746: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
747: }
749: /* Relax residual_fine - F(x_fine) = 0 */
750: for (k=0; k<dmmg[j]->postsmooth; k++) {
751: NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
752: }
754: if (dmmg[j]->monitorall) {
755: /* norm( F(x_fine) - residual_fine ) */
756: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
757: VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
758: VecNorm(dmmg[j]->w,NORM_2,&norm);
759: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
760: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
761: }
762:
763:
764: }
765: theend:
766: return(0);
767: }
769: /*
770: This is function to implement multiplicative FAS with block smoother
773: Options:
774:
775: -dmmg_fas_cycles 1 : FAS v-cycle
776: 2 : FAS w-cycle
779: */
784: PetscErrorCode DMMGSolveFASMCycle9(DMMG *dmmg,PetscInt level,PetscBool * converged)
785: {
787: PetscInt j,k,cycles=1,nlevels=level;//nlevels=dmmg[0]->nlevels-1;
788: // I need to put nlevels=level in order to get grid sequence correctly
789: PetscReal norm;
790: PC_MG **mg;
791: PC pc;
792:
794:
795:
796: PetscOptionsGetInt(PETSC_NULL,"-dmmg_fas_cycles",&cycles,PETSC_NULL);
797:
798: KSPGetPC(dmmg[level]->ksp,&pc);
799: mg = ((PC_MG**)pc->data);
800: // for (j=level; j>0; j--) {
801: j=level;
802: //PetscPrintf(dmmg[level]->comm,"j=%d,nlevels=%d",j,nlevels);
803: if(j) {/* not the coarsest level */
804: /* Relax residual_fine - F(x_fine) = 0 */
805: for (k=0; k<dmmg[j]->presmooth; k++) {
806: NLFRelax_DAAD9(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
807: }
808:
809:
811: /* R*(residual_fine - F(x_fine)) */
812: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
813: VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);
815: if (j == nlevels || dmmg[j]->monitorall) {
816: /* norm( residual_fine - f(x_fine) ) */
817: VecNorm(dmmg[j]->w,NORM_2,&norm);
818:
819: if (j == nlevels) {
820: if (norm < dmmg[level]->abstol) {
821: *converged = PETSC_TRUE;
822: goto theend;
823: }
824: /* if (i == 0) {
825: dmmg[level]->rrtol = norm*dmmg[level]->rtol;
826: } else {*/
827: if (norm < dmmg[level]->rrtol){
828: *converged = PETSC_TRUE;
829: goto theend;
830:
831: }
832: }
833: }
835: if (dmmg[j]->monitorall) {
836: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
837: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
838: }
839: MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);
840:
841: /* F(RI*x_fine) */
842: VecScatterBegin(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
843: VecScatterEnd(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
844: DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);
846: /* residual_coarse = F(RI*x_fine) + R*(residual_fine - F(x_fine)) */
847: VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);
849: /* save RI*x_fine into b (needed when interpolating compute x back up */
850: VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);
851:
852:
853: while (cycles--) {
854: DMMGSolveFASMCycle9(dmmg,level-1,converged);
855: }
856: }
857: else { /* for the coarsest level */
858: for (k=0; k<dmmg[0]->coarsesmooth; k++) {
859: NLFRelax_DAAD9(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);
860: }
861: if (dmmg[0]->monitorall){
862: DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);
863: VecAXPY(dmmg[0]->w,-1.0,dmmg[0]->r);
864: VecNorm(dmmg[0]->w,NORM_2,&norm);
865: for (k=0; k<level+1; k++) {PetscPrintf(dmmg[0]->comm," ");}
866: PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %g\n",norm);
867: }
868: if (j == nlevels || dmmg[j]->monitorall) {
869: /* norm( residual_fine - f(x_fine) ) */
870: VecNorm(dmmg[j]->w,NORM_2,&norm);
871:
872: if (j == nlevels) {
873: if (norm < dmmg[level]->abstol) {
874: *converged = PETSC_TRUE;
875: goto theend;
876: }
877: /* if (i == 0) {
878: dmmg[level]->rrtol = norm*dmmg[level]->rtol;
879: } else {*/
880: if (norm < dmmg[level]->rrtol){
881: *converged = PETSC_TRUE;
882: goto theend;
883:
884: }
885: }
886: }
888:
889: }
890: j=level;
891: if(j) { /* not for the coarsest level */
892: /* x_fine = x_fine + R'*(x_coarse - RI*x_fine) */
893: VecAXPY(dmmg[j-1]->x,-1.0,dmmg[j-1]->b);
894: MatInterpolateAdd(mg[j]->interpolate,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);
896: if (dmmg[j]->monitorall) {
897: /* norm( F(x_fine) - residual_fine ) */
898: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
899: VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
900: VecNorm(dmmg[j]->w,NORM_2,&norm);
901: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
902: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
903: }
905: /* Relax residual_fine - F(x_fine) = 0 */
906: for (k=0; k<dmmg[j]->postsmooth; k++) {
907: NLFRelax_DAAD9(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
908: }
910: if (dmmg[j]->monitorall) {
911: /* norm( F(x_fine) - residual_fine ) */
912: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
913: VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
914: VecNorm(dmmg[j]->w,NORM_2,&norm);
915: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
916: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
917: }
918:
919: /* if(j==nlevels){
920: if (dmmg[level]->monitor){
921: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
922: VecNorm(dmmg[level]->w,NORM_2,&norm);
923:
924: PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);
925:
926: }
927: }*/
928: }
929: theend:
930: return(0);
931: }
933: /*
934: This is function to implement full FAS with block smoother(9 points together)
937: Options:
938:
939: -dmmg_fas_cycles 1 : FAS v-cycle
940: 2 : FAS w-cycle
943: */
948: PetscErrorCode DMMGSolveFASFCycle(DMMG *dmmg,PetscInt l,PetscBool * converged)
949: {
951: PetscInt j;//l = dmmg[0]->nlevels-1;
952: PC_MG **mg;
953: PC pc;
956: KSPGetPC(dmmg[l]->ksp,&pc);
957: mg = ((PC_MG**)pc->data);
958: // restriction all the way down to the coarse level
959: if(l>0) {
960: for(j=l;j>0;j--) {
961:
962: /* R*(residual_fine - F(x_fine)) */
963: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
964: VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);
966: MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);
967:
968: /* F(RI*x_fine) */
969: VecScatterBegin(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
970: VecScatterEnd(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
971: DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);
973: /* residual_coarse = F(RI*x_fine) + R*(residual_fine - F(x_fine)) */
974: VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);
976: /* save RI*x_fine into b (needed when interpolating compute x back up */
977: VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);
979: }
981: // all the way up to the finest level
982: for (j=0; j<l; j++) {
983: DMMGSolveFASMCycle(dmmg,j,PETSC_NULL);
984: /* x_fine = x_fine + R'*(x_coarse - RI*x_fine) */
985: VecAXPY(dmmg[j]->x,-1.0,dmmg[j]->b);
986: MatInterpolateAdd(mg[j+1]->interpolate,dmmg[j]->x,dmmg[j+1]->x,dmmg[j+1]->x);
988: }
989: }
990: DMMGSolveFASMCycle(dmmg,l,converged);
991: return(0);
992: }
996: /*
997: This is function to implement full FAS with block smoother ( 9 points together)
1000: Options:
1001:
1002: -dmmg_fas_cycles 1 : FAS v-cycle
1003: 2 : FAS w-cycle
1006: */
1009: PetscErrorCode DMMGSolveFASFCycle9(DMMG *dmmg,PetscInt l,PetscBool * converged)
1010: {
1012: PetscInt j;//l = dmmg[0]->nlevels-1;
1013: PC_MG **mg;
1014: PC pc;
1017: KSPGetPC(dmmg[l]->ksp,&pc);
1018: mg = ((PC_MG**)pc->data);
1019: // restriction all the way down to the coarse level
1020: if(l>0) {
1021: for(j=l;j>0;j--) {
1022:
1023: /* R*(residual_fine - F(x_fine)) */
1024: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
1025: VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);
1027: MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);
1028:
1029: /* F(RI*x_fine) */
1030: VecScatterBegin(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
1031: VecScatterEnd(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD);
1032: DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);
1034: /* residual_coarse = F(RI*x_fine) + R*(residual_fine - F(x_fine)) */
1035: VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);
1037: /* save RI*x_fine into b (needed when interpolating compute x back up */
1038: VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);
1040: }
1042: // all the way up to the finest level
1043: for (j=0; j<l; j++) {
1044: DMMGSolveFASMCycle9(dmmg,j,PETSC_NULL);
1045: /* x_fine = x_fine + R'*(x_coarse - RI*x_fine) */
1046: VecAXPY(dmmg[j]->x,-1.0,dmmg[j]->b);
1047: MatInterpolateAdd(mg[j+1]->interpolate,dmmg[j]->x,dmmg[j+1]->x,dmmg[j+1]->x);
1049: }
1050: }
1051: DMMGSolveFASMCycle9(dmmg,l,converged);
1052: return(0);
1053: }
1055: /*
1056: This is function is to solve nonlinear system with FAS
1058: Options:
1060: -dmmg_fas_9: using block smoother
1061:
1062: -dmmg_fas_full: using full FAS
1065: */
1068: PetscErrorCode DMMGSolveFASCycle(DMMG *dmmg,PetscInt level)
1069: {
1071: PetscInt i;
1072: PetscBool converged = PETSC_FALSE, flg = PETSC_FALSE,flgb = PETSC_FALSE;
1073: PetscReal norm;
1076: DMMGSolveFASSetUp(dmmg,level);
1077: PetscOptionsGetBool(PETSC_NULL,"-dmmg_fas_9",&flgb,PETSC_NULL);
1078:
1079: if(flgb){
1081: PetscOptionsGetBool(PETSC_NULL,"-dmmg_fas_full",&flg,PETSC_NULL);
1082: if (flg) {
1083: for(i=0;i<1000;i++){
1084: PetscPrintf(dmmg[level]->comm,"%D ",i+1);
1085: DMMGSolveFASFCycle9(dmmg,level,&converged);
1086: if (dmmg[level]->monitor){
1087: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
1088: VecNorm(dmmg[level]->w,NORM_2,&norm);
1089: PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);
1090: }
1091: if (converged) return(0);
1092: }
1093: }
1094: else{
1095: for(i=0;i<1000;i++){
1096: PetscPrintf(dmmg[level]->comm,"%D ",i+1);
1097: DMMGSolveFASMCycle9(dmmg,level,&converged);
1098: if (dmmg[level]->monitor){
1099: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
1100: VecNorm(dmmg[level]->w,NORM_2,&norm);
1101: PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);
1102: }
1104: if (converged) return(0);
1105: }
1106: }
1107: }
1108: else {
1109: flg = PETSC_FALSE;
1110: PetscOptionsGetBool(PETSC_NULL,"-dmmg_fas_full",&flg,PETSC_NULL);
1111: if (flg) {
1112: for(i=0;i<1000;i++){
1113: PetscPrintf(dmmg[level]->comm,"%D ",i+1);
1114: DMMGSolveFASFCycle(dmmg,level,&converged);
1115: if (dmmg[level]->monitor){
1116: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
1117: VecNorm(dmmg[level]->w,NORM_2,&norm);
1118: PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);
1119: }
1120: if (converged) return(0);
1121: }
1122: }
1123: else{
1124: for(i=0;i<1000;i++){
1125: PetscPrintf(dmmg[level]->comm,"%D ",i+1);
1126: DMMGSolveFASMCycle(dmmg,level,&converged);
1127: if (dmmg[level]->monitor){
1128: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
1129: VecNorm(dmmg[level]->w,NORM_2,&norm);
1130: PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);
1131: }
1133: if (converged) return(0);
1134: }
1136: }
1137: }
1138: return(0);
1139: }
1141: /*
1142: This is function is to implement one FAS iteration
1144: Options:
1146: -dmmg_fas_9: using block smoother
1147:
1148: -dmmg_fas_full: using full FAS
1150: */
1153: PetscErrorCode DMMGSolveFASCyclen(DMMG *dmmg,PetscInt level)
1154: {
1156: PetscBool converged = PETSC_FALSE, flg = PETSC_FALSE,flgb = PETSC_FALSE;
1157: PetscReal norm;
1158: // PC_MG **mg;
1159: //PC pc;
1162: DMMGSolveFASSetUp(dmmg,level);
1163: PetscOptionsGetBool(PETSC_NULL,"-dmmg_fas_9",&flgb,PETSC_NULL);
1164:
1165: if(flgb){
1167: PetscOptionsGetBool(PETSC_NULL,"-dmmg_fas_full",&flg,PETSC_NULL);
1168: if (flg) {
1169:
1170: DMMGSolveFASFCycle9(dmmg,level,&converged);
1171: if (dmmg[level]->monitor){
1172: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
1173: VecNorm(dmmg[level]->w,NORM_2,&norm);
1174: PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);
1175: }
1177: }
1178: else{
1179:
1180: DMMGSolveFASMCycle9(dmmg,level,&converged);
1181: if (dmmg[level]->monitor){
1182: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
1183: VecNorm(dmmg[level]->w,NORM_2,&norm);
1184: PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);
1185: }
1188: }
1189: }
1190: else {
1191: flg = PETSC_FALSE;
1192: PetscOptionsGetBool(PETSC_NULL,"-dmmg_fas_full",&flg,PETSC_NULL);
1193: if (flg) {
1194:
1195: DMMGSolveFASFCycle(dmmg,level,&converged);
1196: if (dmmg[level]->monitor){
1197: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
1198: VecNorm(dmmg[level]->w,NORM_2,&norm);
1199: PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);
1200: }
1202: }
1203: else{
1204:
1205: DMMGSolveFASMCycle(dmmg,level,&converged);
1206: if (dmmg[level]->monitor){
1207: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
1208: VecNorm(dmmg[level]->w,NORM_2,&norm);
1209: PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);
1210: }
1214: }
1215: }
1216:
1218: return(0);
1219: }
1222: /*
1224: This is function to implement Nonlinear CG to accelerate FAS
1226: In order to use this acceleration, the option is
1228: -dmmg_fas_NCG
1232: */
1237: PetscErrorCode DMMGSolveFAS_NCG(DMMG *dmmg, PetscInt level)
1238: {
1239: SNES snes = dmmg[level]->snes;
1240: SNES_LS *neP = (SNES_LS*)snes->data;
1242: PetscInt maxits,i,lits;
1243: PetscBool lssucceed;
1244: // MatStructure flg = DIFFERENT_NONZERO_PATTERN;
1245: PetscReal fnorm,gnorm,xnorm,ynorm,betaFR,betaPR,beta,betaHS,betaDY;
1246: Vec Y,X,F,G,W,Gradold,Sk;
1247: //KSP ksp;
1251: VecDuplicate(dmmg[level]->x,&Sk);
1252: PetscObjectReference((PetscObject)dmmg[level]->x);
1253: if (snes->vec_sol) { VecDestroy(snes->vec_sol); }
1254: snes->vec_sol = dmmg[level]->x;
1255: if (!snes->setupcalled) { SNESSetUp(snes); }
1256: if (snes->conv_hist_reset) snes->conv_hist_len = 0;
1257: PetscLogEventBegin(SNES_Solve,snes,0,0,0);
1258: snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
1261: snes->reason = SNES_CONVERGED_ITERATING;
1263: maxits = snes->max_its; /* maximum number of iterations */
1264: X = snes->vec_sol; /* solution vector */
1265: F = snes->vec_func; /* residual vector */
1266: Y = snes->work[0]; /* work vectors */
1267: G = snes->work[1];
1268: W = snes->work[2];
1270: PetscObjectTakeAccess(snes);
1271: snes->iter = 0;
1272: PetscObjectGrantAccess(snes);
1273:
1274: X = dmmg[level]->x;
1275: VecCopy(X,Y);
1276: VecCopy(X,G);
1278: // to get the residual for the F
1279: SNESComputeFunction(snes,X,F);
1280:
1281: VecNorm(F,NORM_2,&fnorm); /* fnorm <- ||F|| */
1282: if (fnorm != fnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1283: PetscObjectTakeAccess(snes);
1284: snes->norm = fnorm;
1285: PetscObjectGrantAccess(snes);
1286: SNESLogConvHistory(snes,fnorm,0);
1287: SNESMonitor(snes,0,fnorm);
1289: if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; return(0);}
1291: /* set parameter for default relative tolerance convergence test */
1292: snes->ttol = fnorm*snes->rtol;
1293: // dmmg[level]->rrtol= snes->ttol;
1295: // set this to store the old grad
1296: Gradold=snes->vec_sol_update;
1297:
1298: // compute the search direction Y
1299: // I need to put Q(x)=x-FAS(x) here
1300: DMMGSolveFASCyclen(dmmg,level);
1301: // F = X - dmmg[level]->x; this is the gradient direction
1302: VecAXPY(Y,-1.0,X);
1303: // copy the gradient to the old
1304: VecCopy(Y,Gradold);
1305: // copy X back
1306: VecCopy(G,X);
1307: VecWAXPY(Sk,-1.0,X,X);
1308: // so far I put X=X_c, F= F(x_c), Gradold= Y=grad(x_c)
1310: // for (i=0; i<maxits; i++) {
1312: for (i=0; i<10000; i++) {
1315: PetscPrintf(PETSC_COMM_WORLD,"iter=%d",i+1);
1316:
1317:
1318: // X=x_c, F=F(x_c),Y search direction; G=F(x_new), W=x_new,
1319: VecNorm(X,NORM_2,&xnorm);
1320: (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);
1321: PetscInfo4(snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);
1322: if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
1323: PetscPrintf(PETSC_COMM_WORLD,"step=%g,oldnorm=%g,norm=%g ",ynorm,fnorm,gnorm);
1324:
1325: fnorm=gnorm; //copy the new function norm; this will effect the line_search
1326: VecWAXPY(Sk,-1.0,X,W);
1327: //update the new solution
1328: ierr=VecCopy(W,X);
1329: ierr=VecCopy(G,F);
1331:
1332: // compute the new search direction G
1333: // I need to put Q(x)=x-FAS(x) here
1335: DMMGSolveFASCyclen(dmmg,level);
1336: //G = X - dmmg[level]->x; G is the new gradient, Y is old gradient
1337:
1338: VecWAXPY(G,-1.0,X,W);
1339: // copy W back to X
1340: VecCopy(W,X);
1341: VecNorm(G,NORM_2,&gnorm);
1342: VecNorm(Gradold,NORM_2,&ynorm);
1343: betaFR = gnorm*gnorm/(ynorm*ynorm); //FR_beta
1344:
1345: VecWAXPY(W,-1.0,Gradold,G);
1346: VecDot(W,G,&gnorm);
1347: // VecNorm(G,NORM_2,&gnorm);
1348: VecNorm(Gradold,NORM_2,&ynorm);
1349: betaPR = gnorm/(ynorm*ynorm); //PR_beta
1350:
1351: if ( betaPR<-betaFR)
1352: {
1353: beta =- betaFR;
1354: }
1355: else {
1356: if (betaPR>betaFR)
1357: {beta=betaFR;}
1358: else{
1360: beta=betaPR;
1361: }
1362: }
1363: // beta=betaFR;
1364: //beta=betaPR;
1366: // try another beta
1367:
1368: VecWAXPY(W,-1.0,Gradold,G);
1369: VecDot(W,G,&betaHS);
1370: VecDot(W,Y,&gnorm);
1371: betaHS=-betaHS/gnorm;
1372: VecDot(G,G,&betaDY);
1373: betaDY=-betaDY/gnorm;
1374: if(betaHS<betaDY)
1375: beta=betaHS;
1376: else
1377: beta=betaDY;
1378: if(beta<0)
1379: beta=0;
1380:
1381: PetscPrintf(PETSC_COMM_WORLD,"betaHS=%g,betaDY=%g\n",betaHS,betaDY);
1382:
1384:
1385: // compute the c_2
1386: VecDot(G,Y,&gnorm);
1387: VecDot(Gradold,Y,&ynorm);
1388: PetscPrintf(PETSC_COMM_WORLD,"beta=%g,c_2=%g\n",beta,fabs(gnorm/ynorm));
1389: VecNorm(G,NORM_2,&gnorm);
1390: VecNorm(Y,NORM_2,&ynorm);
1391: PetscPrintf(PETSC_COMM_WORLD,"size=%g\n",fabs(gnorm/(beta*ynorm)));
1393: // update the direction: Y= G + beta * Y
1394: VecAXPBY(Y,1.0,beta,G);
1395: VecCopy(G,Gradold);
1396: //ierr =VecCopy(G,Y);
1397: snes->iter = i+1;
1398: snes->norm = fnorm;
1399: PetscObjectGrantAccess(snes);
1400: SNESLogConvHistory(snes,fnorm,lits);
1401: SNESMonitor(snes,i+1,fnorm);
1402:
1403: if (!lssucceed) {
1404: PetscBool ismin;
1405: beta=0;
1406: if (++snes->numFailures >= snes->maxFailures) {
1407: snes->reason = SNES_DIVERGED_LINE_SEARCH;
1408: SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);
1409: if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
1410: break;
1411: }
1412: }
1413:
1415: /* Test for convergence */
1416: if (snes->ops->converged) {
1417: VecNorm(X,NORM_2,&xnorm); /* xnorm = || X || */
1418: (*snes->ops->converged)(snes,snes->iter,xnorm,1.0,fnorm,&snes->reason,snes->cnvP);
1419: if (snes->reason) {
1420: break;
1421: }
1422: }
1423: }
1424: if (X != snes->vec_sol) {
1425: VecCopy(X,snes->vec_sol);
1426: }
1427: if (F != snes->vec_func) {
1428: VecCopy(F,snes->vec_func);
1429: }
1430: if (i == maxits) {
1431: PetscInfo1(snes,"SNESSolve_LS: Maximum number of iterations has been reached: %D\n",maxits);
1432: snes->reason = SNES_DIVERGED_MAX_IT;
1433: }
1434: PetscPrintf(PETSC_COMM_WORLD,"reason=%d\n",snes->reason);
1435: // VecView(X,PETSC_VIEWER_STDOUT_WORLD);
1436: return(0);
1437: }
1441: /*
1443: This is function to implement NonGMRES to accelerate FAS
1445: In order to use this acceleration, the option is
1447: -dmmg_fas_NGMRES
1450: Options:
1452: -dmmg_fas_ngmres_m : the number of previous solutions to keep
1454: */
1458: PetscErrorCode DMMGSolveFAS_NGMRES(DMMG *dmmg, PetscInt level)
1459: {
1460: SNES snes = dmmg[level]->snes;
1462: PetscInt maxits=10000,i,k,l,j,subm=3,iter;
1463: PetscOptionsGetInt(PETSC_NULL,"-dmmg_fas_ngmres_m",&subm,PETSC_NULL);
1465: PetscBool restart=PETSC_FALSE, selectA=PETSC_FALSE;
1466: PetscReal fnorm,gnorm,dnorm,dnormtemp,dminnorm,fminnorm,tol=1.e-12,gammaA=2,epsilonB=0.1,deltaB=0.9,gammaC;
1467: Vec X,F,G,W,D,u[subm],res[subm];
1468: PetscScalar H[subm][subm],q[subm][subm],beta[subm],xi[subm],alpha[subm],alphasum,det,Hinv[16];
1472: gammaC=2; if (gammaA>gammaC) gammaC=gammaA;
1473: VecDuplicate(dmmg[level]->x,&X);
1474: VecDuplicate(dmmg[level]->x,&F);
1475: VecDuplicate(dmmg[level]->x,&W);
1476: VecDuplicate(dmmg[level]->x,&G);
1477: VecDuplicate(dmmg[level]->x,&D);
1479: for(i=0;i<subm;i++) {/* get the space for the solution */
1480: VecDuplicate(dmmg[level]->x,&u[i]);
1481: VecDuplicate(dmmg[level]->x,&res[i]);
1482: }
1484: X = dmmg[level]->x;
1485: VecCopy(X,u[0]);
1486:
1487: // to get the residual for the F
1488: SNESComputeFunction(snes,X,F);
1489: VecCopy(F,res[0]);
1490: VecNorm(F,NORM_2,&fnorm); /* fnorm <- ||F|| */
1491: fnorm=fnorm*fnorm;
1492: iter=1;
1493: restartline:
1494:
1495: q[0][0] = fnorm; fminnorm=fnorm;
1498: for (k=1; k<maxits; k++) {
1500:
1501: PetscPrintf(PETSC_COMM_WORLD,"\n k=%d,iter=%d fmin=%g ",k,iter++,sqrt(fminnorm));
1502:
1503: /* compute the X=u^M , F=r, fnorm =||F||*/
1505: DMMGSolveFASCyclen(dmmg,level);
1506: SNESComputeFunction(snes,X,F);
1507: VecNorm(F,NORM_2,&fnorm);
1508: if (fnorm < tol) { return(0);}
1509: fnorm =fnorm*fnorm;
1510: if (fnorm<fminnorm) fminnorm=fnorm;
1511: //PetscPrintf(PETSC_COMM_WORLD,"fmin=%g",sqrt(fminnorm));
1512: /* compute u^A */
1513: //l=min(subm,k)
1514: l=subm;
1515: if (k<l) l=k;
1516:
1517: /* compute the matrix H and RHS */
1518:
1519:
1520: for(i=0;i<l;i++){
1521: VecDot(F,res[i],&xi[i]);
1522: beta[i]=fnorm-xi[i];
1523: }
1524:
1525: for(i=0;i<l;i++){
1526: for(j=0;j<l;j++){
1527: H[i][j]=q[i][j]-xi[i]-xi[j]+fnorm;
1528: }
1529: }
1530: /* Here is special for subm=2 */
1531: if(l==1){
1532: // H[0][0] = q[0][0]-xi[0]-xi[0]+fnorm;
1533: alpha[0]= beta[0]/H[0][0];
1534: }
1535: if(l==2){
1536:
1537: alpha[0]= ( beta[0]*H[1][1] -beta[1]*H[0][1] )/( H[0][0]*H[1][1]- H[0][1]*H[1][0]);
1538: alpha[1]= ( beta[1]*H[0][0] -beta[0]*H[1][0] )/( H[0][0]*H[1][1]- H[0][1]*H[1][0]);
1539:
1540: }
1541: if(l==3) {
1542:
1543: det = H[0][0]*(H[1][1]*H[2][2]-H[1][2]*H[2][1])-H[0][1]*(H[1][0]*H[2][2]-H[1][2]*H[2][0])+H[0][2]*(H[1][0]*H[2][1]-H[1][1]*H[2][0]);
1544: alpha[0]= (beta[0]*(H[1][1]*H[2][2]-H[1][2]*H[2][1])-H[0][1]*(beta[1]*H[2][2]-H[1][2]*beta[2])+H[0][2]*(beta[1]*H[2][1]-H[1][1]*beta[2]))/det;
1545: alpha[1]=(H[0][0]*(beta[1]*H[2][2]-H[1][2]*beta[2])-beta[0]*(H[1][0]*H[2][2]-H[1][2]*H[2][0])+H[0][2]*(H[1][0]*beta[2]-beta[1]*H[2][0]))/det;
1546: alpha[2]=(H[0][0]*(H[1][1]*beta[2]-beta[1]*H[2][1])-H[0][1]*(H[1][0]*beta[2]-beta[1]*H[2][0])+beta[0]*(H[1][0]*H[2][1]-H[1][1]*H[2][0]))/det;
1549: }
1550:
1551: if(l==4){
1552: Hinv[0]=H[0][0];Hinv[1]=H[0][1];Hinv[2]=H[0][2];Hinv[3]=H[0][3];
1553: Hinv[4]=H[1][0];Hinv[5]=H[1][1];Hinv[6]=H[1][2];Hinv[7]=H[1][3];
1554: Hinv[8]=H[2][0];Hinv[9]=H[2][1];Hinv[10]=H[2][2];Hinv[11]=H[2][3];
1555: Hinv[12]=H[3][0];Hinv[13]=H[3][1];Hinv[14]=H[3][2];Hinv[15]=H[3][3];
1556: Kernel_A_gets_inverse_A_4(Hinv,0.0);
1557: for(i=0;i<l;i++)
1558: alpha[i]=Hinv[4*i]*beta[0]+Hinv[4*i+1]*beta[1]+Hinv[4*i+2]*beta[2]+Hinv[4*i+3]*beta[3];
1559:
1560: }
1561: alphasum=0;
1562: for (i=0;i<l;i++)
1563: alphasum=alphasum+alpha[i];
1564:
1565: /* W= u^A */
1566: VecCopy(X,W);
1567: VecAXPBY(W,0.0,1-alphasum,X);
1568:
1569: for(i=0;i<l;i++)
1570: VecAXPY(W,alpha[i],u[i]);
1571:
1572: /* W= F(G) */
1573: SNESComputeFunction(snes,W,G);
1574: VecNorm(G,NORM_2,&gnorm);
1575: gnorm=gnorm*gnorm;
1578:
1579: /* select the uA or uM */
1580: // Criterion A
1581: if(sqrt(gnorm)<gammaA*sqrt(fminnorm)){
1582: //PetscPrintf(PETSC_COMM_WORLD,"Crite A\n");
1583: selectA=PETSC_TRUE;
1584: }
1585: // Criterion B
1586:
1587: ierr=VecCopy(W,D);
1588: ierr=VecAXPY(D,-1,X);
1589: ierr=VecNorm(D,NORM_2,&dnorm);
1590: dminnorm=10000000;
1591: for(i=0;i<l;i++) {
1592: ierr=VecCopy(W,D);
1593: ierr=VecAXPY(D,-1,u[i]);
1594: ierr=VecNorm(D,NORM_2,&dnormtemp);
1595: if(dnormtemp<dminnorm) dminnorm=dnormtemp;
1596: }
1597: if( epsilonB*dnorm<dminnorm || sqrt(gnorm)<deltaB*sqrt(fminnorm))
1598: selectA =PETSC_TRUE;
1599: else
1600: selectA=PETSC_FALSE;
1601:
1602: if(selectA){ /* uA selected */
1603: selectA=PETSC_FALSE;
1604: // PetscPrintf(PETSC_COMM_WORLD,"Select A\n");
1605: VecCopy(W,X);
1606: VecCopy(G,F);
1607: fnorm=gnorm;
1608: }
1609: if (fnorm<fminnorm) fminnorm=fnorm;
1610:
1612: /* Test for convergence */
1613: if (sqrt(fnorm)<tol) {
1614: return(0);
1615: }
1616:
1617: /* Test for restart */
1618: if(sqrt(gnorm)>gammaC*sqrt(fminnorm)) {
1619: PetscPrintf(PETSC_COMM_WORLD,"restart for C ");
1620: restart=PETSC_TRUE;
1621: }
1622: if(epsilonB*dnorm>dminnorm && sqrt(gnorm)>deltaB*sqrt(fminnorm)) {
1623: PetscPrintf(PETSC_COMM_WORLD,"restart for D ");
1624: restart=PETSC_TRUE;
1625: }
1626: /* Prepation for the next iteration */
1627:
1628: //turn off restart
1629: //restart=PETSC_FALSE;
1630: if(restart){
1631: restart=PETSC_FALSE;
1632: goto restartline;
1633: }
1634: else {
1635: j=k%subm;
1636: VecCopy(F,res[j]);
1637: VecCopy(X,u[j]);
1638: for(i=0;i<l;i++){
1639: ierr= VecDot(F,res[i],&q[j][i]);
1640: q[i][j]=q[j][i];
1641: }
1642: if(l<subm)
1643: q[j][j]=fnorm;
1644: }
1645: }
1647: return(0);
1648: }
1650: /*
1651: This is a function to provide the function value:
1653: x-FAS(x)
1655: */
1661: PetscErrorCode DMMGFASFunction(SNES snes,Vec x,Vec f,void *ptr)
1662: {
1663: DMMG* dmmg=(DMMG*)ptr;
1665: Vec temp;
1666: PetscInt level=dmmg[0]->nlevels-1;
1668:
1669: VecDuplicate(dmmg[level]->x,&temp);
1670:
1671:
1672: VecCopy(dmmg[level]->x,temp);
1673: VecCopy(x,dmmg[level]->x);
1674: VecCopy(x,f);
1675:
1676: // I need to put -F(x)=x-FAS(x) here
1677: DMMGSolveFASCyclen(dmmg,level);
1678: VecAXPY(f,-1.0,dmmg[level]->x);
1679: // y = alpha x + y.
1680: ierr=VecCopy(temp,dmmg[level]->x);
1681: //copy W back to X
1682:
1683: return(0);
1684: }
1688: /* this function is to implement Quasi-Newton method with implicit Broyden updating methods(limit memory version)
1690: In order to use this method, the option is -dmmg_fas_QNewton
1692: Options:
1694: -dmmg_fas_QNewton_m: the number of the vectors to keep for inverse of Jacobian
1696: -dmmg_fas_initialJacobian: will use matrix-free GMRES to solve the initial Jacobian
1698: with options -snes_mf -snes_max_it 1 -ksp_max_it n
1701: In this function, does not have line search and nonlinear gmres acceleration
1703: */
1707: PetscErrorCode DMMGSolveFAS_QNewton(DMMG *dmmg, PetscInt level)
1708: {
1710:
1712: SNES snes = dmmg[level]->snes, snes0;
1714: PetscInt maxits=10000,i,k,l,subm=3,subm01;
1715: PetscOptionsGetInt(PETSC_NULL,"-dmmg_fas_QNewton_m",&subm,PETSC_NULL);
1716: subm01=subm-1;
1717: PetscBool flg = PETSC_FALSE;
1718: PetscReal fnorm,gnorm,tol=1.e-12;
1719: Vec X,F,G,W,D,Y,v[subm],w[subm],s0,s1,F0,F1;
1720:
1724: PetscOptionsGetBool(PETSC_NULL,"-dmmg_fas_initialJacobian",&flg,PETSC_NULL);
1725: VecDuplicate(dmmg[level]->x,&X);
1726: VecDuplicate(dmmg[level]->x,&F);
1727: VecDuplicate(dmmg[level]->x,&W);
1728: VecDuplicate(dmmg[level]->x,&G);
1729: VecDuplicate(dmmg[level]->x,&D);
1730: VecDuplicate(dmmg[level]->x,&Y);
1731: VecDuplicate(dmmg[level]->x,&s0);
1732: VecDuplicate(dmmg[level]->x,&s1);
1733: VecDuplicate(dmmg[level]->x,&F0);
1734: VecDuplicate(dmmg[level]->x,&F1);
1735:
1736: // creat a snes for solve the initial Jacobian
1737: SNESCreate(dmmg[level]->comm,&snes0);
1738: SNESSetFunction(snes0,F,DMMGFASFunction,dmmg);
1739: SNESSetFromOptions(snes0);
1741: for(i=0;i<subm;i++) {/* get the space for the solution */
1742: VecDuplicate(dmmg[level]->x,&v[i]);
1743: VecDuplicate(dmmg[level]->x,&w[i]);
1744: }
1746: //We first try B0==I
1747: X = dmmg[level]->x;
1748:
1749: if(flg){
1750: ierr= VecAXPBY(Y,0.0,0.0,X);
1751: ierr= VecCopy(X,s0);
1752: ierr= SNESSolve(snes0,Y,s0);
1753: ierr= VecAXPY(s0,-1.0,X);
1754: }
1755: else{
1756: ierr=VecCopy(X,W);
1757: ierr=VecCopy(X,Y);
1758:
1759: // I need to put -F(x)=x-FAS(x) here
1761: DMMGSolveFASCyclen(dmmg,level);
1762: VecAXPY(Y,-1.0,X);
1763: // y = alpha x + y.
1764: ierr=VecCopy(W,X);
1765: //copy W back to X
1766:
1767: // Y stores the -F(x)
1768: ierr= VecAXPBY(Y,0.0,-1.0,X);
1769: ierr= VecCopy(Y,s0);
1771: }
1772:
1773: VecAXPY(X,1.0,s0);
1774:
1775: for(k=0; k<maxits; k++){
1776:
1777: /* Test for convergence */
1778: SNESComputeFunction(snes,X,F);
1779: VecNorm(F,NORM_2,&fnorm); /* fnorm <- ||F|| */
1780: PetscPrintf(PETSC_COMM_WORLD,"k=%d, fnorm=%g\n",
1781: k,fnorm);
1782: if (fnorm<tol) {
1783: return(0);
1784: }
1785:
1787: if(flg){
1788: // ierr= SNESSolve(snes0,Y,s1);
1789: ierr= VecAXPBY(Y,0.0,0.0,X);
1790: ierr= VecCopy(dmmg[level]->x,s1);
1791: ierr= SNESSolve(snes0,Y,s1);
1792: ierr= VecAXPY(s1,-1.0,X);
1793:
1794: }
1795: else{
1796: ierr=VecCopy(X,W);
1797: ierr=VecCopy(X,Y);
1798:
1799: // I need to put -F(x)=x-FAS(x) here
1801: DMMGSolveFASCyclen(dmmg,level);
1802: VecAXPY(Y,-1.0,X);
1803: // y = alpha x + y.
1804: ierr=VecCopy(W,X);
1805: //copy W back to X
1806:
1807: //So far, I got X=x_k, Y=-F(x_k)
1808: // I should solve the G=-B_0^{-1}F(x_k) first, but I choose B_0^{-1}=I,
1809: ierr=VecCopy(Y,F1);
1810: ierr= VecAXPBY(Y,0.0,-1.0,X);
1811: ierr=VecCopy(Y,s1);
1813: }
1814:
1815: l=subm;
1816: if (k<l) l=k;
1817:
1818:
1819: for (i=0;i<l;i++){
1820: // compute [I+v(i)w(i)^T]*s(k)
1821: ierr= VecDot(w[i],s1,&gnorm);
1822: VecAXPY(s1,gnorm,v[i]);
1823:
1824: }
1825: if(l==subm) {
1826: for(i=0;i<subm01;i++){
1827: ierr= VecCopy(w[i+1],w[i]);
1828: ierr= VecCopy(v[i+1],v[i]);
1829: }
1830: l--;
1831: }
1832: ierr= VecCopy(s0,w[l]);
1833: ierr= VecCopy(s0,Y);
1834: ierr= VecCopy(s1,v[l]);
1835: ierr= VecAXPY(Y,-1.0,s1);
1836: ierr= VecDot(w[l],Y,&gnorm);
1837: ierr= VecAXPBY(v[l],0.0,1.0/gnorm,w[l]);
1839: ierr= VecDot(s1,w[l],&gnorm);
1840: ierr= VecAXPY(s1,gnorm,v[l]);
1841: ierr= VecCopy(s1,s0);
1842: ierr=VecAXPY(X,1.0,s1);
1843: }
1844: return(0);
1845: }