#define PETSCSNES_DLL #include "petscda.h" /*I "petscda.h" I*/ #include "petscmg.h" /*I "petscmg.h" I*/ #include "petscdmmg.h" /*I "petscdmmg.h" I*/ #include "src/inline/ilu.h" #include "src/snes/impls/ls/ls.h" EXTERN_C_BEGIN EXTERN PetscErrorCode PETSCSNES_DLLEXPORT NLFRelax_DAAD(NLF,MatSORType,PetscInt,Vec); EXTERN PetscErrorCode PETSCSNES_DLLEXPORT NLFRelax_DAAD4(NLF,MatSORType,PetscInt,Vec); EXTERN PetscErrorCode PETSCSNES_DLLEXPORT NLFRelax_DAAD9(NLF,MatSORType,PetscInt,Vec); EXTERN PetscErrorCode PETSCSNES_DLLEXPORT NLFRelax_DAADb(NLF,MatSORType,PetscInt,Vec); EXTERN_C_END EXTERN PetscErrorCode DMMGFormFunction(SNES,Vec,Vec,void *); EXTERN PetscErrorCode SNESLSCheckLocalMin_Private(Mat,Vec,Vec,PetscReal,PetscTruth*); #undef __FUNCT__ #define __FUNCT__ "DMMGComputeJacobianWithAdic" /* DMMGComputeJacobianWithAdic - Evaluates the Jacobian via Adic when the user has provided a local function evaluation routine. */ PetscErrorCode DMMGComputeJacobianWithAdic(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr) { DMMG dmmg = (DMMG) ptr; PetscErrorCode ierr; Vec localX; DA da = (DA) dmmg->dm; PetscFunctionBegin; ierr = DAGetLocalVector(da,&localX);CHKERRQ(ierr); ierr = DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); ierr = DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); ierr = DAComputeJacobian1WithAdic(da,localX,*B,dmmg->user);CHKERRQ(ierr); ierr = DARestoreLocalVector(da,&localX);CHKERRQ(ierr); /* Assemble true Jacobian; if it is different */ if (*J != *B) { ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); } ierr = MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR);CHKERRQ(ierr); *flag = SAME_NONZERO_PATTERN; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "SNESDAComputeJacobianWithAdic" /*@ SNESDAComputeJacobianWithAdic - This is a universal Jacobian evaluation routine that may be used with SNESSetJacobian() as long as the user context has a DA as its first record and DASetLocalAdicFunction() has been called. Collective on SNES Input Parameters: + snes - the SNES context . X - input vector . J - Jacobian . B - Jacobian used in preconditioner (usally same as J) . flag - indicates if the matrix changed its structure - ptr - optional user-defined context, as set by SNESSetFunction() Level: intermediate .seealso: DASetLocalFunction(), DASetLocalAdicFunction(), SNESSetFunction(), SNESSetJacobian() @*/ PetscErrorCode PETSCSNES_DLLEXPORT SNESDAComputeJacobianWithAdic(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr) { DA da = *(DA*) ptr; PetscErrorCode ierr; Vec localX; PetscFunctionBegin; ierr = DAGetLocalVector(da,&localX);CHKERRQ(ierr); ierr = DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); ierr = DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); ierr = DAComputeJacobian1WithAdic(da,localX,*B,ptr);CHKERRQ(ierr); ierr = DARestoreLocalVector(da,&localX);CHKERRQ(ierr); /* Assemble true Jacobian; if it is different */ if (*J != *B) { ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); } ierr = MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR);CHKERRQ(ierr); *flag = SAME_NONZERO_PATTERN; PetscFunctionReturn(0); } #include "src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ /* This is pre-beta FAS code. It's design should not be taken seriously! R is the usual multigrid restriction (e.g. the tranpose of peicewise linear interpolation) Q is either a scaled injection or the usual R */ #undef __FUNCT__ #define __FUNCT__ "DMMGSolveFAS" PetscErrorCode DMMGSolveFAS(DMMG *dmmg,PetscInt level) { PetscErrorCode ierr; PetscInt i,j,k; PetscReal norm; PC_MG **mg; PC pc; PetscFunctionBegin; ierr = VecSet(dmmg[level]->r,0.0);CHKERRQ(ierr); for (j=1; j<=level; j++) { if (!dmmg[j]->inject) { ierr = DMGetInjection(dmmg[j-1]->dm,dmmg[j]->dm,&dmmg[j]->inject);CHKERRQ(ierr); } } ierr = KSPGetPC(dmmg[level]->ksp,&pc);CHKERRQ(ierr); mg = ((PC_MG**)pc->data); for (i=0; i<100; i++) { for (j=level; j>0; j--) { /* Relax residual_fine - F(x_fine) = 0 */ for (k=0; kpresmooth; k++) { ierr = NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);CHKERRQ(ierr); } /* R*(residual_fine - F(x_fine)) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); if (j == level || dmmg[j]->monitorall) { /* norm( residual_fine - f(x_fine) ) */ ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); if (j == level) { if (norm < dmmg[level]->abstol) goto theend; if (i == 0) { dmmg[level]->rrtol = norm*dmmg[level]->rtol; } else { if (norm < dmmg[level]->rrtol) goto theend; } } } if (dmmg[j]->monitorall) { for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %G\n",norm);CHKERRQ(ierr); } ierr = MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);CHKERRQ(ierr); /* F(Q*x_fine) */ ierr = VecScatterBegin(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = VecScatterEnd(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);CHKERRQ(ierr); /* residual_coarse = F(Q*x_fine) + R*(residual_fine - F(x_fine)) */ ierr = VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);CHKERRQ(ierr); /* save Q*x_fine into b (needed when interpolating compute x back up */ ierr = VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);CHKERRQ(ierr); } for (j=0; jcoarsesmooth; j++) { ierr = NLFRelax_DAAD(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);CHKERRQ(ierr); } if (dmmg[0]->monitorall){ ierr = DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[0]->w,-1.0,dmmg[0]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[0]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %G\n",norm);CHKERRQ(ierr); } for (j=1; j<=level; j++) { /* x_fine = x_fine + R'*(x_coarse - Q*x_fine) */ ierr = VecAXPY(dmmg[j-1]->x,-1.0,dmmg[j-1]->b);CHKERRQ(ierr); ierr = MatInterpolateAdd(mg[j]->interpolate,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);CHKERRQ(ierr); if (dmmg[j]->monitorall) { /* norm( F(x_fine) - residual_fine ) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %G\n",norm);CHKERRQ(ierr); } /* Relax residual_fine - F(x_fine) = 0 */ for (k=0; kpostsmooth; k++) { ierr = NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);CHKERRQ(ierr); } if (dmmg[j]->monitorall) { /* norm( F(x_fine) - residual_fine ) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %G\n",norm);CHKERRQ(ierr); } } if (dmmg[level]->monitor){ ierr = DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);CHKERRQ(ierr); ierr = VecNorm(dmmg[level]->w,NORM_2,&norm);CHKERRQ(ierr); ierr = PetscPrintf(dmmg[level]->comm,"%D FAS function norm %G\n",i+1,norm);CHKERRQ(ierr); } } theend: PetscFunctionReturn(0); } /* This is the point-block version of FAS */ #undef __FUNCT__ #define __FUNCT__ "DMMGSolveFASb" PetscErrorCode DMMGSolveFASb(DMMG *dmmg,PetscInt level) { PetscErrorCode ierr; PetscInt i,j,k; PetscReal norm; PC_MG **mg; PC pc; PetscFunctionBegin; ierr = VecSet(dmmg[level]->r,0.0);CHKERRQ(ierr); for (j=1; j<=level; j++) { if (!dmmg[j]->inject) { ierr = DMGetInjection(dmmg[j-1]->dm,dmmg[j]->dm,&dmmg[j]->inject);CHKERRQ(ierr); } } ierr = KSPGetPC(dmmg[level]->ksp,&pc);CHKERRQ(ierr); mg = ((PC_MG**)pc->data); for (i=0; i<100; i++) { for (j=level; j>0; j--) { /* Relax residual_fine - F(x_fine) = 0 */ for (k=0; kpresmooth; k++) { ierr = NLFRelax_DAADb(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);CHKERRQ(ierr); } /* R*(residual_fine - F(x_fine)) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); if (j == level || dmmg[j]->monitorall) { /* norm( residual_fine - f(x_fine) ) */ ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); if (j == level) { if (norm < dmmg[level]->abstol) goto theend; if (i == 0) { dmmg[level]->rrtol = norm*dmmg[level]->rtol; } else { if (norm < dmmg[level]->rrtol) goto theend; } } } if (dmmg[j]->monitorall) { for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %G\n",norm);CHKERRQ(ierr); } ierr = MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);CHKERRQ(ierr); /* F(Q*x_fine) */ ierr = VecScatterBegin(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = VecScatterEnd(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);CHKERRQ(ierr); /* residual_coarse = F(Q*x_fine) + R*(residual_fine - F(x_fine)) */ ierr = VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);CHKERRQ(ierr); /* save Q*x_fine into b (needed when interpolating compute x back up */ ierr = VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);CHKERRQ(ierr); } for (j=0; jcoarsesmooth; j++) { ierr = NLFRelax_DAADb(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);CHKERRQ(ierr); } if (dmmg[0]->monitorall){ ierr = DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[0]->w,-1.0,dmmg[0]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[0]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %G\n",norm);CHKERRQ(ierr); } for (j=1; j<=level; j++) { /* x_fine = x_fine + R'*(x_coarse - Q*x_fine) */ ierr = VecAXPY(dmmg[j-1]->x,-1.0,dmmg[j-1]->b);CHKERRQ(ierr); ierr = MatInterpolateAdd(mg[j]->interpolate,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);CHKERRQ(ierr); if (dmmg[j]->monitorall) { /* norm( F(x_fine) - residual_fine ) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %G\n",norm);CHKERRQ(ierr); } /* Relax residual_fine - F(x_fine) = 0 */ for (k=0; kpostsmooth; k++) { ierr = NLFRelax_DAADb(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);CHKERRQ(ierr); } if (dmmg[j]->monitorall) { /* norm( F(x_fine) - residual_fine ) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %G\n",norm);CHKERRQ(ierr); } } if (dmmg[level]->monitor){ ierr = DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);CHKERRQ(ierr); ierr = VecNorm(dmmg[level]->w,NORM_2,&norm);CHKERRQ(ierr); ierr = PetscPrintf(dmmg[level]->comm,"%D FAS function norm %G\n",i+1,norm);CHKERRQ(ierr); } } theend: PetscFunctionReturn(0); } EXTERN_C_BEGIN #include "adic/ad_utils.h" EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "PetscADView" PetscErrorCode PetscADView(PetscInt N,PetscInt nc,double *ptr,PetscViewer viewer) { PetscInt i,j,nlen = PetscADGetDerivTypeSize(); char *cptr = (char*)ptr; double *values; PetscErrorCode ierr; PetscFunctionBegin; for (i=0; ir,zero);CHKERRQ(ierr); for (j=1; j<=level; j++) { if (!dmmg[j]->inject) { ierr = DMGetInjection(dmmg[j-1]->dm,dmmg[j]->dm,&dmmg[j]->inject);CHKERRQ(ierr); } } ierr = KSPGetPC(dmmg[level]->ksp,&pc);CHKERRQ(ierr); mg = ((PC_MG**)pc->data); for (i=0; i<100; i++) { for (j=level; j>0; j--) { ierr = PetscPrintf(PETSC_COMM_WORLD,"I am here");CHKERRQ(ierr); /* Relax residual_fine - F(x_fine) = 0 */ for (k=0; kpresmooth; k++) { ierr = NLFRelax_DAAD4(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);CHKERRQ(ierr); } /* R*(residual_fine - F(x_fine)) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAYPX(dmmg[j]->w,mone,dmmg[j]->r);CHKERRQ(ierr); if (j == level || dmmg[j]->monitorall) { /* norm( residual_fine - f(x_fine) ) */ ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); if (j == level) { if (norm < dmmg[level]->abstol) goto theend; if (i == 0) { dmmg[level]->rrtol = norm*dmmg[level]->rtol; } else { if (norm < dmmg[level]->rrtol) goto theend; } } } ierr = PetscPrintf(PETSC_COMM_WORLD,"I am here");CHKERRQ(ierr); if (dmmg[j]->monitorall) { for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);CHKERRQ(ierr); } ierr = MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);CHKERRQ(ierr); /* F(R*x_fine) */ ierr = VecScatterBegin(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = VecScatterEnd(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);CHKERRQ(ierr); /* residual_coarse = F(R*x_fine) + R*(residual_fine - F(x_fine)) */ ierr = VecAYPX(dmmg[j-1]->r,one,dmmg[j-1]->w);CHKERRQ(ierr); /* save R*x_fine into b (needed when interpolating compute x back up */ ierr = VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);CHKERRQ(ierr); } for (j=0; jpresmooth; j++) { ierr = NLFRelax_DAAD4(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);CHKERRQ(ierr); } if (dmmg[0]->monitorall){ ierr = DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[0]->w,mone,dmmg[0]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[0]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %g\n",norm);CHKERRQ(ierr); } for (j=1; j<=level; j++) { /* x_fine = x_fine + R'*(x_coarse - R*x_fine) */ ierr = VecAXPY(dmmg[j-1]->x,mone,dmmg[j-1]->b);CHKERRQ(ierr); ierr = MatInterpolateAdd(mg[j]->interpolate,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);CHKERRQ(ierr); if (dmmg[j]->monitorall) { /* norm( F(x_fine) - residual_fine ) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[j]->w,mone,dmmg[j]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);CHKERRQ(ierr); } /* Relax residual_fine - F(x_fine) = 0 */ for (k=0; kpostsmooth; k++) { ierr = NLFRelax_DAAD4(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);CHKERRQ(ierr); } if (dmmg[j]->monitorall) { /* norm( F(x_fine) - residual_fine ) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[j]->w,mone,dmmg[j]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);CHKERRQ(ierr); } } if (dmmg[level]->monitor){ ierr = DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);CHKERRQ(ierr); ierr = VecNorm(dmmg[level]->w,NORM_2,&norm);CHKERRQ(ierr); ierr = PetscPrintf(dmmg[level]->comm,"%D FAS function norm %g\n",i,norm);CHKERRQ(ierr); } } theend: PetscFunctionReturn(0); } /* This function provide several FAS v_cycle iteration iter: the number of FAS it run */ #undef __FUNCT__ #define __FUNCT__ "DMMGSolveFASn" PetscErrorCode DMMGSolveFASn(DMMG *dmmg,PetscInt level,PetscInt iter) { PetscErrorCode ierr; PetscInt i,j,k; PetscReal norm; PC_MG **mg; PC pc; PetscFunctionBegin; ierr = VecSet(dmmg[level]->r,0.0);CHKERRQ(ierr); for (j=1; j<=level; j++) { if (!dmmg[j]->inject) { ierr = DMGetInjection(dmmg[j-1]->dm,dmmg[j]->dm,&dmmg[j]->inject);CHKERRQ(ierr); } } ierr = KSPGetPC(dmmg[level]->ksp,&pc);CHKERRQ(ierr); mg = ((PC_MG**)pc->data); for (i=0; i0; j--) { /* Relax residual_fine - F(x_fine) = 0 */ for (k=0; kpresmooth; k++) { ierr = NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);CHKERRQ(ierr); } /* R*(residual_fine - F(x_fine)) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); if (j == level || dmmg[j]->monitorall) { /* norm( residual_fine - f(x_fine) ) */ ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); if (j == level) { if (norm < dmmg[level]->abstol) goto theend; if (i == 0) { dmmg[level]->rrtol = norm*dmmg[level]->rtol; } else { if (norm < dmmg[level]->rrtol) goto theend; } } } if (dmmg[j]->monitorall) { for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);CHKERRQ(ierr); } ierr = MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);CHKERRQ(ierr); /* F(RI*x_fine) */ ierr = VecScatterBegin(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = VecScatterEnd(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);CHKERRQ(ierr); /* residual_coarse = F(RI*x_fine) + R*(residual_fine - F(x_fine)) */ ierr = VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);CHKERRQ(ierr); /* save RI*x_fine into b (needed when interpolating compute x back up */ ierr = VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);CHKERRQ(ierr); } for (j=0; jcoarsesmooth; j++) { ierr = NLFRelax_DAAD(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);CHKERRQ(ierr); } if (dmmg[0]->monitorall){ ierr = DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[0]->w,-1.0,dmmg[0]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[0]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %g\n",norm);CHKERRQ(ierr); } for (j=1; j<=level; j++) { /* x_fine = x_fine + R'*(x_coarse - RI*x_fine) */ ierr = VecAXPY(dmmg[j-1]->x,-1.0,dmmg[j-1]->b);CHKERRQ(ierr); ierr = MatInterpolateAdd(mg[j]->interpolate,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);CHKERRQ(ierr); if (dmmg[j]->monitorall) { /* norm( F(x_fine) - residual_fine ) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);CHKERRQ(ierr); } /* Relax residual_fine - F(x_fine) = 0 */ for (k=0; kpostsmooth; k++) { ierr = NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);CHKERRQ(ierr); } if (dmmg[j]->monitorall) { /* norm( F(x_fine) - residual_fine ) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);CHKERRQ(ierr); } } if (dmmg[level]->monitor){ ierr = DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);CHKERRQ(ierr); ierr = VecNorm(dmmg[level]->w,NORM_2,&norm);CHKERRQ(ierr); ierr = PetscPrintf(dmmg[level]->comm,"%D FAS function norm %g\n",i+1,norm);CHKERRQ(ierr); } } theend: PetscFunctionReturn(0); } /* This is a simple FAS setup function */ #undef __FUNCT__ #define __FUNCT__ "DMMGSolveFASSetUp" PetscErrorCode DMMGSolveFASSetUp(DMMG *dmmg,PetscInt level) { PetscErrorCode ierr; PetscInt j;//,nlevels=dmmg[0]->nlevels-1; //PC pc; PetscFunctionBegin; ierr = VecSet(dmmg[level]->r,0.0);CHKERRQ(ierr); for (j=1; j<=level; j++) { if (!dmmg[j]->inject) { ierr = DMGetInjection(dmmg[j-1]->dm,dmmg[j]->dm,&dmmg[j]->inject);CHKERRQ(ierr); } } ierr = VecSet(dmmg[level]->r,0.0);CHKERRQ(ierr); dmmg[level]->rrtol = 0.0001*dmmg[level]->rtol;//I want to get more precise solution with FAS PetscFunctionReturn(0); } /* This is function to implement multiplicative FAS Options: -dmmg_fas_cycles 1 : FAS v-cycle 2 : FAS w-cycle */ #undef __FUNCT__ #define __FUNCT__ "DMMGSolveFASMCycle" PetscErrorCode DMMGSolveFASMCycle(DMMG *dmmg,PetscInt level,PetscTruth* converged) { PetscErrorCode ierr; PetscInt j,k,cycles=1,nlevels=level;//nlevels=dmmg[0]->nlevels-1; // I need to put nlevels=level in order to get grid sequence correctly PetscReal norm; PC_MG **mg; PC pc; PetscFunctionBegin; ierr = PetscOptionsGetInt(PETSC_NULL,"-dmmg_fas_cycles",&cycles,PETSC_NULL);CHKERRQ(ierr); ierr = KSPGetPC(dmmg[level]->ksp,&pc);CHKERRQ(ierr); mg = ((PC_MG**)pc->data); j=level; if(j) {/* not the coarsest level */ /* Relax residual_fine - F(x_fine) = 0 */ for (k=0; kpresmooth; k++) { ierr = NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);CHKERRQ(ierr); } /* R*(residual_fine - F(x_fine)) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); if (j == nlevels || dmmg[j]->monitorall) { /* norm( residual_fine - f(x_fine) ) */ ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); if (j == nlevels) { if (norm < dmmg[level]->abstol) { *converged = PETSC_TRUE; goto theend; } if (norm < dmmg[level]->rrtol){ *converged = PETSC_TRUE; goto theend; } } } if (dmmg[j]->monitorall) { for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);CHKERRQ(ierr); } ierr = MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);CHKERRQ(ierr); /* F(RI*x_fine) */ ierr = VecScatterBegin(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = VecScatterEnd(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);CHKERRQ(ierr); /* residual_coarse = F(RI*x_fine) + R*(residual_fine - F(x_fine)) */ ierr = VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);CHKERRQ(ierr); /* save RI*x_fine into b (needed when interpolating compute x back up */ ierr = VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);CHKERRQ(ierr); while (cycles--) { ierr = DMMGSolveFASMCycle(dmmg,level-1,converged); } } else { /* for the coarsest level */ for (k=0; kcoarsesmooth; k++) { ierr = NLFRelax_DAAD(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);CHKERRQ(ierr); } if (dmmg[0]->monitorall){ ierr = DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[0]->w,-1.0,dmmg[0]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[0]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %g\n",norm);CHKERRQ(ierr); } if (j == nlevels || dmmg[j]->monitorall) { /* norm( residual_fine - f(x_fine) ) */ ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); if (j == nlevels) { if (norm < dmmg[level]->abstol) { *converged = PETSC_TRUE; goto theend; } if (norm < dmmg[level]->rrtol){ *converged = PETSC_TRUE; goto theend; } } } } j=level; if(j) { /* not for the coarsest level */ /* x_fine = x_fine + R'*(x_coarse - RI*x_fine) */ ierr = VecAXPY(dmmg[j-1]->x,-1.0,dmmg[j-1]->b);CHKERRQ(ierr); ierr = MatInterpolateAdd(mg[j]->interpolate,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);CHKERRQ(ierr); if (dmmg[j]->monitorall) { /* norm( F(x_fine) - residual_fine ) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);CHKERRQ(ierr); } /* Relax residual_fine - F(x_fine) = 0 */ for (k=0; kpostsmooth; k++) { ierr = NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);CHKERRQ(ierr); } if (dmmg[j]->monitorall) { /* norm( F(x_fine) - residual_fine ) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);CHKERRQ(ierr); } } theend: PetscFunctionReturn(0); } /* This is function to implement multiplicative FAS with block smoother Options: -dmmg_fas_cycles 1 : FAS v-cycle 2 : FAS w-cycle */ #undef __FUNCT__ #define __FUNCT__ "DMMGSolveFASMCycle9" PetscErrorCode DMMGSolveFASMCycle9(DMMG *dmmg,PetscInt level,PetscTruth* converged) { PetscErrorCode ierr; PetscInt j,k,cycles=1,nlevels=level;//nlevels=dmmg[0]->nlevels-1; // I need to put nlevels=level in order to get grid sequence correctly PetscReal norm; PC_MG **mg; PC pc; PetscFunctionBegin; ierr = PetscOptionsGetInt(PETSC_NULL,"-dmmg_fas_cycles",&cycles,PETSC_NULL);CHKERRQ(ierr); ierr = KSPGetPC(dmmg[level]->ksp,&pc);CHKERRQ(ierr); mg = ((PC_MG**)pc->data); // for (j=level; j>0; j--) { j=level; //ierr = PetscPrintf(dmmg[level]->comm,"j=%d,nlevels=%d",j,nlevels);CHKERRQ(ierr); if(j) {/* not the coarsest level */ /* Relax residual_fine - F(x_fine) = 0 */ for (k=0; kpresmooth; k++) { ierr = NLFRelax_DAAD9(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);CHKERRQ(ierr); } /* R*(residual_fine - F(x_fine)) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); if (j == nlevels || dmmg[j]->monitorall) { /* norm( residual_fine - f(x_fine) ) */ ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); if (j == nlevels) { if (norm < dmmg[level]->abstol) { *converged = PETSC_TRUE; goto theend; } /* if (i == 0) { dmmg[level]->rrtol = norm*dmmg[level]->rtol; } else {*/ if (norm < dmmg[level]->rrtol){ *converged = PETSC_TRUE; goto theend; } } } if (dmmg[j]->monitorall) { for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);CHKERRQ(ierr); } ierr = MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);CHKERRQ(ierr); /* F(RI*x_fine) */ ierr = VecScatterBegin(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = VecScatterEnd(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);CHKERRQ(ierr); /* residual_coarse = F(RI*x_fine) + R*(residual_fine - F(x_fine)) */ ierr = VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);CHKERRQ(ierr); /* save RI*x_fine into b (needed when interpolating compute x back up */ ierr = VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);CHKERRQ(ierr); while (cycles--) { ierr = DMMGSolveFASMCycle9(dmmg,level-1,converged); } } else { /* for the coarsest level */ for (k=0; kcoarsesmooth; k++) { ierr = NLFRelax_DAAD9(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);CHKERRQ(ierr); } if (dmmg[0]->monitorall){ ierr = DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[0]->w,-1.0,dmmg[0]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[0]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %g\n",norm);CHKERRQ(ierr); } if (j == nlevels || dmmg[j]->monitorall) { /* norm( residual_fine - f(x_fine) ) */ ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); if (j == nlevels) { if (norm < dmmg[level]->abstol) { *converged = PETSC_TRUE; goto theend; } /* if (i == 0) { dmmg[level]->rrtol = norm*dmmg[level]->rtol; } else {*/ if (norm < dmmg[level]->rrtol){ *converged = PETSC_TRUE; goto theend; } } } } j=level; if(j) { /* not for the coarsest level */ /* x_fine = x_fine + R'*(x_coarse - RI*x_fine) */ ierr = VecAXPY(dmmg[j-1]->x,-1.0,dmmg[j-1]->b);CHKERRQ(ierr); ierr = MatInterpolateAdd(mg[j]->interpolate,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);CHKERRQ(ierr); if (dmmg[j]->monitorall) { /* norm( F(x_fine) - residual_fine ) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);CHKERRQ(ierr); } /* Relax residual_fine - F(x_fine) = 0 */ for (k=0; kpostsmooth; k++) { ierr = NLFRelax_DAAD9(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);CHKERRQ(ierr); } if (dmmg[j]->monitorall) { /* norm( F(x_fine) - residual_fine ) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); ierr = VecNorm(dmmg[j]->w,NORM_2,&norm);CHKERRQ(ierr); for (k=0; kcomm," ");CHKERRQ(ierr);} ierr = PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);CHKERRQ(ierr); } /* if(j==nlevels){ if (dmmg[level]->monitor){ ierr = DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);CHKERRQ(ierr); ierr = VecNorm(dmmg[level]->w,NORM_2,&norm);CHKERRQ(ierr); ierr = PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);CHKERRQ(ierr); } }*/ } theend: PetscFunctionReturn(0); } /* This is function to implement full FAS with block smoother(9 points together) Options: -dmmg_fas_cycles 1 : FAS v-cycle 2 : FAS w-cycle */ #undef __FUNCT__ #define __FUNCT__ "DMMGSolveFASFCycle" PetscErrorCode DMMGSolveFASFCycle(DMMG *dmmg,PetscInt l,PetscTruth* converged) { PetscErrorCode ierr; PetscInt j;//l = dmmg[0]->nlevels-1; PC_MG **mg; PC pc; PetscFunctionBegin; ierr = KSPGetPC(dmmg[l]->ksp,&pc);CHKERRQ(ierr); mg = ((PC_MG**)pc->data); // restriction all the way down to the coarse level if(l>0) { for(j=l;j>0;j--) { /* R*(residual_fine - F(x_fine)) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); ierr = MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);CHKERRQ(ierr); /* F(RI*x_fine) */ ierr = VecScatterBegin(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = VecScatterEnd(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);CHKERRQ(ierr); /* residual_coarse = F(RI*x_fine) + R*(residual_fine - F(x_fine)) */ ierr = VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);CHKERRQ(ierr); /* save RI*x_fine into b (needed when interpolating compute x back up */ ierr = VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);CHKERRQ(ierr); } // all the way up to the finest level for (j=0; jx,-1.0,dmmg[j]->b);CHKERRQ(ierr); ierr = MatInterpolateAdd(mg[j+1]->interpolate,dmmg[j]->x,dmmg[j+1]->x,dmmg[j+1]->x);CHKERRQ(ierr); } } ierr = DMMGSolveFASMCycle(dmmg,l,converged);CHKERRQ(ierr); PetscFunctionReturn(0); } /* This is function to implement full FAS with block smoother ( 9 points together) Options: -dmmg_fas_cycles 1 : FAS v-cycle 2 : FAS w-cycle */ #undef __FUNCT__ #define __FUNCT__ "DMMGSolveFASFCycle" PetscErrorCode DMMGSolveFASFCycle9(DMMG *dmmg,PetscInt l,PetscTruth* converged) { PetscErrorCode ierr; PetscInt j;//l = dmmg[0]->nlevels-1; PC_MG **mg; PC pc; PetscFunctionBegin; ierr = KSPGetPC(dmmg[l]->ksp,&pc);CHKERRQ(ierr); mg = ((PC_MG**)pc->data); // restriction all the way down to the coarse level if(l>0) { for(j=l;j>0;j--) { /* R*(residual_fine - F(x_fine)) */ ierr = DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);CHKERRQ(ierr); ierr = VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);CHKERRQ(ierr); ierr = MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);CHKERRQ(ierr); /* F(RI*x_fine) */ ierr = VecScatterBegin(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = VecScatterEnd(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);CHKERRQ(ierr); ierr = DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);CHKERRQ(ierr); /* residual_coarse = F(RI*x_fine) + R*(residual_fine - F(x_fine)) */ ierr = VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);CHKERRQ(ierr); /* save RI*x_fine into b (needed when interpolating compute x back up */ ierr = VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);CHKERRQ(ierr); } // all the way up to the finest level for (j=0; jx,-1.0,dmmg[j]->b);CHKERRQ(ierr); ierr = MatInterpolateAdd(mg[j+1]->interpolate,dmmg[j]->x,dmmg[j+1]->x,dmmg[j+1]->x);CHKERRQ(ierr); } } ierr = DMMGSolveFASMCycle9(dmmg,l,converged);CHKERRQ(ierr); PetscFunctionReturn(0); } /* This is function is to solve nonlinear system with FAS Options: -dmmg_fas_9: using block smoother -dmmg_fas_full: using full FAS */ #undef __FUNCT__ #define __FUNCT__ "DMMGSolveFASCycle" PetscErrorCode DMMGSolveFASCycle(DMMG *dmmg,PetscInt level) { PetscErrorCode ierr; PetscInt i; PetscTruth converged = PETSC_FALSE, flg,flgb; PetscReal norm; PetscFunctionBegin; ierr = DMMGSolveFASSetUp(dmmg,level);CHKERRQ(ierr); ierr = PetscOptionsHasName(PETSC_NULL,"-dmmg_fas_9",&flgb);CHKERRQ(ierr); if(flgb){ ierr = PetscOptionsHasName(PETSC_NULL,"-dmmg_fas_full",&flg);CHKERRQ(ierr); if (flg) { for(i=0;i<1000;i++){ ierr = PetscPrintf(dmmg[level]->comm,"%D ",i+1);CHKERRQ(ierr); ierr = DMMGSolveFASFCycle9(dmmg,level,&converged);CHKERRQ(ierr); if (dmmg[level]->monitor){ ierr = DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);CHKERRQ(ierr); ierr = VecNorm(dmmg[level]->w,NORM_2,&norm);CHKERRQ(ierr); ierr = PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);CHKERRQ(ierr); } if (converged) PetscFunctionReturn(0); } } else{ for(i=0;i<1000;i++){ ierr = PetscPrintf(dmmg[level]->comm,"%D ",i+1);CHKERRQ(ierr); ierr = DMMGSolveFASMCycle9(dmmg,level,&converged);CHKERRQ(ierr); if (dmmg[level]->monitor){ ierr = DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);CHKERRQ(ierr); ierr = VecNorm(dmmg[level]->w,NORM_2,&norm);CHKERRQ(ierr); ierr = PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);CHKERRQ(ierr); } if (converged) PetscFunctionReturn(0); } } } else { ierr = PetscOptionsHasName(PETSC_NULL,"-dmmg_fas_full",&flg);CHKERRQ(ierr); if (flg) { for(i=0;i<1000;i++){ ierr = PetscPrintf(dmmg[level]->comm,"%D ",i+1);CHKERRQ(ierr); ierr = DMMGSolveFASFCycle(dmmg,level,&converged);CHKERRQ(ierr); if (dmmg[level]->monitor){ ierr = DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);CHKERRQ(ierr); ierr = VecNorm(dmmg[level]->w,NORM_2,&norm);CHKERRQ(ierr); ierr = PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);CHKERRQ(ierr); } if (converged) PetscFunctionReturn(0); } } else{ for(i=0;i<1000;i++){ ierr = PetscPrintf(dmmg[level]->comm,"%D ",i+1);CHKERRQ(ierr); ierr = DMMGSolveFASMCycle(dmmg,level,&converged);CHKERRQ(ierr); if (dmmg[level]->monitor){ ierr = DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);CHKERRQ(ierr); ierr = VecNorm(dmmg[level]->w,NORM_2,&norm);CHKERRQ(ierr); ierr = PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);CHKERRQ(ierr); } if (converged) PetscFunctionReturn(0); } } } PetscFunctionReturn(0); } /* This is function is to implement one FAS iteration Options: -dmmg_fas_9: using block smoother -dmmg_fas_full: using full FAS */ #undef __FUNCT__ #define __FUNCT__ "DMMGSolveFASCyclen" PetscErrorCode DMMGSolveFASCyclen(DMMG *dmmg,PetscInt level) { PetscErrorCode ierr; PetscTruth converged = PETSC_FALSE, flg,flgb; PetscReal norm; // PC_MG **mg; //PC pc; PetscFunctionBegin; ierr = DMMGSolveFASSetUp(dmmg,level);CHKERRQ(ierr); ierr = PetscOptionsHasName(PETSC_NULL,"-dmmg_fas_9",&flgb);CHKERRQ(ierr); if(flgb){ ierr = PetscOptionsHasName(PETSC_NULL,"-dmmg_fas_full",&flg);CHKERRQ(ierr); if (flg) { ierr = DMMGSolveFASFCycle9(dmmg,level,&converged);CHKERRQ(ierr); if (dmmg[level]->monitor){ ierr = DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);CHKERRQ(ierr); ierr = VecNorm(dmmg[level]->w,NORM_2,&norm);CHKERRQ(ierr); ierr = PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);CHKERRQ(ierr); } } else{ ierr = DMMGSolveFASMCycle9(dmmg,level,&converged);CHKERRQ(ierr); if (dmmg[level]->monitor){ ierr = DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);CHKERRQ(ierr); ierr = VecNorm(dmmg[level]->w,NORM_2,&norm);CHKERRQ(ierr); ierr = PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);CHKERRQ(ierr); } } } else { ierr = PetscOptionsHasName(PETSC_NULL,"-dmmg_fas_full",&flg);CHKERRQ(ierr); if (flg) { ierr = DMMGSolveFASFCycle(dmmg,level,&converged);CHKERRQ(ierr); if (dmmg[level]->monitor){ ierr = DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);CHKERRQ(ierr); ierr = VecNorm(dmmg[level]->w,NORM_2,&norm);CHKERRQ(ierr); ierr = PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);CHKERRQ(ierr); } } else{ ierr = DMMGSolveFASMCycle(dmmg,level,&converged);CHKERRQ(ierr); if (dmmg[level]->monitor){ ierr = DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);CHKERRQ(ierr); ierr = VecNorm(dmmg[level]->w,NORM_2,&norm);CHKERRQ(ierr); ierr = PetscPrintf(dmmg[level]->comm," FAS function norm %g\n",norm);CHKERRQ(ierr); } } } PetscFunctionReturn(0); } /* This is function to implement Nonlinear CG to accelerate FAS In order to use this acceleration, the option is -dmmg_fas_NCG */ #undef __FUNCT__ #define __FUNCT__ "DMMGSolveFAS_NCG" PetscErrorCode DMMGSolveFAS_NCG(DMMG *dmmg, PetscInt level) { SNES snes = dmmg[level]->snes; SNES_LS *neP = (SNES_LS*)snes->data; PetscErrorCode ierr; PetscInt maxits,i,lits; PetscTruth lssucceed; // MatStructure flg = DIFFERENT_NONZERO_PATTERN; PetscReal fnorm,gnorm,xnorm,ynorm,betaFR,betaPR,beta,betaHS,betaDY; Vec Y,X,F,G,W,Gradold,Sk; //KSP ksp; PetscFunctionBegin; if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); ierr = VecDuplicate(dmmg[level]->x,&Sk);CHKERRQ(ierr); snes->vec_sol = snes->vec_sol_always = dmmg[level]->x; if (!snes->setupcalled) { ierr = SNESSetUp(snes);CHKERRQ(ierr); } if (snes->conv_hist_reset) snes->conv_hist_len = 0; ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; snes->reason = SNES_CONVERGED_ITERATING; maxits = snes->max_its; /* maximum number of iterations */ X = snes->vec_sol; /* solution vector */ F = snes->vec_func; /* residual vector */ Y = snes->work[0]; /* work vectors */ G = snes->work[1]; W = snes->work[2]; ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); snes->iter = 0; ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); X = dmmg[level]->x; ierr = VecCopy(X,Y);CHKERRQ(ierr); ierr = VecCopy(X,G);CHKERRQ(ierr); // to get the residual for the F ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ if (fnorm != fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); snes->norm = fnorm; ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); SNESLogConvHistory(snes,fnorm,0); SNESMonitor(snes,0,fnorm); if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} /* set parameter for default relative tolerance convergence test */ snes->ttol = fnorm*snes->rtol; // dmmg[level]->rrtol= snes->ttol; // set this to store the old grad Gradold=snes->vec_sol_update_always; // compute the search direction Y // I need to put Q(x)=x-FAS(x) here ierr = DMMGSolveFASCyclen(dmmg,level);CHKERRQ(ierr); // F = X - dmmg[level]->x; this is the gradient direction ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr); // copy the gradient to the old ierr = VecCopy(Y,Gradold);CHKERRQ(ierr); // copy X back ierr = VecCopy(G,X);CHKERRQ(ierr); ierr = VecWAXPY(Sk,-1.0,X,X);CHKERRQ(ierr); // so far I put X=X_c, F= F(x_c), Gradold= Y=grad(x_c) // for (i=0; iLineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); ierr = PetscInfo4(snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; ierr = PetscPrintf(PETSC_COMM_WORLD,"step=%g,oldnorm=%g,norm=%g ",ynorm,fnorm,gnorm);CHKERRQ(ierr); fnorm=gnorm; //copy the new function norm; this will effect the line_search ierr = VecWAXPY(Sk,-1.0,X,W);CHKERRQ(ierr); //update the new solution ierr=VecCopy(W,X);CHKERRQ(ierr); ierr=VecCopy(G,F);CHKERRQ(ierr); // compute the new search direction G // I need to put Q(x)=x-FAS(x) here ierr = DMMGSolveFASCyclen(dmmg,level);CHKERRQ(ierr); //G = X - dmmg[level]->x; G is the new gradient, Y is old gradient ierr = VecWAXPY(G,-1.0,X,W);CHKERRQ(ierr); // copy W back to X ierr = VecCopy(W,X);CHKERRQ(ierr); ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); ierr = VecNorm(Gradold,NORM_2,&ynorm);CHKERRQ(ierr); betaFR = gnorm*gnorm/(ynorm*ynorm); //FR_beta ierr = VecWAXPY(W,-1.0,Gradold,G);CHKERRQ(ierr); ierr = VecDot(W,G,&gnorm);CHKERRQ(ierr); // ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); ierr = VecNorm(Gradold,NORM_2,&ynorm);CHKERRQ(ierr); betaPR = gnorm/(ynorm*ynorm); //PR_beta if ( betaPR<-betaFR) { beta =- betaFR; } else { if (betaPR>betaFR) {beta=betaFR;} else{ beta=betaPR; } } // beta=betaFR; //beta=betaPR; // try another beta ierr = VecWAXPY(W,-1.0,Gradold,G);CHKERRQ(ierr); ierr = VecDot(W,G,&betaHS);CHKERRQ(ierr); ierr = VecDot(W,Y,&gnorm);CHKERRQ(ierr); betaHS=-betaHS/gnorm; ierr = VecDot(G,G,&betaDY);CHKERRQ(ierr); betaDY=-betaDY/gnorm; if(betaHSiter = i+1; snes->norm = fnorm; ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); SNESLogConvHistory(snes,fnorm,lits); SNESMonitor(snes,i+1,fnorm); if (!lssucceed) { PetscTruth ismin; beta=0; if (++snes->numFailures >= snes->maxFailures) { snes->reason = SNES_DIVERGED_LS_FAILURE; ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; break; } } /* Test for convergence */ if (snes->converged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ ierr = (*snes->converged)(snes,snes->iter,xnorm,1.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); if (snes->reason) { break; } } } if (X != snes->vec_sol) { ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); } if (F != snes->vec_func) { ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); } snes->vec_sol_always = snes->vec_sol; snes->vec_func_always = snes->vec_func; if (i == maxits) { ierr = PetscInfo1(snes,"SNESSolve_LS: Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); snes->reason = SNES_DIVERGED_MAX_IT; } ierr = PetscPrintf(PETSC_COMM_WORLD,"reason=%d\n",snes->reason);CHKERRQ(ierr); // VecView(X,PETSC_VIEWER_STDOUT_WORLD); PetscFunctionReturn(0); } /* This is function to implement NonGMRES to accelerate FAS In order to use this acceleration, the option is -dmmg_fas_NGMRES Options: -dmmg_fas_ngmres_m : the number of previous solutions to keep */ #undef __FUNCT__ #define __FUNCT__ "DMMGSolveFAS_NGMRES" PetscErrorCode DMMGSolveFAS_NGMRES(DMMG *dmmg, PetscInt level) { SNES snes = dmmg[level]->snes; PetscErrorCode ierr; PetscInt maxits=10000,i,k,l,j,subm=3,iter; ierr = PetscOptionsGetInt(PETSC_NULL,"-dmmg_fas_ngmres_m",&subm,PETSC_NULL);CHKERRQ(ierr); PetscTruth restart=PETSC_FALSE, selectA=PETSC_FALSE; PetscReal fnorm,gnorm,dnorm,dnormtemp,dminnorm,fminnorm,tol=1.e-12,gammaA=2,epsilonB=0.1,deltaB=0.9,gammaC; Vec X,F,G,W,D,u[subm],res[subm]; PetscScalar H[subm][subm],q[subm][subm],beta[subm],xi[subm],alpha[subm],alphasum,det,Hinv[16]; PetscFunctionBegin; gammaC=2; if (gammaA>gammaC) gammaC=gammaA; ierr = VecDuplicate(dmmg[level]->x,&X);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&F);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&W);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&G);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&D);CHKERRQ(ierr); for(i=0;ix,&u[i]);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&res[i]);CHKERRQ(ierr); } X = dmmg[level]->x; ierr = VecCopy(X,u[0]);CHKERRQ(ierr); // to get the residual for the F ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); ierr = VecCopy(F,res[0]);CHKERRQ(ierr); ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ fnorm=fnorm*fnorm; iter=1; restartline: q[0][0] = fnorm; fminnorm=fnorm; for (k=1; kgammaC*sqrt(fminnorm)) { ierr = PetscPrintf(PETSC_COMM_WORLD,"restart for C ");CHKERRQ(ierr); restart=PETSC_TRUE; } if(epsilonB*dnorm>dminnorm && sqrt(gnorm)>deltaB*sqrt(fminnorm)) { ierr = PetscPrintf(PETSC_COMM_WORLD,"restart for D ");CHKERRQ(ierr); restart=PETSC_TRUE; } /* Prepation for the next iteration */ //turn off restart //restart=PETSC_FALSE; if(restart){ restart=PETSC_FALSE; goto restartline; } else { j=k%subm; ierr = VecCopy(F,res[j]);CHKERRQ(ierr); ierr = VecCopy(X,u[j]);CHKERRQ(ierr); for(i=0;inlevels-1; PetscFunctionBegin; ierr = VecDuplicate(dmmg[level]->x,&temp);CHKERRQ(ierr); ierr = VecCopy(dmmg[level]->x,temp);CHKERRQ(ierr); ierr = VecCopy(x,dmmg[level]->x);CHKERRQ(ierr); ierr = VecCopy(x,f);CHKERRQ(ierr); // I need to put -F(x)=x-FAS(x) here ierr = DMMGSolveFASCyclen(dmmg,level);CHKERRQ(ierr); ierr = VecAXPY(f,-1.0,dmmg[level]->x);CHKERRQ(ierr); // y = alpha x + y. ierr=VecCopy(temp,dmmg[level]->x);CHKERRQ(ierr); //copy W back to X PetscFunctionReturn(0); } /* this function is to implement Quasi-Newton method with implicit Broyden updating methods(limit memory version) In order to use this method, the option is -dmmg_fas_QNewton Options: -dmmg_fas_QNewton_m: the number of the vectors to keep for inverse of Jacobian -dmmg_fas_initialJacobian: will use matrix-free GMRES to solve the initial Jacobian with options -snes_mf -snes_max_it 1 -ksp_max_it n In this function, does not have line search and nonlinear gmres acceleration */ #undef __FUNCT__ #define __FUNCT__ "DMMGSolveFAS_QNewton" PetscErrorCode DMMGSolveFAS_QNewton(DMMG *dmmg, PetscInt level) { SNES snes = dmmg[level]->snes, snes0; PetscErrorCode ierr; PetscInt maxits=10000,i,k,l,subm=3,subm01; ierr = PetscOptionsGetInt(PETSC_NULL,"-dmmg_fas_QNewton_m",&subm,PETSC_NULL);CHKERRQ(ierr); subm01=subm-1; PetscTruth flg; PetscReal fnorm,gnorm,tol=1.e-12; Vec X,F,G,W,D,Y,v[subm],w[subm],s0,s1,F0,F1; PetscFunctionBegin; ierr = PetscOptionsHasName(PETSC_NULL,"-dmmg_fas_initialJacobian",&flg);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&X);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&F);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&W);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&G);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&D);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&Y);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&s0);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&s1);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&F0);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&F1);CHKERRQ(ierr); // creat a snes for solve the initial Jacobian ierr = SNESCreate(dmmg[level]->comm,&snes0);CHKERRQ(ierr); ierr = SNESSetFunction(snes0,F,DMMGFASFunction,dmmg);CHKERRQ(ierr); ierr = SNESSetFromOptions(snes0);CHKERRQ(ierr); for(i=0;ix,&v[i]);CHKERRQ(ierr); ierr = VecDuplicate(dmmg[level]->x,&w[i]);CHKERRQ(ierr); } //We first try B0==I X = dmmg[level]->x; if(flg){ ierr= VecAXPBY(Y,0.0,0.0,X);CHKERRQ(ierr); ierr= VecCopy(X,s0);CHKERRQ(ierr); ierr= SNESSolve(snes0,Y,s0);CHKERRQ(ierr); ierr= VecAXPY(s0,-1.0,X);CHKERRQ(ierr); } else{ ierr=VecCopy(X,W);CHKERRQ(ierr); ierr=VecCopy(X,Y);CHKERRQ(ierr); // I need to put -F(x)=x-FAS(x) here ierr = DMMGSolveFASCyclen(dmmg,level);CHKERRQ(ierr); ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr); // y = alpha x + y. ierr=VecCopy(W,X);CHKERRQ(ierr); //copy W back to X // Y stores the -F(x) ierr= VecAXPBY(Y,0.0,-1.0,X);CHKERRQ(ierr); ierr= VecCopy(Y,s0);CHKERRQ(ierr); } ierr = VecAXPY(X,1.0,s0);CHKERRQ(ierr); for(k=0; kx,s1);CHKERRQ(ierr); ierr= SNESSolve(snes0,Y,s1);CHKERRQ(ierr); ierr= VecAXPY(s1,-1.0,X);CHKERRQ(ierr); } else{ ierr=VecCopy(X,W);CHKERRQ(ierr); ierr=VecCopy(X,Y);CHKERRQ(ierr); // I need to put -F(x)=x-FAS(x) here ierr = DMMGSolveFASCyclen(dmmg,level);CHKERRQ(ierr); ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr); // y = alpha x + y. ierr=VecCopy(W,X);CHKERRQ(ierr); //copy W back to X //So far, I got X=x_k, Y=-F(x_k) // I should solve the G=-B_0^{-1}F(x_k) first, but I choose B_0^{-1}=I, ierr=VecCopy(Y,F1);CHKERRQ(ierr); ierr= VecAXPBY(Y,0.0,-1.0,X);CHKERRQ(ierr); ierr=VecCopy(Y,s1);CHKERRQ(ierr); } l=subm; if (k