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