16cab3a1bSJed Brown #include <petscdmda.h> /*I "petscdmda.h" I*/ 2*2961e153SJed Brown #include <private/dmimpl.h> 36cab3a1bSJed Brown #include <private/snesimpl.h> /*I "petscsnes.h" I*/ 46cab3a1bSJed Brown 56cab3a1bSJed Brown /* This structure holds the user-provided DMDA callbacks */ 66cab3a1bSJed Brown typedef struct { 76cab3a1bSJed Brown PetscErrorCode (*residuallocal)(DMDALocalInfo*,void*,void*,void*); 86cab3a1bSJed Brown PetscErrorCode (*jacobianlocal)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*); 96cab3a1bSJed Brown void *residuallocalctx; 106cab3a1bSJed Brown void *jacobianlocalctx; 116cab3a1bSJed Brown } DM_DA_SNES; 126cab3a1bSJed Brown 136cab3a1bSJed Brown #undef __FUNCT__ 146cab3a1bSJed Brown #define __FUNCT__ "SNESDMDestroy_DMDA" 156cab3a1bSJed Brown static PetscErrorCode SNESDMDestroy_DMDA(SNESDM sdm) 166cab3a1bSJed Brown { 176cab3a1bSJed Brown PetscErrorCode ierr; 186cab3a1bSJed Brown 196cab3a1bSJed Brown PetscFunctionBegin; 206cab3a1bSJed Brown ierr = PetscFree(sdm->data);CHKERRQ(ierr); 216cab3a1bSJed Brown PetscFunctionReturn(0); 226cab3a1bSJed Brown } 236cab3a1bSJed Brown 246cab3a1bSJed Brown #undef __FUNCT__ 256cab3a1bSJed Brown #define __FUNCT__ "DMDASNESGetContext" 266cab3a1bSJed Brown static PetscErrorCode DMDASNESGetContext(DM dm,SNESDM sdm,DM_DA_SNES **dmdasnes) 276cab3a1bSJed Brown { 286cab3a1bSJed Brown PetscErrorCode ierr; 296cab3a1bSJed Brown 306cab3a1bSJed Brown PetscFunctionBegin; 316cab3a1bSJed Brown if (!sdm->data) { 326cab3a1bSJed Brown ierr = PetscNewLog(dm,DM_DA_SNES,&sdm->data);CHKERRQ(ierr); 336cab3a1bSJed Brown sdm->destroy = SNESDMDestroy_DMDA; 346cab3a1bSJed Brown } 356cab3a1bSJed Brown *dmdasnes = (DM_DA_SNES*)sdm->data; 366cab3a1bSJed Brown PetscFunctionReturn(0); 376cab3a1bSJed Brown } 386cab3a1bSJed Brown 396cab3a1bSJed Brown #undef __FUNCT__ 406cab3a1bSJed Brown #define __FUNCT__ "SNESComputeFunction_DMDA" 416cab3a1bSJed Brown /* 426cab3a1bSJed Brown This function should eventually replace: 436cab3a1bSJed Brown DMDAComputeFunction() and DMDAComputeFunction1() 446cab3a1bSJed Brown */ 456cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx) 466cab3a1bSJed Brown { 476cab3a1bSJed Brown PetscErrorCode ierr; 486cab3a1bSJed Brown DM dm; 496cab3a1bSJed Brown DM_DA_SNES *dmdasnes = (DM_DA_SNES*)ctx; 506cab3a1bSJed Brown DMDALocalInfo info; 516cab3a1bSJed Brown Vec Xloc; 526cab3a1bSJed Brown void *x,*f; 536cab3a1bSJed Brown 546cab3a1bSJed Brown PetscFunctionBegin; 55*2961e153SJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 56*2961e153SJed Brown PetscValidHeaderSpecific(X,VEC_CLASSID,2); 57*2961e153SJed Brown PetscValidHeaderSpecific(F,VEC_CLASSID,3); 586cab3a1bSJed Brown if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 596cab3a1bSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 606cab3a1bSJed Brown ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 616cab3a1bSJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 626cab3a1bSJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 636cab3a1bSJed Brown ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 646cab3a1bSJed Brown ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 656cab3a1bSJed Brown ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 666cab3a1bSJed Brown CHKMEMQ; 676cab3a1bSJed Brown ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 686cab3a1bSJed Brown CHKMEMQ; 696cab3a1bSJed Brown ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 706cab3a1bSJed Brown ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 716cab3a1bSJed Brown ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 726cab3a1bSJed Brown PetscFunctionReturn(0); 736cab3a1bSJed Brown } 746cab3a1bSJed Brown 756cab3a1bSJed Brown #undef __FUNCT__ 76*2961e153SJed Brown #define __FUNCT__ "SNESComputeJacobian_DMDA" 77*2961e153SJed Brown /* 78*2961e153SJed Brown This function should eventually replace: 79*2961e153SJed Brown DMComputeJacobian() and DMDAComputeJacobian1() 80*2961e153SJed Brown */ 81*2961e153SJed Brown static PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 82*2961e153SJed Brown { 83*2961e153SJed Brown PetscErrorCode ierr; 84*2961e153SJed Brown DM dm; 85*2961e153SJed Brown DM_DA_SNES *dmdasnes = (DM_DA_SNES*)ctx; 86*2961e153SJed Brown DMDALocalInfo info; 87*2961e153SJed Brown Vec Xloc; 88*2961e153SJed Brown void *x; 89*2961e153SJed Brown 90*2961e153SJed Brown PetscFunctionBegin; 91*2961e153SJed Brown if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 92*2961e153SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 93*2961e153SJed Brown 94*2961e153SJed Brown if (dmdasnes->jacobianlocal) { 95*2961e153SJed Brown ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 96*2961e153SJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 97*2961e153SJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 98*2961e153SJed Brown ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 99*2961e153SJed Brown ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 100*2961e153SJed Brown CHKMEMQ; 101*2961e153SJed Brown ierr = (*dmdasnes->jacobianlocal)(&info,x,*A,*B,mstr,dmdasnes->jacobianlocalctx);CHKERRQ(ierr); 102*2961e153SJed Brown CHKMEMQ; 103*2961e153SJed Brown ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 104*2961e153SJed Brown ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 105*2961e153SJed Brown } else { 106*2961e153SJed Brown MatFDColoring fdcoloring; 107*2961e153SJed Brown ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr); 108*2961e153SJed Brown if (!fdcoloring) { 109*2961e153SJed Brown ISColoring coloring; 110*2961e153SJed Brown 111*2961e153SJed Brown ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr); 112*2961e153SJed Brown ierr = MatFDColoringCreate(*B,coloring,&fdcoloring);CHKERRQ(ierr); 113*2961e153SJed Brown ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 114*2961e153SJed Brown switch (dm->coloringtype) { 115*2961e153SJed Brown case IS_COLORING_GLOBAL: 116*2961e153SJed Brown ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode(*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 117*2961e153SJed Brown break; 118*2961e153SJed Brown default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 119*2961e153SJed Brown } 120*2961e153SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr); 121*2961e153SJed Brown ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 122*2961e153SJed Brown ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr); 123*2961e153SJed Brown ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr); 124*2961e153SJed Brown 125*2961e153SJed Brown /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 126*2961e153SJed Brown * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 127*2961e153SJed Brown * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 128*2961e153SJed Brown * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 129*2961e153SJed Brown * taken when the PetscOList for the Vec inside MatFDColoring is destroyed. 130*2961e153SJed Brown */ 131*2961e153SJed Brown ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 132*2961e153SJed Brown } 133*2961e153SJed Brown *mstr = SAME_NONZERO_PATTERN; 134*2961e153SJed Brown ierr = MatFDColoringApply(*B,fdcoloring,X,mstr,snes);CHKERRQ(ierr); 135*2961e153SJed Brown } 136*2961e153SJed Brown /* This will be redundant if the user called both, but it's too common to forget. */ 137*2961e153SJed Brown if (*A != *B) { 138*2961e153SJed Brown ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 139*2961e153SJed Brown ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 140*2961e153SJed Brown } 141*2961e153SJed Brown PetscFunctionReturn(0); 142*2961e153SJed Brown } 143*2961e153SJed Brown 144*2961e153SJed Brown #undef __FUNCT__ 1456cab3a1bSJed Brown #define __FUNCT__ "DMDASNESSetFunctionLocal" 1466cab3a1bSJed Brown /*@C 1476cab3a1bSJed Brown DMDASNESSetFunctionLocal - set a local residual evaluation function 1486cab3a1bSJed Brown 1496cab3a1bSJed Brown Logically Collective 1506cab3a1bSJed Brown 1516cab3a1bSJed Brown Input Arguments: 1526cab3a1bSJed Brown + dm - DM to associate callback with 1536cab3a1bSJed Brown . func - local residual evaluation 1546cab3a1bSJed Brown - ctx - optional context for local residual evaluation 1556cab3a1bSJed Brown 1566cab3a1bSJed Brown Calling sequence for func: 1576cab3a1bSJed Brown + info - DMDALocalInfo defining the subdomain to evaluate the residual on 1586cab3a1bSJed Brown . x - dimensional pointer to state at which to evaluate residual 1596cab3a1bSJed Brown . f - dimensional pointer to residual, write the residual here 1606cab3a1bSJed Brown - ctx - optional context passed above 1616cab3a1bSJed Brown 1626cab3a1bSJed Brown Level: beginner 1636cab3a1bSJed Brown 1646cab3a1bSJed Brown .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 1656cab3a1bSJed Brown @*/ 1666cab3a1bSJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx) 1676cab3a1bSJed Brown { 1686cab3a1bSJed Brown PetscErrorCode ierr; 1696cab3a1bSJed Brown SNESDM sdm; 1706cab3a1bSJed Brown DM_DA_SNES *dmdasnes; 1716cab3a1bSJed Brown 1726cab3a1bSJed Brown PetscFunctionBegin; 1736cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1746cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 1756cab3a1bSJed Brown ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 1766cab3a1bSJed Brown dmdasnes->residuallocal = func; 1776cab3a1bSJed Brown dmdasnes->residuallocalctx = ctx; 1786cab3a1bSJed Brown ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 179*2961e153SJed Brown if (!sdm->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 180*2961e153SJed Brown ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 181*2961e153SJed Brown } 182*2961e153SJed Brown PetscFunctionReturn(0); 183*2961e153SJed Brown } 184*2961e153SJed Brown 185*2961e153SJed Brown #undef __FUNCT__ 186*2961e153SJed Brown #define __FUNCT__ "DMDASNESSetJacobianLocal" 187*2961e153SJed Brown /*@C 188*2961e153SJed Brown DMDASNESSetJacobianLocal - set a local residual evaluation function 189*2961e153SJed Brown 190*2961e153SJed Brown Logically Collective 191*2961e153SJed Brown 192*2961e153SJed Brown Input Arguments: 193*2961e153SJed Brown + dm - DM to associate callback with 194*2961e153SJed Brown . func - local residual evaluation 195*2961e153SJed Brown - ctx - optional context for local residual evaluation 196*2961e153SJed Brown 197*2961e153SJed Brown Calling sequence for func: 198*2961e153SJed Brown + info - DMDALocalInfo defining the subdomain to evaluate the residual on 199*2961e153SJed Brown . x - dimensional pointer to state at which to evaluate residual 200*2961e153SJed Brown . f - dimensional pointer to residual, write the residual here 201*2961e153SJed Brown - ctx - optional context passed above 202*2961e153SJed Brown 203*2961e153SJed Brown Level: beginner 204*2961e153SJed Brown 205*2961e153SJed Brown .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 206*2961e153SJed Brown @*/ 207*2961e153SJed Brown PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*),void *ctx) 208*2961e153SJed Brown { 209*2961e153SJed Brown PetscErrorCode ierr; 210*2961e153SJed Brown SNESDM sdm; 211*2961e153SJed Brown DM_DA_SNES *dmdasnes; 212*2961e153SJed Brown 213*2961e153SJed Brown PetscFunctionBegin; 214*2961e153SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 215*2961e153SJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 216*2961e153SJed Brown ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 217*2961e153SJed Brown dmdasnes->jacobianlocal = func; 218*2961e153SJed Brown dmdasnes->jacobianlocalctx = ctx; 219*2961e153SJed Brown ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 2206cab3a1bSJed Brown PetscFunctionReturn(0); 2216cab3a1bSJed Brown } 222