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*); 86cab3a1bSJed Brown PetscErrorCode (*jacobianlocal)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,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*); 17dcad5f1cSBarry Smith PetscErrorCode (*jacobianplocal)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*); 18dcad5f1cSBarry Smith void *picardlocalctx; 196cab3a1bSJed Brown } DM_DA_SNES; 206cab3a1bSJed Brown 216cab3a1bSJed Brown #undef __FUNCT__ 226cab3a1bSJed Brown #define __FUNCT__ "SNESDMDestroy_DMDA" 236cab3a1bSJed Brown static PetscErrorCode SNESDMDestroy_DMDA(SNESDM 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__ 33e3c0b8a2SPeter Brune #define __FUNCT__ "SNESDMDuplicate_DMDA" 34e3c0b8a2SPeter Brune static PetscErrorCode SNESDMDuplicate_DMDA(SNESDM oldsdm,DM dm) 35e3c0b8a2SPeter Brune { 36e3c0b8a2SPeter Brune PetscErrorCode ierr; 37e3c0b8a2SPeter Brune SNESDM sdm; 38e3c0b8a2SPeter Brune 39e3c0b8a2SPeter Brune PetscFunctionBegin; 40e3c0b8a2SPeter Brune ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 41e3c0b8a2SPeter Brune ierr = PetscNewLog(dm,DM_DA_SNES,&sdm->data);CHKERRQ(ierr); 42dcad5f1cSBarry Smith if (oldsdm->data) { 43dcad5f1cSBarry Smith ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DM_DA_SNES));CHKERRQ(ierr); 44dcad5f1cSBarry Smith } 45e3c0b8a2SPeter Brune PetscFunctionReturn(0); 46e3c0b8a2SPeter Brune } 47e3c0b8a2SPeter Brune 48e3c0b8a2SPeter Brune 49e3c0b8a2SPeter Brune #undef __FUNCT__ 506cab3a1bSJed Brown #define __FUNCT__ "DMDASNESGetContext" 516cab3a1bSJed Brown static PetscErrorCode DMDASNESGetContext(DM dm,SNESDM sdm,DM_DA_SNES **dmdasnes) 526cab3a1bSJed Brown { 536cab3a1bSJed Brown PetscErrorCode ierr; 546cab3a1bSJed Brown 556cab3a1bSJed Brown PetscFunctionBegin; 56b8d0cc33SJed Brown *dmdasnes = PETSC_NULL; 576cab3a1bSJed Brown if (!sdm->data) { 586cab3a1bSJed Brown ierr = PetscNewLog(dm,DM_DA_SNES,&sdm->data);CHKERRQ(ierr); 596cab3a1bSJed Brown sdm->destroy = SNESDMDestroy_DMDA; 60e3c0b8a2SPeter Brune sdm->duplicate = SNESDMDuplicate_DMDA; 616cab3a1bSJed Brown } 626cab3a1bSJed Brown *dmdasnes = (DM_DA_SNES*)sdm->data; 636cab3a1bSJed Brown PetscFunctionReturn(0); 646cab3a1bSJed Brown } 656cab3a1bSJed Brown 666cab3a1bSJed Brown #undef __FUNCT__ 676cab3a1bSJed Brown #define __FUNCT__ "SNESComputeFunction_DMDA" 686cab3a1bSJed Brown /* 696cab3a1bSJed Brown This function should eventually replace: 706cab3a1bSJed Brown DMDAComputeFunction() and DMDAComputeFunction1() 716cab3a1bSJed Brown */ 726cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx) 736cab3a1bSJed Brown { 746cab3a1bSJed Brown PetscErrorCode ierr; 756cab3a1bSJed Brown DM dm; 766cab3a1bSJed Brown DM_DA_SNES *dmdasnes = (DM_DA_SNES*)ctx; 776cab3a1bSJed Brown DMDALocalInfo info; 786cab3a1bSJed Brown Vec Xloc; 796cab3a1bSJed Brown void *x,*f; 806cab3a1bSJed Brown 816cab3a1bSJed Brown PetscFunctionBegin; 822961e153SJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 832961e153SJed Brown PetscValidHeaderSpecific(X,VEC_CLASSID,2); 842961e153SJed Brown PetscValidHeaderSpecific(F,VEC_CLASSID,3); 856cab3a1bSJed Brown if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 866cab3a1bSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 876cab3a1bSJed Brown ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 886cab3a1bSJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 896cab3a1bSJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 906cab3a1bSJed Brown ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 916cab3a1bSJed Brown ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 92709ac0d2SJed Brown switch (dmdasnes->residuallocalimode) { 93709ac0d2SJed Brown case INSERT_VALUES: { 946cab3a1bSJed Brown ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 956cab3a1bSJed Brown CHKMEMQ; 966cab3a1bSJed Brown ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 976cab3a1bSJed Brown CHKMEMQ; 986cab3a1bSJed Brown ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 99709ac0d2SJed Brown } break; 100709ac0d2SJed Brown case ADD_VALUES: { 101709ac0d2SJed Brown Vec Floc; 102709ac0d2SJed Brown ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 103709ac0d2SJed Brown ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 104709ac0d2SJed Brown ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 105709ac0d2SJed Brown CHKMEMQ; 106709ac0d2SJed Brown ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 107709ac0d2SJed Brown CHKMEMQ; 108709ac0d2SJed Brown ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 109709ac0d2SJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 110709ac0d2SJed Brown ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 111709ac0d2SJed Brown ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 112709ac0d2SJed Brown ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 113709ac0d2SJed Brown } break; 114709ac0d2SJed Brown default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 115709ac0d2SJed Brown } 116709ac0d2SJed Brown ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 1176cab3a1bSJed Brown ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 1186cab3a1bSJed Brown PetscFunctionReturn(0); 1196cab3a1bSJed Brown } 1206cab3a1bSJed Brown 1216cab3a1bSJed Brown #undef __FUNCT__ 1222a4ee8f2SPeter Brune #define __FUNCT__ "SNESComputeObjective_DMDA" 1232a4ee8f2SPeter Brune /* 1242a4ee8f2SPeter Brune This function should eventually replace: 1252a4ee8f2SPeter Brune DMDAComputeFunction() and DMDAComputeFunction1() 1262a4ee8f2SPeter Brune */ 1272a4ee8f2SPeter Brune static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx) 1282a4ee8f2SPeter Brune { 1292a4ee8f2SPeter Brune PetscErrorCode ierr; 1302a4ee8f2SPeter Brune DM dm; 1312a4ee8f2SPeter Brune DM_DA_SNES *dmdasnes = (DM_DA_SNES*)ctx; 1322a4ee8f2SPeter Brune DMDALocalInfo info; 1332a4ee8f2SPeter Brune Vec Xloc; 1342a4ee8f2SPeter Brune void *x; 1352a4ee8f2SPeter Brune 1362a4ee8f2SPeter Brune PetscFunctionBegin; 1372a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1382a4ee8f2SPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1392a4ee8f2SPeter Brune PetscValidPointer(ob,3); 1402a4ee8f2SPeter Brune if (!dmdasnes->objectivelocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 1412a4ee8f2SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1422a4ee8f2SPeter Brune ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 1432a4ee8f2SPeter Brune ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1442a4ee8f2SPeter Brune ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1452a4ee8f2SPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 1462a4ee8f2SPeter Brune ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 1472a4ee8f2SPeter Brune CHKMEMQ; 1482a4ee8f2SPeter Brune ierr = (*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx);CHKERRQ(ierr); 1492a4ee8f2SPeter Brune CHKMEMQ; 1502a4ee8f2SPeter Brune ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 1512a4ee8f2SPeter Brune ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 1522a4ee8f2SPeter Brune PetscFunctionReturn(0); 1532a4ee8f2SPeter Brune } 1542a4ee8f2SPeter Brune 1552a4ee8f2SPeter Brune #undef __FUNCT__ 1566427ad5bSPeter Brune #define __FUNCT__ "SNESComputeLocalBlockFunction_DMDA" 1576427ad5bSPeter Brune /* 1586427ad5bSPeter Brune This function should eventually replace: 1596427ad5bSPeter Brune DMDAComputeFunction() and DMDAComputeFunction1() 1606427ad5bSPeter Brune */ 1616427ad5bSPeter Brune PetscErrorCode SNESComputeLocalBlockFunction_DMDA(SNES snes,Vec X,Vec F,void *dmctx) 1626427ad5bSPeter Brune { 1636427ad5bSPeter Brune PetscErrorCode ierr; 1646427ad5bSPeter Brune DM dm = (DM)dmctx; /* the global DMDA context */ 1656427ad5bSPeter Brune DMDALocalInfo info; 1666427ad5bSPeter Brune void *x,*f; 1676427ad5bSPeter Brune DM_DA_SNES *dmdasnes; 1686427ad5bSPeter Brune SNESDM sdm; 1696427ad5bSPeter Brune 1706427ad5bSPeter Brune PetscFunctionBegin; 1716427ad5bSPeter Brune 1726427ad5bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1736427ad5bSPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1746427ad5bSPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 1756427ad5bSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,4); 1766427ad5bSPeter Brune 1776427ad5bSPeter Brune ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 1786427ad5bSPeter Brune ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 179ec31c420SJed Brown if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 1806427ad5bSPeter Brune ierr = DMDAVecGetArray(dm,X,&x);CHKERRQ(ierr); 1816427ad5bSPeter Brune ierr = DMDAGetLocalBlockInfo(dm,&info);CHKERRQ(ierr); 1826427ad5bSPeter Brune switch (dmdasnes->residuallocalimode) { 1836427ad5bSPeter Brune case INSERT_VALUES: { 1846427ad5bSPeter Brune ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 1856427ad5bSPeter Brune CHKMEMQ; 1866427ad5bSPeter Brune ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 1876427ad5bSPeter Brune CHKMEMQ; 1886427ad5bSPeter Brune ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 1896427ad5bSPeter Brune } break; 1906427ad5bSPeter Brune case ADD_VALUES: { 1916427ad5bSPeter Brune Vec Floc; 1926427ad5bSPeter Brune ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 1936427ad5bSPeter Brune ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 1946427ad5bSPeter Brune ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 1956427ad5bSPeter Brune CHKMEMQ; 1966427ad5bSPeter Brune ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 1976427ad5bSPeter Brune CHKMEMQ; 1986427ad5bSPeter Brune ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 1996427ad5bSPeter Brune ierr = VecZeroEntries(F);CHKERRQ(ierr); 2006427ad5bSPeter Brune ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 2016427ad5bSPeter Brune ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 2026427ad5bSPeter Brune ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 2036427ad5bSPeter Brune } break; 2046427ad5bSPeter Brune default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 2056427ad5bSPeter Brune } 2066427ad5bSPeter Brune ierr = DMDAVecRestoreArray(dm,X,&x);CHKERRQ(ierr); 2076427ad5bSPeter Brune PetscFunctionReturn(0); 2086427ad5bSPeter Brune } 2096427ad5bSPeter Brune 2106427ad5bSPeter Brune #undef __FUNCT__ 2112961e153SJed Brown #define __FUNCT__ "SNESComputeJacobian_DMDA" 2122961e153SJed Brown /* 2132961e153SJed Brown This function should eventually replace: 2142961e153SJed Brown DMComputeJacobian() and DMDAComputeJacobian1() 2152961e153SJed Brown */ 2162961e153SJed Brown static PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 2172961e153SJed Brown { 2182961e153SJed Brown PetscErrorCode ierr; 2192961e153SJed Brown DM dm; 2202961e153SJed Brown DM_DA_SNES *dmdasnes = (DM_DA_SNES*)ctx; 2212961e153SJed Brown DMDALocalInfo info; 2222961e153SJed Brown Vec Xloc; 2232961e153SJed Brown void *x; 2242961e153SJed Brown 2252961e153SJed Brown PetscFunctionBegin; 2262961e153SJed Brown if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 2272961e153SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2282961e153SJed Brown 2292961e153SJed Brown if (dmdasnes->jacobianlocal) { 2302961e153SJed Brown ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 2312961e153SJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 2322961e153SJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 2332961e153SJed Brown ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 2342961e153SJed Brown ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 2352961e153SJed Brown CHKMEMQ; 2362961e153SJed Brown ierr = (*dmdasnes->jacobianlocal)(&info,x,*A,*B,mstr,dmdasnes->jacobianlocalctx);CHKERRQ(ierr); 2372961e153SJed Brown CHKMEMQ; 2382961e153SJed Brown ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 2392961e153SJed Brown ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 2402961e153SJed Brown } else { 2412961e153SJed Brown MatFDColoring fdcoloring; 2422961e153SJed Brown ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr); 2432961e153SJed Brown if (!fdcoloring) { 2442961e153SJed Brown ISColoring coloring; 2452961e153SJed Brown 246*c97a5455SBarry Smith ierr = DMCreateColoring(dm,dm->coloringtype,dm->mattype,&coloring);CHKERRQ(ierr); 2472961e153SJed Brown ierr = MatFDColoringCreate(*B,coloring,&fdcoloring);CHKERRQ(ierr); 2482961e153SJed Brown ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 2492961e153SJed Brown switch (dm->coloringtype) { 2502961e153SJed Brown case IS_COLORING_GLOBAL: 2512961e153SJed Brown ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode(*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 2522961e153SJed Brown break; 2532961e153SJed Brown default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 2542961e153SJed Brown } 2552961e153SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr); 2562961e153SJed Brown ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 2572961e153SJed Brown ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr); 2582961e153SJed Brown ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr); 2592961e153SJed Brown 2602961e153SJed Brown /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 2612961e153SJed Brown * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 2622961e153SJed Brown * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 2632961e153SJed Brown * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 2642961e153SJed Brown * taken when the PetscOList for the Vec inside MatFDColoring is destroyed. 2652961e153SJed Brown */ 2662961e153SJed Brown ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 2672961e153SJed Brown } 2682961e153SJed Brown *mstr = SAME_NONZERO_PATTERN; 2692961e153SJed Brown ierr = MatFDColoringApply(*B,fdcoloring,X,mstr,snes);CHKERRQ(ierr); 2702961e153SJed Brown } 2712961e153SJed Brown /* This will be redundant if the user called both, but it's too common to forget. */ 2722961e153SJed Brown if (*A != *B) { 2732961e153SJed Brown ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2742961e153SJed Brown ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2752961e153SJed Brown } 2762961e153SJed Brown PetscFunctionReturn(0); 2772961e153SJed Brown } 2782961e153SJed Brown 2792961e153SJed Brown #undef __FUNCT__ 280dcad5f1cSBarry Smith #define __FUNCT__ "SNESComputeLocalBlockJacobian_DMDA" 2816427ad5bSPeter Brune /* 2826427ad5bSPeter Brune This function should eventually replace: 2836427ad5bSPeter Brune DMComputeJacobian() and DMDAComputeJacobian1() 2846427ad5bSPeter Brune */ 2856427ad5bSPeter Brune PetscErrorCode SNESComputeLocalBlockJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *dmctx) 2866427ad5bSPeter Brune { 2876427ad5bSPeter Brune PetscErrorCode ierr; 2886427ad5bSPeter Brune DM dm = (DM)dmctx; /* the global DMDA context */ 2896427ad5bSPeter Brune DM_DA_SNES *dmdasnes; 2906427ad5bSPeter Brune SNESDM sdm; 2916427ad5bSPeter Brune DMDALocalInfo info; 2926427ad5bSPeter Brune void *x; 2936427ad5bSPeter Brune 2946427ad5bSPeter Brune PetscFunctionBegin; 2956427ad5bSPeter Brune ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 2966427ad5bSPeter Brune ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 297ec31c420SJed Brown if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 2986427ad5bSPeter Brune if (dmdasnes->jacobianlocal) { 2996427ad5bSPeter Brune ierr = DMDAVecGetArray(dm,X,&x);CHKERRQ(ierr); 3006427ad5bSPeter Brune CHKMEMQ; 3016427ad5bSPeter Brune ierr = (*dmdasnes->jacobianlocal)(&info,x,*A,*B,mstr,dmdasnes->jacobianlocalctx);CHKERRQ(ierr); 3026427ad5bSPeter Brune CHKMEMQ; 3036427ad5bSPeter Brune ierr = DMDAVecRestoreArray(dm,X,&x);CHKERRQ(ierr); 3046427ad5bSPeter Brune } 3056427ad5bSPeter Brune /* This will be redundant if the user called both, but it's too common to forget. */ 3066427ad5bSPeter Brune if (*A != *B) { 3076427ad5bSPeter Brune ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3086427ad5bSPeter Brune ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3096427ad5bSPeter Brune } 3106427ad5bSPeter Brune PetscFunctionReturn(0); 3116427ad5bSPeter Brune } 3126427ad5bSPeter Brune 3136427ad5bSPeter Brune #undef __FUNCT__ 3146cab3a1bSJed Brown #define __FUNCT__ "DMDASNESSetFunctionLocal" 3156cab3a1bSJed Brown /*@C 3166cab3a1bSJed Brown DMDASNESSetFunctionLocal - set a local residual evaluation function 3176cab3a1bSJed Brown 3186cab3a1bSJed Brown Logically Collective 3196cab3a1bSJed Brown 3206cab3a1bSJed Brown Input Arguments: 3216cab3a1bSJed Brown + dm - DM to associate callback with 322e3c51ba2SJed Brown . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 3236cab3a1bSJed Brown . func - local residual evaluation 3246cab3a1bSJed Brown - ctx - optional context for local residual evaluation 3256cab3a1bSJed Brown 3266cab3a1bSJed Brown Calling sequence for func: 3276cab3a1bSJed Brown + info - DMDALocalInfo defining the subdomain to evaluate the residual on 3286cab3a1bSJed Brown . x - dimensional pointer to state at which to evaluate residual 3296cab3a1bSJed Brown . f - dimensional pointer to residual, write the residual here 3306cab3a1bSJed Brown - ctx - optional context passed above 3316cab3a1bSJed Brown 3326cab3a1bSJed Brown Level: beginner 3336cab3a1bSJed Brown 3346cab3a1bSJed Brown .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 3356cab3a1bSJed Brown @*/ 336709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx) 3376cab3a1bSJed Brown { 3386cab3a1bSJed Brown PetscErrorCode ierr; 3396cab3a1bSJed Brown SNESDM sdm; 3406cab3a1bSJed Brown DM_DA_SNES *dmdasnes; 3416cab3a1bSJed Brown 3426cab3a1bSJed Brown PetscFunctionBegin; 3436cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 34431092ba1SJed Brown ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 3456cab3a1bSJed Brown ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 346709ac0d2SJed Brown dmdasnes->residuallocalimode = imode; 3476cab3a1bSJed Brown dmdasnes->residuallocal = func; 3486cab3a1bSJed Brown dmdasnes->residuallocalctx = ctx; 3496cab3a1bSJed Brown ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 3506427ad5bSPeter Brune ierr = DMSNESSetBlockFunction(dm,SNESComputeLocalBlockFunction_DMDA,dm);CHKERRQ(ierr); 3512961e153SJed Brown if (!sdm->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 3522961e153SJed Brown ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 3532961e153SJed Brown } 3542961e153SJed Brown PetscFunctionReturn(0); 3552961e153SJed Brown } 3562961e153SJed Brown 3572961e153SJed Brown #undef __FUNCT__ 3582961e153SJed Brown #define __FUNCT__ "DMDASNESSetJacobianLocal" 3592961e153SJed Brown /*@C 3602961e153SJed Brown DMDASNESSetJacobianLocal - set a local residual evaluation function 3612961e153SJed Brown 3622961e153SJed Brown Logically Collective 3632961e153SJed Brown 3642961e153SJed Brown Input Arguments: 3652961e153SJed Brown + dm - DM to associate callback with 3662961e153SJed Brown . func - local residual evaluation 3672961e153SJed Brown - ctx - optional context for local residual evaluation 3682961e153SJed Brown 3692961e153SJed Brown Calling sequence for func: 3702961e153SJed Brown + info - DMDALocalInfo defining the subdomain to evaluate the residual on 3712961e153SJed Brown . x - dimensional pointer to state at which to evaluate residual 3722961e153SJed Brown . f - dimensional pointer to residual, write the residual here 3732961e153SJed Brown - ctx - optional context passed above 3742961e153SJed Brown 3752961e153SJed Brown Level: beginner 3762961e153SJed Brown 3772961e153SJed Brown .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 3782961e153SJed Brown @*/ 3792961e153SJed Brown PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*),void *ctx) 3802961e153SJed Brown { 3812961e153SJed Brown PetscErrorCode ierr; 3822961e153SJed Brown SNESDM sdm; 3832961e153SJed Brown DM_DA_SNES *dmdasnes; 3842961e153SJed Brown 3852961e153SJed Brown PetscFunctionBegin; 3862961e153SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 38731092ba1SJed Brown ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 3882961e153SJed Brown ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 3892961e153SJed Brown dmdasnes->jacobianlocal = func; 3902961e153SJed Brown dmdasnes->jacobianlocalctx = ctx; 3912961e153SJed Brown ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 3926427ad5bSPeter Brune ierr = DMSNESSetBlockJacobian(dm,SNESComputeLocalBlockJacobian_DMDA,dm);CHKERRQ(ierr); 3936cab3a1bSJed Brown PetscFunctionReturn(0); 3946cab3a1bSJed Brown } 3952a4ee8f2SPeter Brune 3962a4ee8f2SPeter Brune 3972a4ee8f2SPeter Brune #undef __FUNCT__ 3982a4ee8f2SPeter Brune #define __FUNCT__ "DMDASNESSetObjectiveLocal" 3992a4ee8f2SPeter Brune /*@C 4002a4ee8f2SPeter Brune DMDASNESSetObjectiveLocal - set a local residual evaluation function 4012a4ee8f2SPeter Brune 4022a4ee8f2SPeter Brune Logically Collective 4032a4ee8f2SPeter Brune 4042a4ee8f2SPeter Brune Input Arguments: 4052a4ee8f2SPeter Brune + dm - DM to associate callback with 4062a4ee8f2SPeter Brune . func - local objective evaluation 4072a4ee8f2SPeter Brune - ctx - optional context for local residual evaluation 4082a4ee8f2SPeter Brune 4092a4ee8f2SPeter Brune Calling sequence for func: 4102a4ee8f2SPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 4112a4ee8f2SPeter Brune . x - dimensional pointer to state at which to evaluate residual 4122a4ee8f2SPeter Brune . ob - eventual objective value 4132a4ee8f2SPeter Brune - ctx - optional context passed above 4142a4ee8f2SPeter Brune 4152a4ee8f2SPeter Brune Level: beginner 4162a4ee8f2SPeter Brune 4172a4ee8f2SPeter Brune .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 4182a4ee8f2SPeter Brune @*/ 4192a4ee8f2SPeter Brune PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx) 4202a4ee8f2SPeter Brune { 4212a4ee8f2SPeter Brune PetscErrorCode ierr; 4222a4ee8f2SPeter Brune SNESDM sdm; 4232a4ee8f2SPeter Brune DM_DA_SNES *dmdasnes; 4242a4ee8f2SPeter Brune 4252a4ee8f2SPeter Brune PetscFunctionBegin; 4262a4ee8f2SPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 42731092ba1SJed Brown ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 4282a4ee8f2SPeter Brune ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 4292a4ee8f2SPeter Brune dmdasnes->objectivelocal = func; 4302a4ee8f2SPeter Brune dmdasnes->objectivelocalctx = ctx; 4312a4ee8f2SPeter Brune ierr = DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes);CHKERRQ(ierr); 4322a4ee8f2SPeter Brune PetscFunctionReturn(0); 4332a4ee8f2SPeter Brune } 434dcad5f1cSBarry Smith 435dcad5f1cSBarry Smith #undef __FUNCT__ 436dcad5f1cSBarry Smith #define __FUNCT__ "SNESComputePicard_DMDA" 437dcad5f1cSBarry Smith static PetscErrorCode SNESComputePicard_DMDA(SNES snes,Vec X,Vec F,void *ctx) 438dcad5f1cSBarry Smith { 439dcad5f1cSBarry Smith PetscErrorCode ierr; 440dcad5f1cSBarry Smith DM dm; 441dcad5f1cSBarry Smith DM_DA_SNES *dmdasnes = (DM_DA_SNES*)ctx; 442dcad5f1cSBarry Smith DMDALocalInfo info; 443dcad5f1cSBarry Smith Vec Xloc; 444dcad5f1cSBarry Smith void *x,*f; 445dcad5f1cSBarry Smith 446dcad5f1cSBarry Smith PetscFunctionBegin; 447dcad5f1cSBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 448dcad5f1cSBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,2); 449dcad5f1cSBarry Smith PetscValidHeaderSpecific(F,VEC_CLASSID,3); 450dcad5f1cSBarry Smith if (!dmdasnes->rhsplocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 451dcad5f1cSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 452dcad5f1cSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 453dcad5f1cSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 454dcad5f1cSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 455dcad5f1cSBarry Smith ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 456dcad5f1cSBarry Smith ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 457dcad5f1cSBarry Smith switch (dmdasnes->residuallocalimode) { 458dcad5f1cSBarry Smith case INSERT_VALUES: { 459dcad5f1cSBarry Smith ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 460dcad5f1cSBarry Smith CHKMEMQ; 461dcad5f1cSBarry Smith ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr); 462dcad5f1cSBarry Smith CHKMEMQ; 463dcad5f1cSBarry Smith ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 464dcad5f1cSBarry Smith } break; 465dcad5f1cSBarry Smith case ADD_VALUES: { 466dcad5f1cSBarry Smith Vec Floc; 467dcad5f1cSBarry Smith ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 468dcad5f1cSBarry Smith ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 469dcad5f1cSBarry Smith ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 470dcad5f1cSBarry Smith CHKMEMQ; 471dcad5f1cSBarry Smith ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr); 472dcad5f1cSBarry Smith CHKMEMQ; 473dcad5f1cSBarry Smith ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 474dcad5f1cSBarry Smith ierr = VecZeroEntries(F);CHKERRQ(ierr); 475dcad5f1cSBarry Smith ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 476dcad5f1cSBarry Smith ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 477dcad5f1cSBarry Smith ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 478dcad5f1cSBarry Smith } break; 479dcad5f1cSBarry Smith default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 480dcad5f1cSBarry Smith } 481dcad5f1cSBarry Smith ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 482dcad5f1cSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 483dcad5f1cSBarry Smith PetscFunctionReturn(0); 484dcad5f1cSBarry Smith } 485dcad5f1cSBarry Smith 486dcad5f1cSBarry Smith #undef __FUNCT__ 487dcad5f1cSBarry Smith #define __FUNCT__ "SNESComputePicardJacobian_DMDA" 488dcad5f1cSBarry Smith static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 489dcad5f1cSBarry Smith { 490dcad5f1cSBarry Smith PetscErrorCode ierr; 491dcad5f1cSBarry Smith DM dm; 492dcad5f1cSBarry Smith DM_DA_SNES *dmdasnes = (DM_DA_SNES*)ctx; 493dcad5f1cSBarry Smith DMDALocalInfo info; 494dcad5f1cSBarry Smith Vec Xloc; 495dcad5f1cSBarry Smith void *x; 496dcad5f1cSBarry Smith 497dcad5f1cSBarry Smith PetscFunctionBegin; 498dcad5f1cSBarry Smith if (!dmdasnes->jacobianplocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 499dcad5f1cSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 500dcad5f1cSBarry Smith 501dcad5f1cSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 502dcad5f1cSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 503dcad5f1cSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 504dcad5f1cSBarry Smith ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 505dcad5f1cSBarry Smith ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 506dcad5f1cSBarry Smith CHKMEMQ; 507dcad5f1cSBarry Smith ierr = (*dmdasnes->jacobianplocal)(&info,x,*A,*B,mstr,dmdasnes->picardlocalctx);CHKERRQ(ierr); 508dcad5f1cSBarry Smith CHKMEMQ; 509dcad5f1cSBarry Smith ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 510dcad5f1cSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 511dcad5f1cSBarry Smith *mstr = SAME_NONZERO_PATTERN; 512dcad5f1cSBarry Smith if (*A != *B) { 513dcad5f1cSBarry Smith ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 514dcad5f1cSBarry Smith ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 515dcad5f1cSBarry Smith } 516dcad5f1cSBarry Smith PetscFunctionReturn(0); 517dcad5f1cSBarry Smith } 518dcad5f1cSBarry Smith 519dcad5f1cSBarry Smith #undef __FUNCT__ 520dcad5f1cSBarry Smith #define __FUNCT__ "DMDASNESSetPicardLocal" 521dcad5f1cSBarry Smith /*@C 522dcad5f1cSBarry Smith DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration 523dcad5f1cSBarry Smith 524dcad5f1cSBarry Smith Logically Collective 525dcad5f1cSBarry Smith 526dcad5f1cSBarry Smith Input Arguments: 527dcad5f1cSBarry Smith + dm - DM to associate callback with 528dcad5f1cSBarry Smith . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 529dcad5f1cSBarry Smith . func - local residual evaluation 530dcad5f1cSBarry Smith - ctx - optional context for local residual evaluation 531dcad5f1cSBarry Smith 532dcad5f1cSBarry Smith Calling sequence for func: 533dcad5f1cSBarry Smith + info - DMDALocalInfo defining the subdomain to evaluate the residual on 534dcad5f1cSBarry Smith . x - dimensional pointer to state at which to evaluate residual 535dcad5f1cSBarry Smith . f - dimensional pointer to residual, write the residual here 536dcad5f1cSBarry Smith - ctx - optional context passed above 537dcad5f1cSBarry Smith 538dcad5f1cSBarry Smith Notes: The user must use 539dcad5f1cSBarry Smith extern PetscErrorCode SNESPicardComputeFunction(SNES,Vec,Vec,void*); 540dcad5f1cSBarry Smith extern PetscErrorCode SNESPicardComputeJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 541dcad5f1cSBarry Smith ierr = SNESSetFunction(snes,PETSC_NULL,SNESPicardComputeFunction,&user);CHKERRQ(ierr); 542dcad5f1cSBarry Smith in their code before calling this routine. 543dcad5f1cSBarry Smith 544dcad5f1cSBarry Smith 545dcad5f1cSBarry Smith Level: beginner 546dcad5f1cSBarry Smith 547dcad5f1cSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 548dcad5f1cSBarry Smith @*/ 549dcad5f1cSBarry Smith PetscErrorCode DMDASNESSetPicardLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*), 550dcad5f1cSBarry Smith PetscErrorCode (*jac)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*),void *ctx) 551dcad5f1cSBarry Smith { 552dcad5f1cSBarry Smith PetscErrorCode ierr; 553dcad5f1cSBarry Smith SNESDM sdm; 554dcad5f1cSBarry Smith DM_DA_SNES *dmdasnes; 555dcad5f1cSBarry Smith 556dcad5f1cSBarry Smith PetscFunctionBegin; 557dcad5f1cSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 558dcad5f1cSBarry Smith ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 559dcad5f1cSBarry Smith ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 560dcad5f1cSBarry Smith dmdasnes->residuallocalimode = imode; 561dcad5f1cSBarry Smith dmdasnes->rhsplocal = func; 562dcad5f1cSBarry Smith dmdasnes->jacobianplocal = jac; 563dcad5f1cSBarry Smith dmdasnes->picardlocalctx = ctx; 564dcad5f1cSBarry Smith ierr = DMSNESSetPicard(dm,SNESComputePicard_DMDA,SNESComputePicardJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 565dcad5f1cSBarry Smith PetscFunctionReturn(0); 566dcad5f1cSBarry Smith } 567