16cab3a1bSJed Brown #include <petscdmda.h> /*I "petscdmda.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 3af0996ceSBarry Smith #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 21942e3340SBarry Smith static PetscErrorCode DMSNESDestroy_DMDA(DMSNES sdm) 226cab3a1bSJed Brown { 236cab3a1bSJed Brown PetscErrorCode ierr; 246cab3a1bSJed Brown 256cab3a1bSJed Brown PetscFunctionBegin; 266cab3a1bSJed Brown ierr = PetscFree(sdm->data);CHKERRQ(ierr); 276cab3a1bSJed Brown PetscFunctionReturn(0); 286cab3a1bSJed Brown } 296cab3a1bSJed Brown 3022c6f798SBarry Smith static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm,DMSNES sdm) 31e3c0b8a2SPeter Brune { 32e3c0b8a2SPeter Brune PetscErrorCode ierr; 33e3c0b8a2SPeter Brune 34e3c0b8a2SPeter Brune PetscFunctionBegin; 35b00a9115SJed Brown ierr = PetscNewLog(sdm,(DMSNES_DA**)&sdm->data);CHKERRQ(ierr); 36dcad5f1cSBarry Smith if (oldsdm->data) { 37942e3340SBarry Smith ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_DA));CHKERRQ(ierr); 38dcad5f1cSBarry Smith } 39e3c0b8a2SPeter Brune PetscFunctionReturn(0); 40e3c0b8a2SPeter Brune } 41e3c0b8a2SPeter Brune 42e3c0b8a2SPeter Brune 43942e3340SBarry Smith static PetscErrorCode DMDASNESGetContext(DM dm,DMSNES sdm,DMSNES_DA **dmdasnes) 446cab3a1bSJed Brown { 456cab3a1bSJed Brown PetscErrorCode ierr; 466cab3a1bSJed Brown 476cab3a1bSJed Brown PetscFunctionBegin; 480298fd71SBarry Smith *dmdasnes = NULL; 496cab3a1bSJed Brown if (!sdm->data) { 50b00a9115SJed Brown ierr = PetscNewLog(dm,(DMSNES_DA**)&sdm->data);CHKERRQ(ierr); 5122c6f798SBarry Smith sdm->ops->destroy = DMSNESDestroy_DMDA; 5222c6f798SBarry Smith sdm->ops->duplicate = DMSNESDuplicate_DMDA; 536cab3a1bSJed Brown } 54942e3340SBarry Smith *dmdasnes = (DMSNES_DA*)sdm->data; 556cab3a1bSJed Brown PetscFunctionReturn(0); 566cab3a1bSJed Brown } 576cab3a1bSJed Brown 586cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx) 596cab3a1bSJed Brown { 606cab3a1bSJed Brown PetscErrorCode ierr; 616cab3a1bSJed Brown DM dm; 62942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 636cab3a1bSJed Brown DMDALocalInfo info; 646cab3a1bSJed Brown Vec Xloc; 656cab3a1bSJed Brown void *x,*f; 666cab3a1bSJed Brown 676cab3a1bSJed Brown PetscFunctionBegin; 682961e153SJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 692961e153SJed Brown PetscValidHeaderSpecific(X,VEC_CLASSID,2); 702961e153SJed Brown PetscValidHeaderSpecific(F,VEC_CLASSID,3); 71ce94432eSBarry Smith if (!dmdasnes->residuallocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 726cab3a1bSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 736cab3a1bSJed Brown ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 746cab3a1bSJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 756cab3a1bSJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 766cab3a1bSJed Brown ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 776cab3a1bSJed Brown ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 78709ac0d2SJed Brown switch (dmdasnes->residuallocalimode) { 79709ac0d2SJed Brown case INSERT_VALUES: { 806cab3a1bSJed Brown ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 8162cd874fSBarry Smith ierr = PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0);CHKERRQ(ierr); 826cab3a1bSJed Brown CHKMEMQ; 836cab3a1bSJed Brown ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 846cab3a1bSJed Brown CHKMEMQ; 8562cd874fSBarry Smith ierr = PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0);CHKERRQ(ierr); 866cab3a1bSJed Brown ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 87709ac0d2SJed Brown } break; 88709ac0d2SJed Brown case ADD_VALUES: { 89709ac0d2SJed Brown Vec Floc; 90709ac0d2SJed Brown ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 91709ac0d2SJed Brown ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 92709ac0d2SJed Brown ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 9362cd874fSBarry Smith ierr = PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0);CHKERRQ(ierr); 94709ac0d2SJed Brown CHKMEMQ; 95709ac0d2SJed Brown ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 96709ac0d2SJed Brown CHKMEMQ; 9762cd874fSBarry Smith ierr = PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0);CHKERRQ(ierr); 98709ac0d2SJed Brown ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 99709ac0d2SJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 100709ac0d2SJed Brown ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 101709ac0d2SJed Brown ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 102709ac0d2SJed Brown ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 103709ac0d2SJed Brown } break; 104ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 105709ac0d2SJed Brown } 106709ac0d2SJed Brown ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 1076cab3a1bSJed Brown ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 108422a814eSBarry Smith if (snes->domainerror) { 109422a814eSBarry Smith ierr = VecSetInf(F);CHKERRQ(ierr); 110422a814eSBarry Smith } 1116cab3a1bSJed Brown PetscFunctionReturn(0); 1126cab3a1bSJed Brown } 1136cab3a1bSJed Brown 1142a4ee8f2SPeter Brune static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx) 1152a4ee8f2SPeter Brune { 1162a4ee8f2SPeter Brune PetscErrorCode ierr; 1172a4ee8f2SPeter Brune DM dm; 118942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 1192a4ee8f2SPeter Brune DMDALocalInfo info; 1202a4ee8f2SPeter Brune Vec Xloc; 1212a4ee8f2SPeter Brune void *x; 1222a4ee8f2SPeter Brune 1232a4ee8f2SPeter Brune PetscFunctionBegin; 1242a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1252a4ee8f2SPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1262a4ee8f2SPeter Brune PetscValidPointer(ob,3); 127ce94432eSBarry Smith if (!dmdasnes->objectivelocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 1282a4ee8f2SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1292a4ee8f2SPeter Brune ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 1302a4ee8f2SPeter Brune ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1312a4ee8f2SPeter Brune ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1322a4ee8f2SPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 1332a4ee8f2SPeter Brune ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 1342a4ee8f2SPeter Brune CHKMEMQ; 1352a4ee8f2SPeter Brune ierr = (*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx);CHKERRQ(ierr); 1362a4ee8f2SPeter Brune CHKMEMQ; 1372a4ee8f2SPeter Brune ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 1382a4ee8f2SPeter Brune ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 1392a4ee8f2SPeter Brune PetscFunctionReturn(0); 1402a4ee8f2SPeter Brune } 1412a4ee8f2SPeter Brune 1426427ad5bSPeter Brune 143cc2e6a90SBarry Smith /* Routine is called by example, hence must be labeled PETSC_EXTERN */ 144cc2e6a90SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx) 1452961e153SJed Brown { 1462961e153SJed Brown PetscErrorCode ierr; 1472961e153SJed Brown DM dm; 148942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 1492961e153SJed Brown DMDALocalInfo info; 1502961e153SJed Brown Vec Xloc; 1512961e153SJed Brown void *x; 1522961e153SJed Brown 1532961e153SJed Brown PetscFunctionBegin; 154ce94432eSBarry Smith if (!dmdasnes->residuallocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 1552961e153SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1562961e153SJed Brown 1572961e153SJed Brown if (dmdasnes->jacobianlocal) { 1582961e153SJed Brown ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 1592961e153SJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1602961e153SJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1612961e153SJed Brown ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 1622961e153SJed Brown ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 1632961e153SJed Brown CHKMEMQ; 164d1e9a80fSBarry Smith ierr = (*dmdasnes->jacobianlocal)(&info,x,A,B,dmdasnes->jacobianlocalctx);CHKERRQ(ierr); 1652961e153SJed Brown CHKMEMQ; 1662961e153SJed Brown ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 1672961e153SJed Brown ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 1682961e153SJed Brown } else { 1692961e153SJed Brown MatFDColoring fdcoloring; 1702961e153SJed Brown ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr); 1712961e153SJed Brown if (!fdcoloring) { 1722961e153SJed Brown ISColoring coloring; 1732961e153SJed Brown 174b412c318SBarry Smith ierr = DMCreateColoring(dm,dm->coloringtype,&coloring);CHKERRQ(ierr); 17594ab13aaSBarry Smith ierr = MatFDColoringCreate(B,coloring,&fdcoloring);CHKERRQ(ierr); 1762961e153SJed Brown switch (dm->coloringtype) { 1772961e153SJed Brown case IS_COLORING_GLOBAL: 1782961e153SJed Brown ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 1792961e153SJed Brown break; 180ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 1812961e153SJed Brown } 1822961e153SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr); 1832961e153SJed Brown ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 18494ab13aaSBarry Smith ierr = MatFDColoringSetUp(B,coloring,fdcoloring);CHKERRQ(ierr); 185f86b9fbaSHong Zhang ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 1862961e153SJed Brown ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr); 1872961e153SJed Brown ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr); 1882961e153SJed Brown 1892961e153SJed Brown /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 1902961e153SJed Brown * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 1912961e153SJed Brown * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 1922961e153SJed Brown * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 193140e18c1SBarry Smith * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 1942961e153SJed Brown */ 1952961e153SJed Brown ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 1962961e153SJed Brown } 197d1e9a80fSBarry Smith ierr = MatFDColoringApply(B,fdcoloring,X,snes);CHKERRQ(ierr); 1982961e153SJed Brown } 1992961e153SJed Brown /* This will be redundant if the user called both, but it's too common to forget. */ 20094ab13aaSBarry Smith if (A != B) { 20194ab13aaSBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20294ab13aaSBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2032961e153SJed Brown } 2042961e153SJed Brown PetscFunctionReturn(0); 2052961e153SJed Brown } 2062961e153SJed Brown 2076cab3a1bSJed Brown /*@C 2086cab3a1bSJed Brown DMDASNESSetFunctionLocal - set a local residual evaluation function 2096cab3a1bSJed Brown 2106cab3a1bSJed Brown Logically Collective 2116cab3a1bSJed Brown 2126cab3a1bSJed Brown Input Arguments: 2136cab3a1bSJed Brown + dm - DM to associate callback with 214e3c51ba2SJed Brown . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 2156cab3a1bSJed Brown . func - local residual evaluation 2166cab3a1bSJed Brown - ctx - optional context for local residual evaluation 2176cab3a1bSJed Brown 2180038884cSBarry Smith Calling sequence: 2190038884cSBarry Smith For PetscErrorCode (*func)(DMDALocalInfo *info,void *x, void *f, void *ctx), 2206cab3a1bSJed Brown + info - DMDALocalInfo defining the subdomain to evaluate the residual on 2210038884cSBarry Smith . x - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x) 2220038884cSBarry Smith . f - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f) 2236cab3a1bSJed Brown - ctx - optional context passed above 2246cab3a1bSJed Brown 2256cab3a1bSJed Brown Level: beginner 2266cab3a1bSJed Brown 2270038884cSBarry Smith .seealso: DMDASNESSetJacobianLocal(), DMSNESSetFunction(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 2286cab3a1bSJed Brown @*/ 229709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx) 2306cab3a1bSJed Brown { 2316cab3a1bSJed Brown PetscErrorCode ierr; 232942e3340SBarry Smith DMSNES sdm; 233942e3340SBarry Smith DMSNES_DA *dmdasnes; 2346cab3a1bSJed Brown 2356cab3a1bSJed Brown PetscFunctionBegin; 2366cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 237942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 2386cab3a1bSJed Brown ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 2391aa26658SKarl Rupp 240709ac0d2SJed Brown dmdasnes->residuallocalimode = imode; 2416cab3a1bSJed Brown dmdasnes->residuallocal = func; 2426cab3a1bSJed Brown dmdasnes->residuallocalctx = ctx; 2431aa26658SKarl Rupp 2446cab3a1bSJed Brown ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 24522c6f798SBarry Smith if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 2462961e153SJed Brown ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 2472961e153SJed Brown } 2482961e153SJed Brown PetscFunctionReturn(0); 2492961e153SJed Brown } 2502961e153SJed Brown 2512961e153SJed Brown /*@C 252e79ed307SJed Brown DMDASNESSetJacobianLocal - set a local Jacobian evaluation function 2532961e153SJed Brown 2542961e153SJed Brown Logically Collective 2552961e153SJed Brown 2562961e153SJed Brown Input Arguments: 2572961e153SJed Brown + dm - DM to associate callback with 2586b7eb5ceSMatthew G. Knepley . func - local Jacobian evaluation 2596b7eb5ceSMatthew G. Knepley - ctx - optional context for local Jacobian evaluation 2602961e153SJed Brown 2610038884cSBarry Smith Calling sequence: 2620038884cSBarry Smith For PetscErrorCode (*func)(DMDALocalInfo *info,void *x,Mat J,Mat M,void *ctx), 2636b7eb5ceSMatthew G. Knepley + info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at 2640038884cSBarry Smith . x - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x) 2656b7eb5ceSMatthew G. Knepley . J - Mat object for the Jacobian 2666b7eb5ceSMatthew G. Knepley . M - Mat object for the Jacobian preconditioner matrix 2672961e153SJed Brown - ctx - optional context passed above 2682961e153SJed Brown 2692961e153SJed Brown Level: beginner 2702961e153SJed Brown 2710038884cSBarry Smith .seealso: DMDASNESSetFunctionLocal(), DMSNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 2722961e153SJed Brown @*/ 273d1e9a80fSBarry Smith PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx) 2742961e153SJed Brown { 2752961e153SJed Brown PetscErrorCode ierr; 276942e3340SBarry Smith DMSNES sdm; 277942e3340SBarry Smith DMSNES_DA *dmdasnes; 2782961e153SJed Brown 2792961e153SJed Brown PetscFunctionBegin; 2802961e153SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 281942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 2822961e153SJed Brown ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 2831aa26658SKarl Rupp 2842961e153SJed Brown dmdasnes->jacobianlocal = func; 2852961e153SJed Brown dmdasnes->jacobianlocalctx = ctx; 2861aa26658SKarl Rupp 2872961e153SJed Brown ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 2886cab3a1bSJed Brown PetscFunctionReturn(0); 2896cab3a1bSJed Brown } 2902a4ee8f2SPeter Brune 2912a4ee8f2SPeter Brune 2922a4ee8f2SPeter Brune /*@C 2932a4ee8f2SPeter Brune DMDASNESSetObjectiveLocal - set a local residual evaluation function 2942a4ee8f2SPeter Brune 2952a4ee8f2SPeter Brune Logically Collective 2962a4ee8f2SPeter Brune 2972a4ee8f2SPeter Brune Input Arguments: 2982a4ee8f2SPeter Brune + dm - DM to associate callback with 2992a4ee8f2SPeter Brune . func - local objective evaluation 3002a4ee8f2SPeter Brune - ctx - optional context for local residual evaluation 3012a4ee8f2SPeter Brune 3022a4ee8f2SPeter Brune Calling sequence for func: 3032a4ee8f2SPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 3042a4ee8f2SPeter Brune . x - dimensional pointer to state at which to evaluate residual 3052a4ee8f2SPeter Brune . ob - eventual objective value 3062a4ee8f2SPeter Brune - ctx - optional context passed above 3072a4ee8f2SPeter Brune 3082a4ee8f2SPeter Brune Level: beginner 3092a4ee8f2SPeter Brune 3100038884cSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 3112a4ee8f2SPeter Brune @*/ 3122a4ee8f2SPeter Brune PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx) 3132a4ee8f2SPeter Brune { 3142a4ee8f2SPeter Brune PetscErrorCode ierr; 315942e3340SBarry Smith DMSNES sdm; 316942e3340SBarry Smith DMSNES_DA *dmdasnes; 3172a4ee8f2SPeter Brune 3182a4ee8f2SPeter Brune PetscFunctionBegin; 3192a4ee8f2SPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 320942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 3212a4ee8f2SPeter Brune ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 3221aa26658SKarl Rupp 3232a4ee8f2SPeter Brune dmdasnes->objectivelocal = func; 3242a4ee8f2SPeter Brune dmdasnes->objectivelocalctx = ctx; 3251aa26658SKarl Rupp 3262a4ee8f2SPeter Brune ierr = DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes);CHKERRQ(ierr); 3272a4ee8f2SPeter Brune PetscFunctionReturn(0); 3282a4ee8f2SPeter Brune } 329dcad5f1cSBarry Smith 330dcad5f1cSBarry Smith static PetscErrorCode SNESComputePicard_DMDA(SNES snes,Vec X,Vec F,void *ctx) 331dcad5f1cSBarry Smith { 332dcad5f1cSBarry Smith PetscErrorCode ierr; 333dcad5f1cSBarry Smith DM dm; 334942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 335dcad5f1cSBarry Smith DMDALocalInfo info; 336dcad5f1cSBarry Smith Vec Xloc; 337dcad5f1cSBarry Smith void *x,*f; 338dcad5f1cSBarry Smith 339dcad5f1cSBarry Smith PetscFunctionBegin; 340dcad5f1cSBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 341dcad5f1cSBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,2); 342dcad5f1cSBarry Smith PetscValidHeaderSpecific(F,VEC_CLASSID,3); 343ce94432eSBarry Smith if (!dmdasnes->rhsplocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 344dcad5f1cSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 345dcad5f1cSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 346dcad5f1cSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 347dcad5f1cSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 348dcad5f1cSBarry Smith ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 349dcad5f1cSBarry Smith ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 350dcad5f1cSBarry Smith switch (dmdasnes->residuallocalimode) { 351dcad5f1cSBarry Smith case INSERT_VALUES: { 352dcad5f1cSBarry Smith ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 353dcad5f1cSBarry Smith CHKMEMQ; 354dcad5f1cSBarry Smith ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr); 355dcad5f1cSBarry Smith CHKMEMQ; 356dcad5f1cSBarry Smith ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 357dcad5f1cSBarry Smith } break; 358dcad5f1cSBarry Smith case ADD_VALUES: { 359dcad5f1cSBarry Smith Vec Floc; 360dcad5f1cSBarry Smith ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 361dcad5f1cSBarry Smith ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 362dcad5f1cSBarry Smith ierr = DMDAVecGetArray(dm,Floc,&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,Floc,&f);CHKERRQ(ierr); 367dcad5f1cSBarry Smith ierr = VecZeroEntries(F);CHKERRQ(ierr); 368dcad5f1cSBarry Smith ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 369dcad5f1cSBarry Smith ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 370dcad5f1cSBarry Smith ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 371dcad5f1cSBarry Smith } break; 372ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 373dcad5f1cSBarry Smith } 374dcad5f1cSBarry Smith ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 375dcad5f1cSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 376dcad5f1cSBarry Smith PetscFunctionReturn(0); 377dcad5f1cSBarry Smith } 378dcad5f1cSBarry Smith 379d1e9a80fSBarry Smith static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx) 380dcad5f1cSBarry Smith { 381dcad5f1cSBarry Smith PetscErrorCode ierr; 382dcad5f1cSBarry Smith DM dm; 383942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 384dcad5f1cSBarry Smith DMDALocalInfo info; 385dcad5f1cSBarry Smith Vec Xloc; 386dcad5f1cSBarry Smith void *x; 387dcad5f1cSBarry Smith 388dcad5f1cSBarry Smith PetscFunctionBegin; 389ce94432eSBarry Smith if (!dmdasnes->jacobianplocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 390dcad5f1cSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 391dcad5f1cSBarry Smith 392dcad5f1cSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 393dcad5f1cSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 394dcad5f1cSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 395dcad5f1cSBarry Smith ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 396dcad5f1cSBarry Smith ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 397dcad5f1cSBarry Smith CHKMEMQ; 398d1e9a80fSBarry Smith ierr = (*dmdasnes->jacobianplocal)(&info,x,A,B,dmdasnes->picardlocalctx);CHKERRQ(ierr); 399dcad5f1cSBarry Smith CHKMEMQ; 400dcad5f1cSBarry Smith ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 401dcad5f1cSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 40294ab13aaSBarry Smith if (A != B) { 40394ab13aaSBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 40494ab13aaSBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 405dcad5f1cSBarry Smith } 406dcad5f1cSBarry Smith PetscFunctionReturn(0); 407dcad5f1cSBarry Smith } 408dcad5f1cSBarry Smith 409dcad5f1cSBarry Smith /*@C 410dcad5f1cSBarry Smith DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration 411dcad5f1cSBarry Smith 412dcad5f1cSBarry Smith Logically Collective 413dcad5f1cSBarry Smith 414dcad5f1cSBarry Smith Input Arguments: 415dcad5f1cSBarry Smith + dm - DM to associate callback with 416dcad5f1cSBarry Smith . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 417dcad5f1cSBarry Smith . func - local residual evaluation 418dcad5f1cSBarry Smith - ctx - optional context for local residual evaluation 419dcad5f1cSBarry Smith 420dcad5f1cSBarry Smith Calling sequence for func: 421dcad5f1cSBarry Smith + info - DMDALocalInfo defining the subdomain to evaluate the residual on 422dcad5f1cSBarry Smith . x - dimensional pointer to state at which to evaluate residual 423dcad5f1cSBarry Smith . f - dimensional pointer to residual, write the residual here 424dcad5f1cSBarry Smith - ctx - optional context passed above 425dcad5f1cSBarry Smith 42695452b02SPatrick Sanan Notes: 42795452b02SPatrick Sanan The user must use 4280298fd71SBarry Smith ierr = SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user);CHKERRQ(ierr); 429dcad5f1cSBarry Smith in their code before calling this routine. 430dcad5f1cSBarry Smith 431dcad5f1cSBarry Smith 432dcad5f1cSBarry Smith Level: beginner 433dcad5f1cSBarry Smith 434dcad5f1cSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 435dcad5f1cSBarry Smith @*/ 436dcad5f1cSBarry Smith PetscErrorCode DMDASNESSetPicardLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*), 437d1e9a80fSBarry Smith PetscErrorCode (*jac)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx) 438dcad5f1cSBarry Smith { 439dcad5f1cSBarry Smith PetscErrorCode ierr; 440942e3340SBarry Smith DMSNES sdm; 441942e3340SBarry Smith DMSNES_DA *dmdasnes; 442dcad5f1cSBarry Smith 443dcad5f1cSBarry Smith PetscFunctionBegin; 444dcad5f1cSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 445942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 446dcad5f1cSBarry Smith ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 4471aa26658SKarl Rupp 448dcad5f1cSBarry Smith dmdasnes->residuallocalimode = imode; 449dcad5f1cSBarry Smith dmdasnes->rhsplocal = func; 450dcad5f1cSBarry Smith dmdasnes->jacobianplocal = jac; 451dcad5f1cSBarry Smith dmdasnes->picardlocalctx = ctx; 4521aa26658SKarl Rupp 453dcad5f1cSBarry Smith ierr = DMSNESSetPicard(dm,SNESComputePicard_DMDA,SNESComputePicardJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 454*bbc1464cSBarry Smith ierr = DMSNESSetMFFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 455dcad5f1cSBarry Smith PetscFunctionReturn(0); 456dcad5f1cSBarry Smith } 457