16cab3a1bSJed Brown #include <petscdmda.h> /*I "petscdmda.h" I*/ 2b45d2f2cSJed Brown #include <petsc-private/dmimpl.h> 3b45d2f2cSJed Brown #include <petsc-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*); 8d1e9a80fSBarry Smith PetscErrorCode (*jacobianlocal)(DMDALocalInfo*,void*,Mat,Mat,void*); 92a4ee8f2SPeter Brune PetscErrorCode (*objectivelocal)(DMDALocalInfo*,void*,PetscReal*,void*); 106cab3a1bSJed Brown void *residuallocalctx; 116cab3a1bSJed Brown void *jacobianlocalctx; 122a4ee8f2SPeter Brune void *objectivelocalctx; 13709ac0d2SJed Brown InsertMode residuallocalimode; 14dcad5f1cSBarry Smith 15dcad5f1cSBarry Smith /* For Picard iteration defined locally */ 16dcad5f1cSBarry Smith PetscErrorCode (*rhsplocal)(DMDALocalInfo*,void*,void*,void*); 17d1e9a80fSBarry Smith PetscErrorCode (*jacobianplocal)(DMDALocalInfo*,void*,Mat,Mat,void*); 18dcad5f1cSBarry Smith void *picardlocalctx; 19942e3340SBarry Smith } DMSNES_DA; 206cab3a1bSJed Brown 216cab3a1bSJed Brown #undef __FUNCT__ 22942e3340SBarry Smith #define __FUNCT__ "DMSNESDestroy_DMDA" 23942e3340SBarry Smith static PetscErrorCode DMSNESDestroy_DMDA(DMSNES sdm) 246cab3a1bSJed Brown { 256cab3a1bSJed Brown PetscErrorCode ierr; 266cab3a1bSJed Brown 276cab3a1bSJed Brown PetscFunctionBegin; 286cab3a1bSJed Brown ierr = PetscFree(sdm->data);CHKERRQ(ierr); 296cab3a1bSJed Brown PetscFunctionReturn(0); 306cab3a1bSJed Brown } 316cab3a1bSJed Brown 326cab3a1bSJed Brown #undef __FUNCT__ 33942e3340SBarry Smith #define __FUNCT__ "DMSNESDuplicate_DMDA" 3422c6f798SBarry Smith static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm,DMSNES sdm) 35e3c0b8a2SPeter Brune { 36e3c0b8a2SPeter Brune PetscErrorCode ierr; 37e3c0b8a2SPeter Brune 38e3c0b8a2SPeter Brune PetscFunctionBegin; 39b00a9115SJed Brown ierr = PetscNewLog(sdm,(DMSNES_DA**)&sdm->data);CHKERRQ(ierr); 40dcad5f1cSBarry Smith if (oldsdm->data) { 41942e3340SBarry Smith ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_DA));CHKERRQ(ierr); 42dcad5f1cSBarry Smith } 43e3c0b8a2SPeter Brune PetscFunctionReturn(0); 44e3c0b8a2SPeter Brune } 45e3c0b8a2SPeter Brune 46e3c0b8a2SPeter Brune 47e3c0b8a2SPeter Brune #undef __FUNCT__ 486cab3a1bSJed Brown #define __FUNCT__ "DMDASNESGetContext" 49942e3340SBarry Smith static PetscErrorCode DMDASNESGetContext(DM dm,DMSNES sdm,DMSNES_DA **dmdasnes) 506cab3a1bSJed Brown { 516cab3a1bSJed Brown PetscErrorCode ierr; 526cab3a1bSJed Brown 536cab3a1bSJed Brown PetscFunctionBegin; 540298fd71SBarry Smith *dmdasnes = NULL; 556cab3a1bSJed Brown if (!sdm->data) { 56b00a9115SJed Brown ierr = PetscNewLog(dm,(DMSNES_DA**)&sdm->data);CHKERRQ(ierr); 5722c6f798SBarry Smith sdm->ops->destroy = DMSNESDestroy_DMDA; 5822c6f798SBarry Smith sdm->ops->duplicate = DMSNESDuplicate_DMDA; 596cab3a1bSJed Brown } 60942e3340SBarry Smith *dmdasnes = (DMSNES_DA*)sdm->data; 616cab3a1bSJed Brown PetscFunctionReturn(0); 626cab3a1bSJed Brown } 636cab3a1bSJed Brown 646cab3a1bSJed Brown #undef __FUNCT__ 656cab3a1bSJed Brown #define __FUNCT__ "SNESComputeFunction_DMDA" 666cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx) 676cab3a1bSJed Brown { 686cab3a1bSJed Brown PetscErrorCode ierr; 696cab3a1bSJed Brown DM dm; 70942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 716cab3a1bSJed Brown DMDALocalInfo info; 726cab3a1bSJed Brown Vec Xloc; 736cab3a1bSJed Brown void *x,*f; 746cab3a1bSJed Brown 756cab3a1bSJed Brown PetscFunctionBegin; 762961e153SJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 772961e153SJed Brown PetscValidHeaderSpecific(X,VEC_CLASSID,2); 782961e153SJed Brown PetscValidHeaderSpecific(F,VEC_CLASSID,3); 79ce94432eSBarry Smith if (!dmdasnes->residuallocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 806cab3a1bSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 816cab3a1bSJed Brown ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 826cab3a1bSJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 836cab3a1bSJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 846cab3a1bSJed Brown ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 856cab3a1bSJed Brown ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 86709ac0d2SJed Brown switch (dmdasnes->residuallocalimode) { 87709ac0d2SJed Brown case INSERT_VALUES: { 886cab3a1bSJed Brown ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 896cab3a1bSJed Brown CHKMEMQ; 906cab3a1bSJed Brown ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 916cab3a1bSJed Brown CHKMEMQ; 926cab3a1bSJed Brown ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 93709ac0d2SJed Brown } break; 94709ac0d2SJed Brown case ADD_VALUES: { 95709ac0d2SJed Brown Vec Floc; 96709ac0d2SJed Brown ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 97709ac0d2SJed Brown ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 98709ac0d2SJed Brown ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 99709ac0d2SJed Brown CHKMEMQ; 100709ac0d2SJed Brown ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 101709ac0d2SJed Brown CHKMEMQ; 102709ac0d2SJed Brown ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 103709ac0d2SJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 104709ac0d2SJed Brown ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 105709ac0d2SJed Brown ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 106709ac0d2SJed Brown ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 107709ac0d2SJed Brown } break; 108ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 109709ac0d2SJed Brown } 110709ac0d2SJed Brown ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 1116cab3a1bSJed Brown ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 1126cab3a1bSJed Brown PetscFunctionReturn(0); 1136cab3a1bSJed Brown } 1146cab3a1bSJed Brown 1156cab3a1bSJed Brown #undef __FUNCT__ 1162a4ee8f2SPeter Brune #define __FUNCT__ "SNESComputeObjective_DMDA" 1172a4ee8f2SPeter Brune static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx) 1182a4ee8f2SPeter Brune { 1192a4ee8f2SPeter Brune PetscErrorCode ierr; 1202a4ee8f2SPeter Brune DM dm; 121942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 1222a4ee8f2SPeter Brune DMDALocalInfo info; 1232a4ee8f2SPeter Brune Vec Xloc; 1242a4ee8f2SPeter Brune void *x; 1252a4ee8f2SPeter Brune 1262a4ee8f2SPeter Brune PetscFunctionBegin; 1272a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1282a4ee8f2SPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1292a4ee8f2SPeter Brune PetscValidPointer(ob,3); 130ce94432eSBarry Smith if (!dmdasnes->objectivelocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 1312a4ee8f2SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1322a4ee8f2SPeter Brune ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 1332a4ee8f2SPeter Brune ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1342a4ee8f2SPeter Brune ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1352a4ee8f2SPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 1362a4ee8f2SPeter Brune ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 1372a4ee8f2SPeter Brune CHKMEMQ; 1382a4ee8f2SPeter Brune ierr = (*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx);CHKERRQ(ierr); 1392a4ee8f2SPeter Brune CHKMEMQ; 1402a4ee8f2SPeter Brune ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 1412a4ee8f2SPeter Brune ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 1422a4ee8f2SPeter Brune PetscFunctionReturn(0); 1432a4ee8f2SPeter Brune } 1442a4ee8f2SPeter Brune 1456427ad5bSPeter Brune 1466427ad5bSPeter Brune #undef __FUNCT__ 1472961e153SJed Brown #define __FUNCT__ "SNESComputeJacobian_DMDA" 148d1e9a80fSBarry Smith static PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx) 1492961e153SJed Brown { 1502961e153SJed Brown PetscErrorCode ierr; 1512961e153SJed Brown DM dm; 152942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 1532961e153SJed Brown DMDALocalInfo info; 1542961e153SJed Brown Vec Xloc; 1552961e153SJed Brown void *x; 1562961e153SJed Brown 1572961e153SJed Brown PetscFunctionBegin; 158ce94432eSBarry Smith if (!dmdasnes->residuallocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 1592961e153SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1602961e153SJed Brown 1612961e153SJed Brown if (dmdasnes->jacobianlocal) { 1622961e153SJed Brown ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 1632961e153SJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1642961e153SJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1652961e153SJed Brown ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 1662961e153SJed Brown ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 1672961e153SJed Brown CHKMEMQ; 168d1e9a80fSBarry Smith ierr = (*dmdasnes->jacobianlocal)(&info,x,A,B,dmdasnes->jacobianlocalctx);CHKERRQ(ierr); 1692961e153SJed Brown CHKMEMQ; 1702961e153SJed Brown ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 1712961e153SJed Brown ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 1722961e153SJed Brown } else { 1732961e153SJed Brown MatFDColoring fdcoloring; 1742961e153SJed Brown ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr); 1752961e153SJed Brown if (!fdcoloring) { 1762961e153SJed Brown ISColoring coloring; 1772961e153SJed Brown 178b412c318SBarry Smith ierr = DMCreateColoring(dm,dm->coloringtype,&coloring);CHKERRQ(ierr); 17994ab13aaSBarry Smith ierr = MatFDColoringCreate(B,coloring,&fdcoloring);CHKERRQ(ierr); 1802961e153SJed Brown switch (dm->coloringtype) { 1812961e153SJed Brown case IS_COLORING_GLOBAL: 1822961e153SJed Brown ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 1832961e153SJed Brown break; 184ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 1852961e153SJed Brown } 1862961e153SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr); 1872961e153SJed Brown ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 18894ab13aaSBarry Smith ierr = MatFDColoringSetUp(B,coloring,fdcoloring);CHKERRQ(ierr); 189f86b9fbaSHong Zhang ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 1902961e153SJed Brown ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr); 1912961e153SJed Brown ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr); 1922961e153SJed Brown 1932961e153SJed Brown /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 1942961e153SJed Brown * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 1952961e153SJed Brown * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 1962961e153SJed Brown * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 197140e18c1SBarry Smith * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 1982961e153SJed Brown */ 1992961e153SJed Brown ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 2002961e153SJed Brown } 201d1e9a80fSBarry Smith ierr = MatFDColoringApply(B,fdcoloring,X,snes);CHKERRQ(ierr); 2022961e153SJed Brown } 2032961e153SJed Brown /* This will be redundant if the user called both, but it's too common to forget. */ 20494ab13aaSBarry Smith if (A != B) { 20594ab13aaSBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20694ab13aaSBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2072961e153SJed Brown } 2082961e153SJed Brown PetscFunctionReturn(0); 2092961e153SJed Brown } 2102961e153SJed Brown 2112961e153SJed Brown #undef __FUNCT__ 2126cab3a1bSJed Brown #define __FUNCT__ "DMDASNESSetFunctionLocal" 2136cab3a1bSJed Brown /*@C 2146cab3a1bSJed Brown DMDASNESSetFunctionLocal - set a local residual evaluation function 2156cab3a1bSJed Brown 2166cab3a1bSJed Brown Logically Collective 2176cab3a1bSJed Brown 2186cab3a1bSJed Brown Input Arguments: 2196cab3a1bSJed Brown + dm - DM to associate callback with 220e3c51ba2SJed Brown . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 2216cab3a1bSJed Brown . func - local residual evaluation 2226cab3a1bSJed Brown - ctx - optional context for local residual evaluation 2236cab3a1bSJed Brown 2246cab3a1bSJed Brown Calling sequence for func: 2256cab3a1bSJed Brown + info - DMDALocalInfo defining the subdomain to evaluate the residual on 2266cab3a1bSJed Brown . x - dimensional pointer to state at which to evaluate residual 2276cab3a1bSJed Brown . f - dimensional pointer to residual, write the residual here 2286cab3a1bSJed Brown - ctx - optional context passed above 2296cab3a1bSJed Brown 2306cab3a1bSJed Brown Level: beginner 2316cab3a1bSJed Brown 2326cab3a1bSJed Brown .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 2336cab3a1bSJed Brown @*/ 234709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx) 2356cab3a1bSJed Brown { 2366cab3a1bSJed Brown PetscErrorCode ierr; 237942e3340SBarry Smith DMSNES sdm; 238942e3340SBarry Smith DMSNES_DA *dmdasnes; 2396cab3a1bSJed Brown 2406cab3a1bSJed Brown PetscFunctionBegin; 2416cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 242942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 2436cab3a1bSJed Brown ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 2441aa26658SKarl Rupp 245709ac0d2SJed Brown dmdasnes->residuallocalimode = imode; 2466cab3a1bSJed Brown dmdasnes->residuallocal = func; 2476cab3a1bSJed Brown dmdasnes->residuallocalctx = ctx; 2481aa26658SKarl Rupp 2496cab3a1bSJed Brown ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 25022c6f798SBarry Smith if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 2512961e153SJed Brown ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 2522961e153SJed Brown } 2532961e153SJed Brown PetscFunctionReturn(0); 2542961e153SJed Brown } 2552961e153SJed Brown 2562961e153SJed Brown #undef __FUNCT__ 2572961e153SJed Brown #define __FUNCT__ "DMDASNESSetJacobianLocal" 2582961e153SJed Brown /*@C 259e79ed307SJed Brown DMDASNESSetJacobianLocal - set a local Jacobian evaluation function 2602961e153SJed Brown 2612961e153SJed Brown Logically Collective 2622961e153SJed Brown 2632961e153SJed Brown Input Arguments: 2642961e153SJed Brown + dm - DM to associate callback with 265*6b7eb5ceSMatthew G. Knepley . func - local Jacobian evaluation 266*6b7eb5ceSMatthew G. Knepley - ctx - optional context for local Jacobian evaluation 2672961e153SJed Brown 2682961e153SJed Brown Calling sequence for func: 269*6b7eb5ceSMatthew G. Knepley + info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at 270*6b7eb5ceSMatthew G. Knepley . x - dimensional pointer to state at which to evaluate Jacobian 271*6b7eb5ceSMatthew G. Knepley . J - Mat object for the Jacobian 272*6b7eb5ceSMatthew G. Knepley . M - Mat object for the Jacobian preconditioner matrix 2732961e153SJed Brown - ctx - optional context passed above 2742961e153SJed Brown 2752961e153SJed Brown Level: beginner 2762961e153SJed Brown 2772961e153SJed Brown .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 2782961e153SJed Brown @*/ 279d1e9a80fSBarry Smith PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx) 2802961e153SJed Brown { 2812961e153SJed Brown PetscErrorCode ierr; 282942e3340SBarry Smith DMSNES sdm; 283942e3340SBarry Smith DMSNES_DA *dmdasnes; 2842961e153SJed Brown 2852961e153SJed Brown PetscFunctionBegin; 2862961e153SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 287942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 2882961e153SJed Brown ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 2891aa26658SKarl Rupp 2902961e153SJed Brown dmdasnes->jacobianlocal = func; 2912961e153SJed Brown dmdasnes->jacobianlocalctx = ctx; 2921aa26658SKarl Rupp 2932961e153SJed Brown ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 2946cab3a1bSJed Brown PetscFunctionReturn(0); 2956cab3a1bSJed Brown } 2962a4ee8f2SPeter Brune 2972a4ee8f2SPeter Brune 2982a4ee8f2SPeter Brune #undef __FUNCT__ 2992a4ee8f2SPeter Brune #define __FUNCT__ "DMDASNESSetObjectiveLocal" 3002a4ee8f2SPeter Brune /*@C 3012a4ee8f2SPeter Brune DMDASNESSetObjectiveLocal - set a local residual evaluation function 3022a4ee8f2SPeter Brune 3032a4ee8f2SPeter Brune Logically Collective 3042a4ee8f2SPeter Brune 3052a4ee8f2SPeter Brune Input Arguments: 3062a4ee8f2SPeter Brune + dm - DM to associate callback with 3072a4ee8f2SPeter Brune . func - local objective evaluation 3082a4ee8f2SPeter Brune - ctx - optional context for local residual evaluation 3092a4ee8f2SPeter Brune 3102a4ee8f2SPeter Brune Calling sequence for func: 3112a4ee8f2SPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 3122a4ee8f2SPeter Brune . x - dimensional pointer to state at which to evaluate residual 3132a4ee8f2SPeter Brune . ob - eventual objective value 3142a4ee8f2SPeter Brune - ctx - optional context passed above 3152a4ee8f2SPeter Brune 3162a4ee8f2SPeter Brune Level: beginner 3172a4ee8f2SPeter Brune 3182a4ee8f2SPeter Brune .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 3192a4ee8f2SPeter Brune @*/ 3202a4ee8f2SPeter Brune PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx) 3212a4ee8f2SPeter Brune { 3222a4ee8f2SPeter Brune PetscErrorCode ierr; 323942e3340SBarry Smith DMSNES sdm; 324942e3340SBarry Smith DMSNES_DA *dmdasnes; 3252a4ee8f2SPeter Brune 3262a4ee8f2SPeter Brune PetscFunctionBegin; 3272a4ee8f2SPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 328942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 3292a4ee8f2SPeter Brune ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 3301aa26658SKarl Rupp 3312a4ee8f2SPeter Brune dmdasnes->objectivelocal = func; 3322a4ee8f2SPeter Brune dmdasnes->objectivelocalctx = ctx; 3331aa26658SKarl Rupp 3342a4ee8f2SPeter Brune ierr = DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes);CHKERRQ(ierr); 3352a4ee8f2SPeter Brune PetscFunctionReturn(0); 3362a4ee8f2SPeter Brune } 337dcad5f1cSBarry Smith 338dcad5f1cSBarry Smith #undef __FUNCT__ 339dcad5f1cSBarry Smith #define __FUNCT__ "SNESComputePicard_DMDA" 340dcad5f1cSBarry Smith static PetscErrorCode SNESComputePicard_DMDA(SNES snes,Vec X,Vec F,void *ctx) 341dcad5f1cSBarry Smith { 342dcad5f1cSBarry Smith PetscErrorCode ierr; 343dcad5f1cSBarry Smith DM dm; 344942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 345dcad5f1cSBarry Smith DMDALocalInfo info; 346dcad5f1cSBarry Smith Vec Xloc; 347dcad5f1cSBarry Smith void *x,*f; 348dcad5f1cSBarry Smith 349dcad5f1cSBarry Smith PetscFunctionBegin; 350dcad5f1cSBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 351dcad5f1cSBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,2); 352dcad5f1cSBarry Smith PetscValidHeaderSpecific(F,VEC_CLASSID,3); 353ce94432eSBarry Smith if (!dmdasnes->rhsplocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 354dcad5f1cSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 355dcad5f1cSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 356dcad5f1cSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 357dcad5f1cSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 358dcad5f1cSBarry Smith ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 359dcad5f1cSBarry Smith ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 360dcad5f1cSBarry Smith switch (dmdasnes->residuallocalimode) { 361dcad5f1cSBarry Smith case INSERT_VALUES: { 362dcad5f1cSBarry Smith ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 363dcad5f1cSBarry Smith CHKMEMQ; 364dcad5f1cSBarry Smith ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr); 365dcad5f1cSBarry Smith CHKMEMQ; 366dcad5f1cSBarry Smith ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 367dcad5f1cSBarry Smith } break; 368dcad5f1cSBarry Smith case ADD_VALUES: { 369dcad5f1cSBarry Smith Vec Floc; 370dcad5f1cSBarry Smith ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 371dcad5f1cSBarry Smith ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 372dcad5f1cSBarry Smith ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 373dcad5f1cSBarry Smith CHKMEMQ; 374dcad5f1cSBarry Smith ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr); 375dcad5f1cSBarry Smith CHKMEMQ; 376dcad5f1cSBarry Smith ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 377dcad5f1cSBarry Smith ierr = VecZeroEntries(F);CHKERRQ(ierr); 378dcad5f1cSBarry Smith ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 379dcad5f1cSBarry Smith ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 380dcad5f1cSBarry Smith ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 381dcad5f1cSBarry Smith } break; 382ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 383dcad5f1cSBarry Smith } 384dcad5f1cSBarry Smith ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 385dcad5f1cSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 386dcad5f1cSBarry Smith PetscFunctionReturn(0); 387dcad5f1cSBarry Smith } 388dcad5f1cSBarry Smith 389dcad5f1cSBarry Smith #undef __FUNCT__ 390dcad5f1cSBarry Smith #define __FUNCT__ "SNESComputePicardJacobian_DMDA" 391d1e9a80fSBarry Smith static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx) 392dcad5f1cSBarry Smith { 393dcad5f1cSBarry Smith PetscErrorCode ierr; 394dcad5f1cSBarry Smith DM dm; 395942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 396dcad5f1cSBarry Smith DMDALocalInfo info; 397dcad5f1cSBarry Smith Vec Xloc; 398dcad5f1cSBarry Smith void *x; 399dcad5f1cSBarry Smith 400dcad5f1cSBarry Smith PetscFunctionBegin; 401ce94432eSBarry Smith if (!dmdasnes->jacobianplocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 402dcad5f1cSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 403dcad5f1cSBarry Smith 404dcad5f1cSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 405dcad5f1cSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 406dcad5f1cSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 407dcad5f1cSBarry Smith ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 408dcad5f1cSBarry Smith ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 409dcad5f1cSBarry Smith CHKMEMQ; 410d1e9a80fSBarry Smith ierr = (*dmdasnes->jacobianplocal)(&info,x,A,B,dmdasnes->picardlocalctx);CHKERRQ(ierr); 411dcad5f1cSBarry Smith CHKMEMQ; 412dcad5f1cSBarry Smith ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 413dcad5f1cSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 41494ab13aaSBarry Smith if (A != B) { 41594ab13aaSBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 41694ab13aaSBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 417dcad5f1cSBarry Smith } 418dcad5f1cSBarry Smith PetscFunctionReturn(0); 419dcad5f1cSBarry Smith } 420dcad5f1cSBarry Smith 421dcad5f1cSBarry Smith #undef __FUNCT__ 422dcad5f1cSBarry Smith #define __FUNCT__ "DMDASNESSetPicardLocal" 423dcad5f1cSBarry Smith /*@C 424dcad5f1cSBarry Smith DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration 425dcad5f1cSBarry Smith 426dcad5f1cSBarry Smith Logically Collective 427dcad5f1cSBarry Smith 428dcad5f1cSBarry Smith Input Arguments: 429dcad5f1cSBarry Smith + dm - DM to associate callback with 430dcad5f1cSBarry Smith . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 431dcad5f1cSBarry Smith . func - local residual evaluation 432dcad5f1cSBarry Smith - ctx - optional context for local residual evaluation 433dcad5f1cSBarry Smith 434dcad5f1cSBarry Smith Calling sequence for func: 435dcad5f1cSBarry Smith + info - DMDALocalInfo defining the subdomain to evaluate the residual on 436dcad5f1cSBarry Smith . x - dimensional pointer to state at which to evaluate residual 437dcad5f1cSBarry Smith . f - dimensional pointer to residual, write the residual here 438dcad5f1cSBarry Smith - ctx - optional context passed above 439dcad5f1cSBarry Smith 440dcad5f1cSBarry Smith Notes: The user must use 441dcad5f1cSBarry Smith extern PetscErrorCode SNESPicardComputeFunction(SNES,Vec,Vec,void*); 44294ab13aaSBarry Smith extern PetscErrorCode SNESPicardComputeJacobian(SNES,Vec,Mat,Mat,MatStructure*,void*); 4430298fd71SBarry Smith ierr = SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user);CHKERRQ(ierr); 444dcad5f1cSBarry Smith in their code before calling this routine. 445dcad5f1cSBarry Smith 446dcad5f1cSBarry Smith 447dcad5f1cSBarry Smith Level: beginner 448dcad5f1cSBarry Smith 449dcad5f1cSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 450dcad5f1cSBarry Smith @*/ 451dcad5f1cSBarry Smith PetscErrorCode DMDASNESSetPicardLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*), 452d1e9a80fSBarry Smith PetscErrorCode (*jac)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx) 453dcad5f1cSBarry Smith { 454dcad5f1cSBarry Smith PetscErrorCode ierr; 455942e3340SBarry Smith DMSNES sdm; 456942e3340SBarry Smith DMSNES_DA *dmdasnes; 457dcad5f1cSBarry Smith 458dcad5f1cSBarry Smith PetscFunctionBegin; 459dcad5f1cSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 460942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 461dcad5f1cSBarry Smith ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 4621aa26658SKarl Rupp 463dcad5f1cSBarry Smith dmdasnes->residuallocalimode = imode; 464dcad5f1cSBarry Smith dmdasnes->rhsplocal = func; 465dcad5f1cSBarry Smith dmdasnes->jacobianplocal = jac; 466dcad5f1cSBarry Smith dmdasnes->picardlocalctx = ctx; 4671aa26658SKarl Rupp 468dcad5f1cSBarry Smith ierr = DMSNESSetPicard(dm,SNESComputePicard_DMDA,SNESComputePicardJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 469dcad5f1cSBarry Smith PetscFunctionReturn(0); 470dcad5f1cSBarry Smith } 471