1 #include <../src/snes/impls/gs/gsimpl.h> 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "GSDestroy_Private" 5 PetscErrorCode GSDestroy_Private(ISColoring coloring) 6 { 7 PetscErrorCode ierr; 8 9 PetscFunctionBegin; 10 ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 11 PetscFunctionReturn(0); 12 13 } 14 15 #undef __FUNCT__ 16 #define __FUNCT__ "SNESComputeGSDefaultSecant" 17 PETSC_EXTERN PetscErrorCode SNESComputeGSDefaultSecant(SNES snes,Vec X,Vec B,void *ctx) 18 { 19 PetscErrorCode ierr; 20 SNES_GS *gs = (SNES_GS*)snes->data; 21 PetscInt i,j,k,ncolors; 22 DM dm; 23 PetscBool flg; 24 ISColoring coloring; 25 MatColoring mc; 26 Vec W,G,F; 27 PetscScalar h=gs->h; 28 IS *coloris; 29 PetscScalar f,g,x,w,d; 30 PetscReal dxt,xt,ft,ft1=0; 31 const PetscInt *idx; 32 PetscInt size,s; 33 PetscReal atol,rtol,stol; 34 PetscInt its; 35 PetscErrorCode (*func)(SNES,Vec,Vec,void*); 36 void *fctx; 37 PetscContainer colorcontainer; 38 PetscBool mat = gs->secant_mat,equal,isdone,alldone; 39 PetscScalar *xa,*fa,*wa,*ga; 40 41 PetscFunctionBegin; 42 if (snes->nwork < 3) { 43 ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr); 44 } 45 W = snes->work[0]; 46 G = snes->work[1]; 47 F = snes->work[2]; 48 ierr = VecGetOwnershipRange(X,&s,NULL);CHKERRQ(ierr); 49 ierr = SNESGSGetTolerances(snes,&atol,&rtol,&stol,&its);CHKERRQ(ierr); 50 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 51 ierr = SNESGetFunction(snes,NULL,&func,&fctx);CHKERRQ(ierr); 52 ierr = PetscObjectQuery((PetscObject)snes,"SNESGSColoring",(PetscObject*)&colorcontainer);CHKERRQ(ierr); 53 if (!colorcontainer) { 54 /* create the coloring */ 55 ierr = DMHasColoring(dm,&flg);CHKERRQ(ierr); 56 if (flg && !mat) { 57 ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,&coloring);CHKERRQ(ierr); 58 } else { 59 if (!snes->jacobian) {ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);} 60 ierr = MatColoringCreate(snes->jacobian,&mc);CHKERRQ(ierr); 61 ierr = MatColoringSetDistance(mc,1);CHKERRQ(ierr); 62 ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr); 63 ierr = MatColoringApply(mc,&coloring);CHKERRQ(ierr); 64 ierr = MatColoringDestroy(&mc);CHKERRQ(ierr); 65 } 66 ierr = PetscContainerCreate(PetscObjectComm((PetscObject)snes),&colorcontainer);CHKERRQ(ierr); 67 ierr = PetscContainerSetPointer(colorcontainer,(void *)coloring);CHKERRQ(ierr); 68 ierr = PetscContainerSetUserDestroy(colorcontainer,(PetscErrorCode (*)(void *))GSDestroy_Private);CHKERRQ(ierr); 69 ierr = PetscObjectCompose((PetscObject)snes,"SNESGSColoring",(PetscObject)colorcontainer);CHKERRQ(ierr); 70 ierr = PetscContainerDestroy(&colorcontainer);CHKERRQ(ierr); 71 } else { 72 ierr = PetscContainerGetPointer(colorcontainer,(void **)&coloring);CHKERRQ(ierr); 73 } 74 ierr = ISColoringGetIS(coloring,&ncolors,&coloris);CHKERRQ(ierr); 75 ierr = VecEqual(X,snes->vec_sol,&equal);CHKERRQ(ierr); 76 if (equal && snes->normschedule == SNES_NORM_ALWAYS) { 77 /* assume that the function is already computed */ 78 ierr = VecCopy(snes->vec_func,F);CHKERRQ(ierr); 79 } else { 80 ierr = PetscLogEventBegin(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr); 81 ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr); 82 ierr = PetscLogEventEnd(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr); 83 if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);} 84 } 85 ierr = VecGetArray(X,&xa);CHKERRQ(ierr); 86 ierr = VecGetArray(F,&fa);CHKERRQ(ierr); 87 ierr = VecGetArray(G,&ga);CHKERRQ(ierr); 88 ierr = VecGetArray(W,&wa);CHKERRQ(ierr); 89 for (i=0;i<ncolors;i++) { 90 ierr = ISGetIndices(coloris[i],&idx);CHKERRQ(ierr); 91 ierr = ISGetLocalSize(coloris[i],&size);CHKERRQ(ierr); 92 ierr = VecCopy(X,W);CHKERRQ(ierr); 93 for (j=0;j<size;j++) { 94 wa[idx[j]-s] += h; 95 } 96 ierr = PetscLogEventBegin(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr); 97 ierr = (*func)(snes,W,G,fctx);CHKERRQ(ierr); 98 ierr = PetscLogEventEnd(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr); 99 if (B) {ierr = VecAXPY(G,-1.0,B);CHKERRQ(ierr);} 100 for (k=0;k<its;k++) { 101 dxt = 0.; 102 xt = 0.; 103 ft = 0.; 104 for (j=0;j<size;j++) { 105 f = fa[idx[j]-s]; 106 x = xa[idx[j]-s]; 107 g = ga[idx[j]-s]; 108 w = wa[idx[j]-s]; 109 if (PetscAbsScalar(g-f) > atol) { 110 d = (x*g-w*f) / PetscRealPart(g-f); 111 } else { 112 d = x; 113 } 114 dxt += PetscRealPart(PetscSqr(d-x)); 115 xt += PetscRealPart(PetscSqr(x)); 116 ft += PetscRealPart(PetscSqr(f)); 117 xa[idx[j]-s] = d; 118 } 119 120 if (k == 0) ft1 = PetscSqrtReal(ft); 121 if (k<its-1) { 122 isdone = PETSC_FALSE; 123 if (stol*PetscSqrtReal(xt) > PetscSqrtReal(dxt)) isdone = PETSC_TRUE; 124 if (PetscSqrtReal(ft) < atol) isdone = PETSC_TRUE; 125 if (rtol*ft1 > PetscSqrtReal(ft)) isdone = PETSC_TRUE; 126 ierr = MPI_Allreduce(&isdone,&alldone,1,MPIU_BOOL,MPI_BAND,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 127 if (alldone) break; 128 } 129 if (i < ncolors-1 || k < its-1) { 130 ierr = PetscLogEventBegin(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr); 131 ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr); 132 ierr = PetscLogEventEnd(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr); 133 if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);} 134 } 135 if (k<its-1) { 136 ierr = VecSwap(X,W);CHKERRQ(ierr); 137 ierr = VecSwap(F,G);CHKERRQ(ierr); 138 } 139 } 140 } 141 ierr = VecRestoreArray(X,&xa);CHKERRQ(ierr); 142 ierr = VecRestoreArray(F,&fa);CHKERRQ(ierr); 143 ierr = VecRestoreArray(G,&ga);CHKERRQ(ierr); 144 ierr = VecRestoreArray(W,&wa);CHKERRQ(ierr); 145 ierr = ISColoringRestoreIS(coloring,&coloris);CHKERRQ(ierr); 146 PetscFunctionReturn(0); 147 } 148