xref: /petsc/src/ts/tutorials/ex29.c (revision 9566063d113dddea24716c546802770db7481bc0)
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*9566063dSJacob Faibussowitsch   PetscCall(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*9566063dSJacob Faibussowitsch   PetscCall(DMDACreateCompatibleDMDA(dmc,2,&cdmc));
49*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)dmc,"coefficientdm",(PetscObject)cdmc));
50c4762a1bSJed Brown 
51*9566063dSJacob Faibussowitsch   PetscCall(DMGetNamedGlobalVector(cdm,"coefficient",&c));
52*9566063dSJacob Faibussowitsch   PetscCall(DMGetNamedGlobalVector(cdmc,"coefficient",&cc));
53*9566063dSJacob Faibussowitsch   PetscCall(DMGetNamedLocalVector(cdmc,"coefficient",&ccl));
54c4762a1bSJed Brown 
55*9566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(cdmc,cdm,&J,&vscale));
56*9566063dSJacob Faibussowitsch   PetscCall(MatRestrict(J,c,cc));
57*9566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(cc,vscale,cc));
58c4762a1bSJed Brown 
59*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
60*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vscale));
61c4762a1bSJed Brown 
62*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(cdmc,cc,INSERT_VALUES,ccl));
63*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(cdmc,cc,INSERT_VALUES,ccl));
64c4762a1bSJed Brown 
65*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreNamedGlobalVector(cdm,"coefficient",&c));
66*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreNamedGlobalVector(cdmc,"coefficient",&cc));
67*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreNamedLocalVector(cdmc,"coefficient",&ccl));
68c4762a1bSJed Brown 
69*9566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dmc,CoefficientCoarsenHook,NULL,NULL));
70*9566063dSJacob Faibussowitsch   PetscCall(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*9566063dSJacob Faibussowitsch   PetscCall(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*9566063dSJacob Faibussowitsch   PetscCall(DMDACreateCompatibleDMDA(subdm,2,&csubdm));
88*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)subdm,"coefficientdm",(PetscObject)csubdm));
89c4762a1bSJed Brown 
90*9566063dSJacob Faibussowitsch   PetscCall(DMGetNamedGlobalVector(cdm,"coefficient",&c));
91*9566063dSJacob Faibussowitsch   PetscCall(DMGetNamedLocalVector(csubdm,"coefficient",&cc));
92c4762a1bSJed Brown 
93*9566063dSJacob Faibussowitsch   PetscCall(DMCreateDomainDecompositionScatters(cdm,1,&csubdm,&iscat,&oscat,&gscat));
94c4762a1bSJed Brown 
95*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(*gscat,c,cc,INSERT_VALUES,SCATTER_FORWARD));
96*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(*gscat,c,cc,INSERT_VALUES,SCATTER_FORWARD));
97c4762a1bSJed Brown 
98*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(iscat));
99*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(oscat));
100*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(gscat));
101*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(iscat));
102*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(oscat));
103*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(gscat));
104c4762a1bSJed Brown 
105*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreNamedGlobalVector(cdm,"coefficient",&c));
106*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreNamedLocalVector(csubdm,"coefficient",&cc));
107c4762a1bSJed Brown 
108*9566063dSJacob Faibussowitsch   PetscCall(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*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
120*9566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
121*9566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSARKIMEX));
122*9566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
123*9566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
124*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
125*9566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
126*9566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
127c4762a1bSJed Brown 
128*9566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da,0,"u"));
129*9566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&x));
130c4762a1bSJed Brown 
131*9566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
132c4762a1bSJed Brown 
133*9566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(da,NULL,x));
134*9566063dSJacob Faibussowitsch   PetscCall(DMDATSSetIFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*))FormIFunctionLocal,NULL));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* set up the coefficient */
137c4762a1bSJed Brown 
138*9566063dSJacob Faibussowitsch   PetscCall(DMDACreateCompatibleDMDA(da,2,&cda));
139*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)da,"coefficientdm",(PetscObject)cda));
140c4762a1bSJed Brown 
141*9566063dSJacob Faibussowitsch   PetscCall(DMGetNamedGlobalVector(cda,"coefficient",&c));
142*9566063dSJacob Faibussowitsch   PetscCall(DMGetNamedLocalVector(cda,"coefficient",&clocal));
143c4762a1bSJed Brown 
144*9566063dSJacob Faibussowitsch   PetscCall(FormDiffusionCoefficient(cda,NULL,c));
145c4762a1bSJed Brown 
146*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(cda,c,INSERT_VALUES,clocal));
147*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(cda,c,INSERT_VALUES,clocal));
148c4762a1bSJed Brown 
149*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreNamedLocalVector(cda,"coefficient",&clocal));
150*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreNamedGlobalVector(cda,"coefficient",&c));
151c4762a1bSJed Brown 
152*9566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(da,CoefficientCoarsenHook,NULL,NULL));
153*9566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(da,CoefficientSubDomainRestrictHook,NULL,NULL));
154c4762a1bSJed Brown 
155*9566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts,10000));
156*9566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,10000.0));
157*9566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
158*9566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,0.05));
159*9566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,x));
160*9566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
161c4762a1bSJed Brown 
162*9566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,x));
163c4762a1bSJed Brown 
164*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
165*9566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
166*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
167*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&cda));
168c4762a1bSJed Brown 
169*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
170b122ec5aSJacob 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;
182*9566063dSJacob Faibussowitsch   PetscCall(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 
184*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,X,&x));
185*9566063dSJacob Faibussowitsch   PetscCall(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 
195*9566063dSJacob Faibussowitsch   PetscCall(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;
207*9566063dSJacob Faibussowitsch   PetscCall(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);
211*9566063dSJacob Faibussowitsch   PetscCall(VecMin(X,NULL,&min));
212c4762a1bSJed Brown    */
213c4762a1bSJed Brown 
214*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,X,&x));
215*9566063dSJacob Faibussowitsch   PetscCall(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 
227*9566063dSJacob Faibussowitsch   PetscCall(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;
243*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)info->da,"coefficientdm",(PetscObject*)&cdm));
244*9566063dSJacob Faibussowitsch   PetscCall(DMGetNamedLocalVector(cdm,"coefficient",&C));
245*9566063dSJacob Faibussowitsch   PetscCall(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   }
283*9566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(11.*info->ym*info->xm));
284c4762a1bSJed Brown 
285*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cdm,C,&c));
286*9566063dSJacob Faibussowitsch   PetscCall(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