1c4762a1bSJed Brown static char help[] = "Time-Dependent Allan-Cahn example in 2D with Varying Coefficients"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown This example is mainly here to show how to transfer coefficients between subdomains and levels in 5c4762a1bSJed Brown multigrid and domain decomposition. 6c4762a1bSJed Brown */ 7c4762a1bSJed Brown 8c4762a1bSJed Brown /*T 9c4762a1bSJed Brown Concepts: Alan-Cahn 10c4762a1bSJed Brown Concepts: DMDA with timestepping 11c4762a1bSJed Brown Concepts: Variable Coefficient 12c4762a1bSJed Brown Concepts: Nonlinear Domain Decomposition and Multigrid with Coefficients 13c4762a1bSJed Brown Processors: n 14c4762a1bSJed Brown T*/ 15c4762a1bSJed Brown 16c4762a1bSJed Brown #include <petscdm.h> 17c4762a1bSJed Brown #include <petscdmda.h> 18c4762a1bSJed Brown #include <petscsnes.h> 19c4762a1bSJed Brown #include <petscts.h> 20c4762a1bSJed Brown 21c4762a1bSJed Brown typedef struct { 22c4762a1bSJed Brown PetscScalar epsilon; 23c4762a1bSJed Brown PetscScalar beta; 24c4762a1bSJed Brown } Coeff; 25c4762a1bSJed Brown 26c4762a1bSJed Brown typedef struct { 27c4762a1bSJed Brown PetscScalar u; 28c4762a1bSJed Brown } Field; 29c4762a1bSJed Brown 30c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(DM da,void *ctx,Vec X); 31c4762a1bSJed Brown extern PetscErrorCode FormDiffusionCoefficient(DM da,void *ctx,Vec X); 32c4762a1bSJed Brown extern PetscErrorCode FormIFunctionLocal(DMDALocalInfo*,PetscReal,Field**,Field**,Field**,void*); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* hooks */ 35c4762a1bSJed Brown 36c4762a1bSJed Brown static PetscErrorCode CoefficientCoarsenHook(DM dm, DM dmc,void *ctx) 37c4762a1bSJed Brown { 38c4762a1bSJed Brown Vec c,cc,ccl; 39c4762a1bSJed Brown PetscErrorCode ierr; 40c4762a1bSJed Brown Mat J; 41c4762a1bSJed Brown Vec vscale; 42c4762a1bSJed Brown DM cdm,cdmc; 43c4762a1bSJed Brown 44c4762a1bSJed Brown PetscFunctionBegin; 45c4762a1bSJed Brown ierr = PetscObjectQuery((PetscObject)dm,"coefficientdm",(PetscObject*)&cdm);CHKERRQ(ierr); 46c4762a1bSJed Brown 47*3c633725SBarry Smith PetscCheck(cdm,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"The coefficient DM needs to be set up!"); 48c4762a1bSJed Brown 49c4762a1bSJed Brown ierr = DMDACreateCompatibleDMDA(dmc,2,&cdmc);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject)dmc,"coefficientdm",(PetscObject)cdmc);CHKERRQ(ierr); 51c4762a1bSJed Brown 52c4762a1bSJed Brown ierr = DMGetNamedGlobalVector(cdm,"coefficient",&c);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = DMGetNamedGlobalVector(cdmc,"coefficient",&cc);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = DMGetNamedLocalVector(cdmc,"coefficient",&ccl);CHKERRQ(ierr); 55c4762a1bSJed Brown 56c4762a1bSJed Brown ierr = DMCreateInterpolation(cdmc,cdm,&J,&vscale);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = MatRestrict(J,c,cc);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = VecPointwiseMult(cc,vscale,cc);CHKERRQ(ierr); 59c4762a1bSJed Brown 60c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = VecDestroy(&vscale);CHKERRQ(ierr); 62c4762a1bSJed Brown 63c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(cdmc,cc,INSERT_VALUES,ccl);CHKERRQ(ierr); 64c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(cdmc,cc,INSERT_VALUES,ccl);CHKERRQ(ierr); 65c4762a1bSJed Brown 66c4762a1bSJed Brown ierr = DMRestoreNamedGlobalVector(cdm,"coefficient",&c);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = DMRestoreNamedGlobalVector(cdmc,"coefficient",&cc);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = DMRestoreNamedLocalVector(cdmc,"coefficient",&ccl);CHKERRQ(ierr); 69c4762a1bSJed Brown 70c4762a1bSJed Brown ierr = DMCoarsenHookAdd(dmc,CoefficientCoarsenHook,NULL,NULL);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = DMDestroy(&cdmc);CHKERRQ(ierr); 72c4762a1bSJed Brown PetscFunctionReturn(0); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* This could restrict auxiliary information to the coarse level. 76c4762a1bSJed Brown */ 77c4762a1bSJed Brown static PetscErrorCode CoefficientSubDomainRestrictHook(DM dm,DM subdm,void *ctx) 78c4762a1bSJed Brown { 79c4762a1bSJed Brown Vec c,cc; 80c4762a1bSJed Brown DM cdm,csubdm; 81c4762a1bSJed Brown PetscErrorCode ierr; 82c4762a1bSJed Brown VecScatter *iscat,*oscat,*gscat; 83c4762a1bSJed Brown 84c4762a1bSJed Brown PetscFunctionBegin; 85c4762a1bSJed Brown ierr = PetscObjectQuery((PetscObject)dm,"coefficientdm",(PetscObject*)&cdm);CHKERRQ(ierr); 86c4762a1bSJed Brown 87*3c633725SBarry Smith PetscCheck(cdm,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"The coefficient DM needs to be set up!"); 88c4762a1bSJed Brown 89c4762a1bSJed Brown ierr = DMDACreateCompatibleDMDA(subdm,2,&csubdm);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject)subdm,"coefficientdm",(PetscObject)csubdm);CHKERRQ(ierr); 91c4762a1bSJed Brown 92c4762a1bSJed Brown ierr = DMGetNamedGlobalVector(cdm,"coefficient",&c);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = DMGetNamedLocalVector(csubdm,"coefficient",&cc);CHKERRQ(ierr); 94c4762a1bSJed Brown 95c4762a1bSJed Brown ierr = DMCreateDomainDecompositionScatters(cdm,1,&csubdm,&iscat,&oscat,&gscat);CHKERRQ(ierr); 96c4762a1bSJed Brown 97c4762a1bSJed Brown ierr = VecScatterBegin(*gscat,c,cc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = VecScatterEnd(*gscat,c,cc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 99c4762a1bSJed Brown 100c4762a1bSJed Brown ierr = VecScatterDestroy(iscat);CHKERRQ(ierr); 101c4762a1bSJed Brown ierr = VecScatterDestroy(oscat);CHKERRQ(ierr); 102c4762a1bSJed Brown ierr = VecScatterDestroy(gscat);CHKERRQ(ierr); 103c4762a1bSJed Brown ierr = PetscFree(iscat);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = PetscFree(oscat);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = PetscFree(gscat);CHKERRQ(ierr); 106c4762a1bSJed Brown 107c4762a1bSJed Brown ierr = DMRestoreNamedGlobalVector(cdm,"coefficient",&c);CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = DMRestoreNamedLocalVector(csubdm,"coefficient",&cc);CHKERRQ(ierr); 109c4762a1bSJed Brown 110c4762a1bSJed Brown ierr = DMDestroy(&csubdm);CHKERRQ(ierr); 111c4762a1bSJed Brown PetscFunctionReturn(0); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114c4762a1bSJed Brown int main(int argc,char **argv) 115c4762a1bSJed Brown 116c4762a1bSJed Brown { 117c4762a1bSJed Brown TS ts; 118c4762a1bSJed Brown Vec x,c,clocal; 119c4762a1bSJed Brown PetscErrorCode ierr; 120c4762a1bSJed Brown DM da,cda; 121c4762a1bSJed Brown 122c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 123c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr); 127c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 129c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr); 130c4762a1bSJed Brown 131c4762a1bSJed Brown ierr = DMDASetFieldName(da,0,"u");CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr); 133c4762a1bSJed Brown 134c4762a1bSJed Brown ierr = TSSetDM(ts, da);CHKERRQ(ierr); 135c4762a1bSJed Brown 136c4762a1bSJed Brown ierr = FormInitialGuess(da,NULL,x);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = DMDATSSetIFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*))FormIFunctionLocal,NULL);CHKERRQ(ierr); 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* set up the coefficient */ 140c4762a1bSJed Brown 141c4762a1bSJed Brown ierr = DMDACreateCompatibleDMDA(da,2,&cda);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject)da,"coefficientdm",(PetscObject)cda);CHKERRQ(ierr); 143c4762a1bSJed Brown 144c4762a1bSJed Brown ierr = DMGetNamedGlobalVector(cda,"coefficient",&c);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = DMGetNamedLocalVector(cda,"coefficient",&clocal);CHKERRQ(ierr); 146c4762a1bSJed Brown 147c4762a1bSJed Brown ierr = FormDiffusionCoefficient(cda,NULL,c);CHKERRQ(ierr); 148c4762a1bSJed Brown 149c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(cda,c,INSERT_VALUES,clocal);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(cda,c,INSERT_VALUES,clocal);CHKERRQ(ierr); 151c4762a1bSJed Brown 152c4762a1bSJed Brown ierr = DMRestoreNamedLocalVector(cda,"coefficient",&clocal);CHKERRQ(ierr); 153c4762a1bSJed Brown ierr = DMRestoreNamedGlobalVector(cda,"coefficient",&c);CHKERRQ(ierr); 154c4762a1bSJed Brown 155c4762a1bSJed Brown ierr = DMCoarsenHookAdd(da,CoefficientCoarsenHook,NULL,NULL);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = DMSubDomainHookAdd(da,CoefficientSubDomainRestrictHook,NULL,NULL);CHKERRQ(ierr); 157c4762a1bSJed Brown 158c4762a1bSJed Brown ierr = TSSetMaxSteps(ts,10000);CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = TSSetMaxTime(ts,10000.0);CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.05);CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 164c4762a1bSJed Brown 165c4762a1bSJed Brown ierr = TSSolve(ts,x);CHKERRQ(ierr); 166c4762a1bSJed Brown 167c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = DMDestroy(&cda);CHKERRQ(ierr); 171c4762a1bSJed Brown 172c4762a1bSJed Brown ierr = PetscFinalize(); 173c4762a1bSJed Brown return ierr; 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 177c4762a1bSJed Brown 178c4762a1bSJed Brown PetscErrorCode FormInitialGuess(DM da,void *ctx,Vec X) 179c4762a1bSJed Brown { 180c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 181c4762a1bSJed Brown PetscErrorCode ierr; 182c4762a1bSJed Brown Field **x; 183c4762a1bSJed Brown PetscReal x0,x1; 184c4762a1bSJed Brown 185c4762a1bSJed Brown PetscFunctionBeginUser; 186c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr); 187c4762a1bSJed Brown 188c4762a1bSJed Brown ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr); 189c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 190c4762a1bSJed Brown 191c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 192c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 193c4762a1bSJed Brown x0 = 10.0*(i - 0.5*(Mx-1)) / (Mx-1); 194c4762a1bSJed Brown x1 = 10.0*(j - 0.5*(Mx-1)) / (My-1); 195c4762a1bSJed Brown x[j][i].u = PetscCosReal(2.0*PetscSqrtReal(x1*x1 + x0*x0)); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr); 200c4762a1bSJed Brown PetscFunctionReturn(0); 201c4762a1bSJed Brown 202c4762a1bSJed Brown } 203c4762a1bSJed Brown 204c4762a1bSJed Brown PetscErrorCode FormDiffusionCoefficient(DM da,void *ctx,Vec X) 205c4762a1bSJed Brown { 206c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 207c4762a1bSJed Brown PetscErrorCode ierr; 208c4762a1bSJed Brown Coeff **x; 209c4762a1bSJed Brown PetscReal x1,x0; 210c4762a1bSJed Brown 211c4762a1bSJed Brown PetscFunctionBeginUser; 212c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr); 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* 215c4762a1bSJed Brown ierr = VecSetRandom(X,NULL); 216c4762a1bSJed Brown ierr = VecMin(X,NULL,&min);CHKERRQ(ierr); 217c4762a1bSJed Brown */ 218c4762a1bSJed Brown 219c4762a1bSJed Brown ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 221c4762a1bSJed Brown 222c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 223c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 224c4762a1bSJed Brown x0 = 10.0*(i - 0.5*(Mx-1)) / (Mx-1); 225c4762a1bSJed Brown x1 = 10.0*(j - 0.5*(My-1)) / (My-1); 226c4762a1bSJed Brown 227c4762a1bSJed Brown x[j][i].epsilon = 0.0; 228c4762a1bSJed Brown x[j][i].beta = 0.05+0.05*PetscSqrtReal(x0*x0+x1*x1); 229c4762a1bSJed Brown } 230c4762a1bSJed Brown } 231c4762a1bSJed Brown 232c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr); 233c4762a1bSJed Brown PetscFunctionReturn(0); 234c4762a1bSJed Brown 235c4762a1bSJed Brown } 236c4762a1bSJed Brown 237c4762a1bSJed Brown PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info,PetscReal ptime,Field **x,Field **xt,Field **f,void *ctx) 238c4762a1bSJed Brown { 239c4762a1bSJed Brown PetscErrorCode ierr; 240c4762a1bSJed Brown PetscInt i,j; 241c4762a1bSJed Brown PetscReal hx,hy,dhx,dhy,hxdhy,hydhx,scale; 242c4762a1bSJed Brown PetscScalar u,uxx,uyy; 243c4762a1bSJed Brown PetscScalar ux,uy,bx,by; 244c4762a1bSJed Brown Vec C; 245c4762a1bSJed Brown Coeff **c; 246c4762a1bSJed Brown DM cdm; 247c4762a1bSJed Brown 248c4762a1bSJed Brown PetscFunctionBeginUser; 249c4762a1bSJed Brown ierr = PetscObjectQuery((PetscObject)info->da,"coefficientdm",(PetscObject*)&cdm);CHKERRQ(ierr); 250c4762a1bSJed Brown ierr = DMGetNamedLocalVector(cdm,"coefficient",&C);CHKERRQ(ierr); 251c4762a1bSJed Brown ierr = DMDAVecGetArray(cdm,C,&c);CHKERRQ(ierr); 252c4762a1bSJed Brown 253c4762a1bSJed Brown hx = 10.0/((PetscReal)(info->mx-1)); 254c4762a1bSJed Brown hy = 10.0/((PetscReal)(info->my-1)); 255c4762a1bSJed Brown 256c4762a1bSJed Brown dhx = 1. / hx; 257c4762a1bSJed Brown dhy = 1. / hy; 258c4762a1bSJed Brown 259c4762a1bSJed Brown hxdhy = hx/hy; 260c4762a1bSJed Brown hydhx = hy/hx; 261c4762a1bSJed Brown scale = hx*hy; 262c4762a1bSJed Brown 263c4762a1bSJed Brown for (j=info->ys; j<info->ys+info->ym; j++) { 264c4762a1bSJed Brown for (i=info->xs; i<info->xs+info->xm; i++) { 265c4762a1bSJed Brown f[j][i].u = xt[j][i].u*scale; 266c4762a1bSJed Brown 267c4762a1bSJed Brown u = x[j][i].u; 268c4762a1bSJed Brown 269c4762a1bSJed Brown f[j][i].u += scale*(u*u - 1.)*u; 270c4762a1bSJed Brown 271c4762a1bSJed Brown if (i == 0) f[j][i].u += (x[j][i].u - x[j][i+1].u)*dhx; 272c4762a1bSJed Brown else if (i == info->mx-1) f[j][i].u += (x[j][i].u - x[j][i-1].u)*dhx; 273c4762a1bSJed Brown else if (j == 0) f[j][i].u += (x[j][i].u - x[j+1][i].u)*dhy; 274c4762a1bSJed Brown else if (j == info->my-1) f[j][i].u += (x[j][i].u - x[j-1][i].u)*dhy; 275c4762a1bSJed Brown else { 276c4762a1bSJed Brown uyy = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy; 277c4762a1bSJed Brown uxx = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx; 278c4762a1bSJed Brown 279c4762a1bSJed Brown bx = 0.5*(c[j][i+1].beta - c[j][i-1].beta)*dhx; 280c4762a1bSJed Brown by = 0.5*(c[j+1][i].beta - c[j-1][i].beta)*dhy; 281c4762a1bSJed Brown 282c4762a1bSJed Brown ux = 0.5*(x[j][i+1].u - x[j][i-1].u)*dhx; 283c4762a1bSJed Brown uy = 0.5*(x[j+1][i].u - x[j-1][i].u)*dhy; 284c4762a1bSJed Brown 285c4762a1bSJed Brown f[j][i].u += c[j][i].beta*(uxx + uyy) + scale*(bx*ux + by*uy); 286c4762a1bSJed Brown } 287c4762a1bSJed Brown } 288c4762a1bSJed Brown } 289c4762a1bSJed Brown ierr = PetscLogFlops(11.*info->ym*info->xm);CHKERRQ(ierr); 290c4762a1bSJed Brown 291c4762a1bSJed Brown ierr = DMDAVecRestoreArray(cdm,C,&c);CHKERRQ(ierr); 292c4762a1bSJed Brown ierr = DMRestoreNamedLocalVector(cdm,"coefficient",&C);CHKERRQ(ierr); 293c4762a1bSJed Brown PetscFunctionReturn(0); 294c4762a1bSJed Brown } 295c4762a1bSJed Brown 296c4762a1bSJed Brown /*TEST 297c4762a1bSJed Brown 298c4762a1bSJed Brown test: 299c4762a1bSJed Brown args: -da_refine 4 -ts_max_steps 10 -ts_rtol 1e-3 -ts_atol 1e-3 -ts_type arkimex -ts_monitor -snes_monitor -snes_type ngmres -npc_snes_type nasm -npc_snes_nasm_type restrict -da_overlap 4 300c4762a1bSJed Brown nsize: 16 301c4762a1bSJed Brown requires: !single 302c4762a1bSJed Brown output_file: output/ex29.out 303c4762a1bSJed Brown 304c4762a1bSJed Brown TEST*/ 305