xref: /petsc/src/ts/tutorials/ex29.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
445f80ce2aSJacob 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 
485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreateCompatibleDMDA(dmc,2,&cdmc));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectCompose((PetscObject)dmc,"coefficientdm",(PetscObject)cdmc));
50c4762a1bSJed Brown 
515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetNamedGlobalVector(cdm,"coefficient",&c));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetNamedGlobalVector(cdmc,"coefficient",&cc));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetNamedLocalVector(cdmc,"coefficient",&ccl));
54c4762a1bSJed Brown 
555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateInterpolation(cdmc,cdm,&J,&vscale));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatRestrict(J,c,cc));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(cc,vscale,cc));
58c4762a1bSJed Brown 
595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&vscale));
61c4762a1bSJed Brown 
625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(cdmc,cc,INSERT_VALUES,ccl));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(cdmc,cc,INSERT_VALUES,ccl));
64c4762a1bSJed Brown 
655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreNamedGlobalVector(cdm,"coefficient",&c));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreNamedGlobalVector(cdmc,"coefficient",&cc));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreNamedLocalVector(cdmc,"coefficient",&ccl));
68c4762a1bSJed Brown 
695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCoarsenHookAdd(dmc,CoefficientCoarsenHook,NULL,NULL));
705f80ce2aSJacob 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;
835f80ce2aSJacob 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 
875f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreateCompatibleDMDA(subdm,2,&csubdm));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectCompose((PetscObject)subdm,"coefficientdm",(PetscObject)csubdm));
89c4762a1bSJed Brown 
905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetNamedGlobalVector(cdm,"coefficient",&c));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetNamedLocalVector(csubdm,"coefficient",&cc));
92c4762a1bSJed Brown 
935f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDomainDecompositionScatters(cdm,1,&csubdm,&iscat,&oscat,&gscat));
94c4762a1bSJed Brown 
955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(*gscat,c,cc,INSERT_VALUES,SCATTER_FORWARD));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(*gscat,c,cc,INSERT_VALUES,SCATTER_FORWARD));
97c4762a1bSJed Brown 
985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(iscat));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(oscat));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(gscat));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(iscat));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(oscat));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(gscat));
104c4762a1bSJed Brown 
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreNamedGlobalVector(cdm,"coefficient",&c));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreNamedLocalVector(csubdm,"coefficient",&cc));
107c4762a1bSJed Brown 
1085f80ce2aSJacob 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   DM             da,cda;
118c4762a1bSJed Brown 
119*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD, &ts));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSARKIMEX));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
1235f80ce2aSJacob 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));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
127c4762a1bSJed Brown 
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,0,"u"));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&x));
130c4762a1bSJed Brown 
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts, da));
132c4762a1bSJed Brown 
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialGuess(da,NULL,x));
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDATSSetIFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*))FormIFunctionLocal,NULL));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* set up the coefficient */
137c4762a1bSJed Brown 
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreateCompatibleDMDA(da,2,&cda));
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectCompose((PetscObject)da,"coefficientdm",(PetscObject)cda));
140c4762a1bSJed Brown 
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetNamedGlobalVector(cda,"coefficient",&c));
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetNamedLocalVector(cda,"coefficient",&clocal));
143c4762a1bSJed Brown 
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(FormDiffusionCoefficient(cda,NULL,c));
145c4762a1bSJed Brown 
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(cda,c,INSERT_VALUES,clocal));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(cda,c,INSERT_VALUES,clocal));
148c4762a1bSJed Brown 
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreNamedLocalVector(cda,"coefficient",&clocal));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreNamedGlobalVector(cda,"coefficient",&c));
151c4762a1bSJed Brown 
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCoarsenHookAdd(da,CoefficientCoarsenHook,NULL,NULL));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSubDomainHookAdd(da,CoefficientSubDomainRestrictHook,NULL,NULL));
154c4762a1bSJed Brown 
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxSteps(ts,10000));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,10000.0));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,0.05));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,x));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
161c4762a1bSJed Brown 
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,x));
163c4762a1bSJed Brown 
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&cda));
168c4762a1bSJed Brown 
169*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
170*b122ec5aSJacob Faibussowitsch   return 0;
171c4762a1bSJed Brown }
172c4762a1bSJed Brown 
173c4762a1bSJed Brown /* ------------------------------------------------------------------- */
174c4762a1bSJed Brown 
175c4762a1bSJed Brown PetscErrorCode FormInitialGuess(DM da,void *ctx,Vec X)
176c4762a1bSJed Brown {
177c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
178c4762a1bSJed Brown   Field          **x;
179c4762a1bSJed Brown   PetscReal      x0,x1;
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   PetscFunctionBeginUser;
1825f80ce2aSJacob 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));
183c4762a1bSJed Brown 
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,X,&x));
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
188c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
189c4762a1bSJed Brown       x0        = 10.0*(i - 0.5*(Mx-1)) / (Mx-1);
190c4762a1bSJed Brown       x1        = 10.0*(j - 0.5*(Mx-1)) / (My-1);
191c4762a1bSJed Brown       x[j][i].u = PetscCosReal(2.0*PetscSqrtReal(x1*x1 + x0*x0));
192c4762a1bSJed Brown     }
193c4762a1bSJed Brown   }
194c4762a1bSJed Brown 
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,X,&x));
196c4762a1bSJed Brown   PetscFunctionReturn(0);
197c4762a1bSJed Brown 
198c4762a1bSJed Brown }
199c4762a1bSJed Brown 
200c4762a1bSJed Brown PetscErrorCode FormDiffusionCoefficient(DM da,void *ctx,Vec X)
201c4762a1bSJed Brown {
202c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
203c4762a1bSJed Brown   Coeff          **x;
204c4762a1bSJed Brown   PetscReal      x1,x0;
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   PetscFunctionBeginUser;
2075f80ce2aSJacob 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));
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /*
210c4762a1bSJed Brown   ierr = VecSetRandom(X,NULL);
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecMin(X,NULL,&min));
212c4762a1bSJed Brown    */
213c4762a1bSJed Brown 
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,X,&x));
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
218c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
219c4762a1bSJed Brown       x0 = 10.0*(i - 0.5*(Mx-1)) / (Mx-1);
220c4762a1bSJed Brown       x1 = 10.0*(j - 0.5*(My-1)) / (My-1);
221c4762a1bSJed Brown 
222c4762a1bSJed Brown       x[j][i].epsilon = 0.0;
223c4762a1bSJed Brown       x[j][i].beta    = 0.05+0.05*PetscSqrtReal(x0*x0+x1*x1);
224c4762a1bSJed Brown     }
225c4762a1bSJed Brown   }
226c4762a1bSJed Brown 
2275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,X,&x));
228c4762a1bSJed Brown   PetscFunctionReturn(0);
229c4762a1bSJed Brown 
230c4762a1bSJed Brown }
231c4762a1bSJed Brown 
232c4762a1bSJed Brown PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info,PetscReal ptime,Field **x,Field **xt,Field **f,void *ctx)
233c4762a1bSJed Brown {
234c4762a1bSJed Brown   PetscInt       i,j;
235c4762a1bSJed Brown   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx,scale;
236c4762a1bSJed Brown   PetscScalar    u,uxx,uyy;
237c4762a1bSJed Brown   PetscScalar    ux,uy,bx,by;
238c4762a1bSJed Brown   Vec            C;
239c4762a1bSJed Brown   Coeff          **c;
240c4762a1bSJed Brown   DM             cdm;
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   PetscFunctionBeginUser;
2435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectQuery((PetscObject)info->da,"coefficientdm",(PetscObject*)&cdm));
2445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetNamedLocalVector(cdm,"coefficient",&C));
2455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(cdm,C,&c));
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   hx = 10.0/((PetscReal)(info->mx-1));
248c4762a1bSJed Brown   hy = 10.0/((PetscReal)(info->my-1));
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   dhx = 1. / hx;
251c4762a1bSJed Brown   dhy = 1. / hy;
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   hxdhy =  hx/hy;
254c4762a1bSJed Brown   hydhx =  hy/hx;
255c4762a1bSJed Brown   scale =   hx*hy;
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
258c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
259c4762a1bSJed Brown       f[j][i].u = xt[j][i].u*scale;
260c4762a1bSJed Brown 
261c4762a1bSJed Brown       u = x[j][i].u;
262c4762a1bSJed Brown 
263c4762a1bSJed Brown       f[j][i].u += scale*(u*u - 1.)*u;
264c4762a1bSJed Brown 
265c4762a1bSJed Brown       if (i == 0)               f[j][i].u += (x[j][i].u - x[j][i+1].u)*dhx;
266c4762a1bSJed Brown       else if (i == info->mx-1) f[j][i].u += (x[j][i].u - x[j][i-1].u)*dhx;
267c4762a1bSJed Brown       else if (j == 0)          f[j][i].u += (x[j][i].u - x[j+1][i].u)*dhy;
268c4762a1bSJed Brown       else if (j == info->my-1) f[j][i].u += (x[j][i].u - x[j-1][i].u)*dhy;
269c4762a1bSJed Brown       else {
270c4762a1bSJed Brown         uyy     = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
271c4762a1bSJed Brown         uxx     = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
272c4762a1bSJed Brown 
273c4762a1bSJed Brown         bx      = 0.5*(c[j][i+1].beta - c[j][i-1].beta)*dhx;
274c4762a1bSJed Brown         by      = 0.5*(c[j+1][i].beta - c[j-1][i].beta)*dhy;
275c4762a1bSJed Brown 
276c4762a1bSJed Brown         ux      = 0.5*(x[j][i+1].u - x[j][i-1].u)*dhx;
277c4762a1bSJed Brown         uy      = 0.5*(x[j+1][i].u - x[j-1][i].u)*dhy;
278c4762a1bSJed Brown 
279c4762a1bSJed Brown         f[j][i].u += c[j][i].beta*(uxx + uyy) + scale*(bx*ux + by*uy);
280c4762a1bSJed Brown        }
281c4762a1bSJed Brown     }
282c4762a1bSJed Brown   }
2835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(11.*info->ym*info->xm));
284c4762a1bSJed Brown 
2855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(cdm,C,&c));
2865f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreNamedLocalVector(cdm,"coefficient",&C));
287c4762a1bSJed Brown   PetscFunctionReturn(0);
288c4762a1bSJed Brown }
289c4762a1bSJed Brown 
290c4762a1bSJed Brown /*TEST
291c4762a1bSJed Brown 
292c4762a1bSJed Brown     test:
293c4762a1bSJed 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
294c4762a1bSJed Brown       nsize: 16
295c4762a1bSJed Brown       requires: !single
296c4762a1bSJed Brown       output_file: output/ex29.out
297c4762a1bSJed Brown 
298c4762a1bSJed Brown TEST*/
299