xref: /petsc/src/ts/tutorials/ex48.c (revision 4e278199b78715991f5c71ebbd945c1489263e6c)
1 static char help[] = "Evolution of magnetic islands.\n\
2 The aim of this model is to self-consistently study the interaction between the tearing mode and small scale drift-wave turbulence.\n\n\n";
3 
4 /*F
5 This is a three field model for the density $\tilde n$, vorticity $\tilde\Omega$, and magnetic flux $\tilde\psi$, using auxiliary variables potential $\tilde\phi$ and current $j_z$.
6 \begin{equation}
7   \begin{aligned}
8     \partial_t \tilde n       &= \left\{ \tilde n, \tilde\phi \right\} + \beta \left\{ j_z, \tilde\psi \right\} + \left\{ \ln n_0, \tilde\phi \right\} + \mu \nabla^2_\perp \tilde n \\
9   \partial_t \tilde\Omega   &= \left\{ \tilde\Omega, \tilde\phi \right\} + \beta \left\{ j_z, \tilde\psi \right\} + \mu \nabla^2_\perp \tilde\Omega \\
10   \partial_t \tilde\psi     &= \left\{ \psi_0 + \tilde\psi, \tilde\phi - \tilde n \right\} - \left\{ \ln n_0, \tilde\psi \right\} + \frac{\eta}{\beta} \nabla^2_\perp \tilde\psi \\
11   \nabla^2_\perp\tilde\phi        &= \tilde\Omega \\
12   j_z  &= -\nabla^2_\perp  \left(\tilde\psi + \psi_0  \right)\\
13   \end{aligned}
14 \end{equation}
15 F*/
16 
17 #include <petscdmplex.h>
18 #include <petscts.h>
19 #include <petscds.h>
20 
21 typedef struct {
22   PetscInt       debug;             /* The debugging level */
23   PetscBool      plotRef;           /* Plot the reference fields */
24   PetscReal      lower[3], upper[3];
25   /* Problem definition */
26   PetscErrorCode (**initialFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
27   PetscReal      mu, eta, beta;
28   PetscReal      a,b,Jo,Jop,m,ke,kx,ky,DeltaPrime,eps;
29   /* solver */
30   PetscBool      implicit;
31 } AppCtx;
32 
33 static AppCtx *s_ctx;
34 
35 static PetscScalar poissonBracket(PetscInt dim, const PetscScalar df[], const PetscScalar dg[])
36 {
37   PetscScalar ret = df[0]*dg[1] - df[1]*dg[0];
38   return ret;
39 }
40 
41 enum field_idx {DENSITY,OMEGA,PSI,PHI,JZ};
42 
43 static void f0_n(PetscInt dim, PetscInt Nf, PetscInt NfAux,
44                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
45                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
46                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
47 {
48   const PetscScalar *pnDer   = &u_x[uOff_x[DENSITY]];
49   const PetscScalar *ppsiDer = &u_x[uOff_x[PSI]];
50   const PetscScalar *pphiDer = &u_x[uOff_x[PHI]];
51   const PetscScalar *jzDer   = &u_x[uOff_x[JZ]];
52   const PetscScalar *logRefDenDer = &a_x[aOff_x[DENSITY]];
53   f0[0] += - poissonBracket(dim,pnDer, pphiDer) - s_ctx->beta*poissonBracket(dim,jzDer, ppsiDer) - poissonBracket(dim,logRefDenDer, pphiDer);
54   if (u_t) f0[0] += u_t[DENSITY];
55 }
56 
57 static void f1_n(PetscInt dim, PetscInt Nf, PetscInt NfAux,
58                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
59                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
60                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
61 {
62   const PetscScalar *pnDer = &u_x[uOff_x[DENSITY]];
63   PetscInt           d;
64 
65   for (d = 0; d < dim-1; ++d) f1[d] = -s_ctx->mu*pnDer[d];
66 }
67 
68 static void f0_Omega(PetscInt dim, PetscInt Nf, PetscInt NfAux,
69                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
70                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
71                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
72 {
73   const PetscScalar *pOmegaDer = &u_x[uOff_x[OMEGA]];
74   const PetscScalar *ppsiDer   = &u_x[uOff_x[PSI]];
75   const PetscScalar *pphiDer   = &u_x[uOff_x[PHI]];
76   const PetscScalar *jzDer     = &u_x[uOff_x[JZ]];
77 
78   f0[0] += - poissonBracket(dim,pOmegaDer, pphiDer) - s_ctx->beta*poissonBracket(dim,jzDer, ppsiDer);
79   if (u_t) f0[0] += u_t[OMEGA];
80 }
81 
82 static void f1_Omega(PetscInt dim, PetscInt Nf, PetscInt NfAux,
83                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
84                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
85                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
86 {
87   const PetscScalar *pOmegaDer = &u_x[uOff_x[OMEGA]];
88   PetscInt           d;
89 
90   for (d = 0; d < dim-1; ++d) f1[d] = -s_ctx->mu*pOmegaDer[d];
91 }
92 
93 static void f0_psi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
94                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
95                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
96                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
97 {
98   const PetscScalar *pnDer     = &u_x[uOff_x[DENSITY]];
99   const PetscScalar *ppsiDer   = &u_x[uOff_x[PSI]];
100   const PetscScalar *pphiDer   = &u_x[uOff_x[PHI]];
101   const PetscScalar *refPsiDer = &a_x[aOff_x[PSI]];
102   const PetscScalar *logRefDenDer= &a_x[aOff_x[DENSITY]];
103   PetscScalar       psiDer[3];
104   PetscScalar       phi_n_Der[3];
105   PetscInt          d;
106   if (dim < 2) {MPI_Abort(MPI_COMM_WORLD,1); return;} /* this is needed so that the clang static analyzer does not generate a warning about variables used by not set */
107   for (d = 0; d < dim; ++d) {
108     psiDer[d]    = refPsiDer[d] + ppsiDer[d];
109     phi_n_Der[d] = pphiDer[d]   - pnDer[d];
110   }
111   f0[0] = - poissonBracket(dim,psiDer, phi_n_Der) + poissonBracket(dim,logRefDenDer, ppsiDer);
112   if (u_t) f0[0] += u_t[PSI];
113 }
114 
115 static void f1_psi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
116                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
117                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
118                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
119 {
120   const PetscScalar *ppsi = &u_x[uOff_x[PSI]];
121   PetscInt           d;
122 
123   for (d = 0; d < dim-1; ++d) f1[d] = -(s_ctx->eta/s_ctx->beta)*ppsi[d];
124 }
125 
126 static void f0_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
127                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
128                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
129                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
130 {
131   f0[0] = -u[uOff[OMEGA]];
132 }
133 
134 static void f1_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
135                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
136                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
137                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
138 {
139   const PetscScalar *pphi = &u_x[uOff_x[PHI]];
140   PetscInt           d;
141 
142   for (d = 0; d < dim-1; ++d) f1[d] = pphi[d];
143 }
144 
145 static void f0_jz(PetscInt dim, PetscInt Nf, PetscInt NfAux,
146                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
147                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
148                   PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
149 {
150   f0[0] = u[uOff[JZ]];
151 }
152 
153 static void f1_jz(PetscInt dim, PetscInt Nf, PetscInt NfAux,
154                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
155                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
156                   PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
157 {
158   const PetscScalar *ppsi = &u_x[uOff_x[PSI]];
159   const PetscScalar *refPsiDer = &a_x[aOff_x[PSI]]; /* aOff_x[PSI] == 2*PSI */
160   PetscInt           d;
161 
162   for (d = 0; d < dim-1; ++d) f1[d] = ppsi[d] + refPsiDer[d];
163 }
164 
165 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
166 {
167   PetscErrorCode ierr;
168 
169   PetscFunctionBeginUser;
170   options->debug    = 1;
171   options->plotRef  = PETSC_FALSE;
172   options->implicit = PETSC_FALSE;
173   options->mu       = 0;
174   options->eta      = 0;
175   options->beta     = 1;
176   options->a        = 1;
177   options->b        = PETSC_PI;
178   options->Jop      = 0;
179   options->m        = 1;
180   options->eps      = 1.e-6;
181 
182   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
183   ierr = PetscOptionsInt("-debug", "The debugging level", "ex48.c", options->debug, &options->debug, NULL);CHKERRQ(ierr);
184   ierr = PetscOptionsBool("-plot_ref", "Plot the reference fields", "ex48.c", options->plotRef, &options->plotRef, NULL);CHKERRQ(ierr);
185   ierr = PetscOptionsReal("-mu", "mu", "ex48.c", options->mu, &options->mu, NULL);CHKERRQ(ierr);
186   ierr = PetscOptionsReal("-eta", "eta", "ex48.c", options->eta, &options->eta, NULL);CHKERRQ(ierr);
187   ierr = PetscOptionsReal("-beta", "beta", "ex48.c", options->beta, &options->beta, NULL);CHKERRQ(ierr);
188   ierr = PetscOptionsReal("-Jop", "Jop", "ex48.c", options->Jop, &options->Jop, NULL);CHKERRQ(ierr);
189   ierr = PetscOptionsReal("-m", "m", "ex48.c", options->m, &options->m, NULL);CHKERRQ(ierr);
190   ierr = PetscOptionsReal("-eps", "eps", "ex48.c", options->eps, &options->eps, NULL);CHKERRQ(ierr);
191   ierr = PetscOptionsBool("-implicit", "Use implicit time integrator", "ex48.c", options->implicit, &options->implicit, NULL);CHKERRQ(ierr);
192   ierr = PetscOptionsEnd();CHKERRQ(ierr);
193   options->ke = PetscSqrtScalar(options->Jop);
194   if (options->Jop==0.0) {
195     options->Jo = 1.0/PetscPowScalar(options->a,2);
196   } else {
197     options->Jo = options->Jop*PetscCosReal(options->ke*options->a)/(1.0-PetscCosReal(options->ke*options->a));
198   }
199   options->ky = PETSC_PI*options->m/options->b;
200   if (PetscPowReal(options->ky, 2) < options->Jop) {
201     options->kx = PetscSqrtScalar(options->Jop-PetscPowScalar(options->ky,2));
202     options->DeltaPrime = -2.0*options->kx*options->a*PetscCosReal(options->kx*options->a)/PetscSinReal(options->kx*options->a);
203   } else if (PetscPowReal(options->ky, 2) > options->Jop) {
204     options->kx = PetscSqrtScalar(PetscPowScalar(options->ky,2)-options->Jop);
205     options->DeltaPrime = -2.0*options->kx*options->a*PetscCoshReal(options->kx*options->a)/PetscSinhReal(options->kx*options->a);
206   } else { /*they're equal (or there's a NaN), lim(x*cot(x))_x->0=1*/
207     options->kx = 0;
208     options->DeltaPrime = -2.0;
209   }
210   ierr = PetscPrintf(comm, "DeltaPrime=%g\n",options->DeltaPrime);CHKERRQ(ierr);
211 
212   PetscFunctionReturn(0);
213 }
214 
215 static void f_n(PetscInt dim, PetscInt Nf, PetscInt NfAux,
216                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
217                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
218                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
219 {
220   const PetscScalar *pn = &u[uOff[DENSITY]];
221   *f0 = *pn;
222 }
223 
224 static PetscErrorCode PostStep(TS ts)
225 {
226   PetscErrorCode    ierr;
227   DM                dm;
228   AppCtx            *ctx;
229   PetscInt          stepi,num;
230   Vec               X;
231   PetscFunctionBegin;
232   ierr = TSGetApplicationContext(ts, &ctx);CHKERRQ(ierr);
233   if (ctx->debug<1) PetscFunctionReturn(0);
234   ierr = TSGetSolution(ts, &X);CHKERRQ(ierr);
235   ierr = VecGetDM(X, &dm);CHKERRQ(ierr);
236   ierr = TSGetStepNumber(ts, &stepi);CHKERRQ(ierr);
237   ierr = DMGetOutputSequenceNumber(dm, &num, NULL);CHKERRQ(ierr);
238   if (num < 0) {ierr = DMSetOutputSequenceNumber(dm, 0, 0.0);CHKERRQ(ierr);}
239   ierr = PetscObjectSetName((PetscObject) X, "u");CHKERRQ(ierr);
240   ierr = VecViewFromOptions(X, NULL, "-vec_view");CHKERRQ(ierr);
241   /* print integrals */
242   {
243     PetscDS          prob;
244     DM               plex;
245     PetscScalar den, tt[5];
246     ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
247     ierr = DMGetDS(plex, &prob);CHKERRQ(ierr);
248     ierr = PetscDSSetObjective(prob, 0, &f_n);CHKERRQ(ierr);
249     ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
250     den = tt[0];
251     ierr = DMDestroy(&plex);CHKERRQ(ierr);
252     PetscPrintf(PetscObjectComm((PetscObject)dm), "%D) total perturbed mass = %g\n", stepi, (double) PetscRealPart(den));CHKERRQ(ierr);
253   }
254   PetscFunctionReturn(0);
255 }
256 
257 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
258 {
259   PetscErrorCode ierr;
260 
261   PetscFunctionBeginUser;
262   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
263   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
264   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
265   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
266 
267   ierr = DMGetBoundingBox(*dm, ctx->lower, ctx->upper);CHKERRQ(ierr);
268   ctx->a = (ctx->upper[0] - ctx->lower[0])/2.0;
269   ctx->b = (ctx->upper[1] - ctx->lower[1])/2.0;
270   PetscFunctionReturn(0);
271 }
272 
273 static PetscErrorCode log_n_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
274 {
275   AppCtx *lctx = (AppCtx*)ctx;
276   u[0] = 2.*lctx->a + x[0];
277   return 0;
278 }
279 
280 static PetscErrorCode Omega_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
281 {
282   u[0] = 0.0;
283   return 0;
284 }
285 
286 static PetscErrorCode psi_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
287 {
288   AppCtx *lctx = (AppCtx*)ctx;
289   /* This sets up a symmetrix By flux aroound the mid point in x, which represents a current density flux along z.  The stability
290      is analytically known and reported in ProcessOptions. */
291   if (lctx->ke!=0.0) {
292     u[0] = (PetscCosReal(lctx->ke*(x[0]-lctx->a))-PetscCosReal(lctx->ke*lctx->a))/(1.0-PetscCosReal(lctx->ke*lctx->a));
293   } else {
294     u[0] = 1.0-PetscPowScalar((x[0]-lctx->a)/lctx->a,2);
295   }
296   return 0;
297 }
298 
299 static PetscErrorCode initialSolution_n(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
300 {
301   u[0] = 0.0;
302   return 0;
303 }
304 
305 static PetscErrorCode initialSolution_Omega(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
306 {
307   u[0] = 0.0;
308   return 0;
309 }
310 
311 static PetscErrorCode initialSolution_psi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *a_ctx)
312 {
313   AppCtx *ctx = (AppCtx*)a_ctx;
314   PetscScalar r = ctx->eps*(PetscScalar) (rand()) / (PetscScalar) (RAND_MAX);
315   if (x[0] == ctx->lower[0] || x[0] == ctx->upper[0]) r = 0;
316   u[0] = r;
317   /* PetscPrintf(PETSC_COMM_WORLD, "rand psi %lf\n",u[0]); */
318   return 0;
319 }
320 
321 static PetscErrorCode initialSolution_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
322 {
323   u[0] = 0.0;
324   return 0;
325 }
326 
327 static PetscErrorCode initialSolution_jz(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
328 {
329   u[0] = 0.0;
330   return 0;
331 }
332 
333 static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
334 {
335   PetscDS        ds;
336   DMLabel        label;
337   const PetscInt id = 1;
338   PetscErrorCode ierr, f;
339 
340   PetscFunctionBeginUser;
341   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
342   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
343   ierr = PetscDSSetResidual(ds, 0, f0_n,     f1_n);CHKERRQ(ierr);
344   ierr = PetscDSSetResidual(ds, 1, f0_Omega, f1_Omega);CHKERRQ(ierr);
345   ierr = PetscDSSetResidual(ds, 2, f0_psi,   f1_psi);CHKERRQ(ierr);
346   ierr = PetscDSSetResidual(ds, 3, f0_phi,   f1_phi);CHKERRQ(ierr);
347   ierr = PetscDSSetResidual(ds, 4, f0_jz,    f1_jz);CHKERRQ(ierr);
348   ctx->initialFuncs[0] = initialSolution_n;
349   ctx->initialFuncs[1] = initialSolution_Omega;
350   ctx->initialFuncs[2] = initialSolution_psi;
351   ctx->initialFuncs[3] = initialSolution_phi;
352   ctx->initialFuncs[4] = initialSolution_jz;
353   for (f = 0; f < 5; ++f) {
354     ierr = PetscDSSetImplicit(ds, f, ctx->implicit);CHKERRQ(ierr);
355     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, f, 0, NULL, (void (*)(void)) ctx->initialFuncs[f], NULL, ctx, NULL);CHKERRQ(ierr);
356   }
357   ierr = PetscDSSetContext(ds, 0, ctx);CHKERRQ(ierr);
358   PetscFunctionReturn(0);
359 }
360 
361 static PetscErrorCode SetupEquilibriumFields(DM dm, DM dmAux, AppCtx *ctx)
362 {
363   PetscErrorCode (*eqFuncs[3])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar [], void *) = {log_n_0, Omega_0, psi_0};
364   Vec            eq;
365   PetscErrorCode ierr;
366   AppCtx *ctxarr[3];
367 
368   ctxarr[0] = ctxarr[1] = ctxarr[2] = ctx; /* each variable could have a different context */
369   PetscFunctionBegin;
370   ierr = DMCreateLocalVector(dmAux, &eq);CHKERRQ(ierr);
371   ierr = DMProjectFunctionLocal(dmAux, 0.0, eqFuncs, (void **)ctxarr, INSERT_ALL_VALUES, eq);CHKERRQ(ierr);
372   ierr = DMSetAuxiliaryVec(dm, NULL, 0, eq);CHKERRQ(ierr);
373   if (ctx->plotRef) {  /* plot reference functions */
374     PetscViewer       viewer = NULL;
375     PetscBool         isHDF5,isVTK;
376     char              buf[256];
377     Vec               global;
378     PetscInt          dim;
379 
380     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
381     ierr = DMCreateGlobalVector(dmAux,&global);CHKERRQ(ierr);
382     ierr = VecSet(global,.0);CHKERRQ(ierr); /* BCs! */
383     ierr = DMLocalToGlobalBegin(dmAux,eq,INSERT_VALUES,global);CHKERRQ(ierr);
384     ierr = DMLocalToGlobalEnd(dmAux,eq,INSERT_VALUES,global);CHKERRQ(ierr);
385     ierr = PetscViewerCreate(PetscObjectComm((PetscObject)dmAux),&viewer);CHKERRQ(ierr);
386 #ifdef PETSC_HAVE_HDF5
387     ierr = PetscViewerSetType(viewer,PETSCVIEWERHDF5);CHKERRQ(ierr);
388 #else
389     ierr = PetscViewerSetType(viewer,PETSCVIEWERVTK);CHKERRQ(ierr);
390 #endif
391     ierr = PetscViewerSetFromOptions(viewer);CHKERRQ(ierr);
392     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5);CHKERRQ(ierr);
393     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK);CHKERRQ(ierr);
394     if (isHDF5) {
395       ierr = PetscSNPrintf(buf, 256, "uEquilibrium-%dD.h5", dim);CHKERRQ(ierr);
396     } else if (isVTK) {
397       ierr = PetscSNPrintf(buf, 256, "uEquilibrium-%dD.vtu", dim);CHKERRQ(ierr);
398       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);CHKERRQ(ierr);
399     }
400     ierr = PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
401     ierr = PetscViewerFileSetName(viewer,buf);CHKERRQ(ierr);
402     if (isHDF5) {ierr = DMView(dmAux,viewer);CHKERRQ(ierr);}
403     /* view equilibrium fields, this will overwrite fine grids with coarse grids! */
404     ierr = PetscObjectSetName((PetscObject) global, "u0");CHKERRQ(ierr);
405     ierr = VecView(global,viewer);CHKERRQ(ierr);
406     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
407     ierr = VecDestroy(&global);CHKERRQ(ierr);
408   }
409   ierr = VecDestroy(&eq);CHKERRQ(ierr);
410   PetscFunctionReturn(0);
411 }
412 
413 static PetscErrorCode SetupAuxDM(DM dm, PetscInt NfAux, PetscFE feAux[], AppCtx *user)
414 {
415   DM             dmAux, coordDM;
416   PetscInt       f;
417   PetscErrorCode ierr;
418 
419   PetscFunctionBegin;
420   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
421   ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
422   if (!feAux) PetscFunctionReturn(0);
423   ierr = DMClone(dm, &dmAux);CHKERRQ(ierr);
424   ierr = DMSetCoordinateDM(dmAux, coordDM);CHKERRQ(ierr);
425   for (f = 0; f < NfAux; ++f) {ierr = DMSetField(dmAux, f, NULL, (PetscObject) feAux[f]);CHKERRQ(ierr);}
426   ierr = DMCreateDS(dmAux);CHKERRQ(ierr);
427   ierr = SetupEquilibriumFields(dm, dmAux, user);CHKERRQ(ierr);
428   ierr = DMDestroy(&dmAux);CHKERRQ(ierr);
429   PetscFunctionReturn(0);
430 }
431 
432 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
433 {
434   DM              cdm = dm;
435   PetscFE         fe[5], feAux[3];
436   PetscInt        dim, Nf = 5, NfAux = 3, f;
437   PetscBool       simplex;
438   MPI_Comm        comm;
439   PetscErrorCode  ierr;
440 
441   PetscFunctionBeginUser;
442   /* Create finite element */
443   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
444   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
445   ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr);
446   ierr = PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[0]);CHKERRQ(ierr);
447   ierr = PetscObjectSetName((PetscObject) fe[0], "density");CHKERRQ(ierr);
448   ierr = PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[1]);CHKERRQ(ierr);
449   ierr = PetscObjectSetName((PetscObject) fe[1], "vorticity");CHKERRQ(ierr);
450   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
451   ierr = PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[2]);CHKERRQ(ierr);
452   ierr = PetscObjectSetName((PetscObject) fe[2], "flux");CHKERRQ(ierr);
453   ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr);
454   ierr = PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[3]);CHKERRQ(ierr);
455   ierr = PetscObjectSetName((PetscObject) fe[3], "potential");CHKERRQ(ierr);
456   ierr = PetscFECopyQuadrature(fe[0], fe[3]);CHKERRQ(ierr);
457   ierr = PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[4]);CHKERRQ(ierr);
458   ierr = PetscObjectSetName((PetscObject) fe[4], "current");CHKERRQ(ierr);
459   ierr = PetscFECopyQuadrature(fe[0], fe[4]);CHKERRQ(ierr);
460 
461   ierr = PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &feAux[0]);CHKERRQ(ierr);
462   ierr = PetscObjectSetName((PetscObject) feAux[0], "n_0");CHKERRQ(ierr);
463   ierr = PetscFECopyQuadrature(fe[0], feAux[0]);CHKERRQ(ierr);
464   ierr = PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &feAux[1]);CHKERRQ(ierr);
465   ierr = PetscObjectSetName((PetscObject) feAux[1], "vorticity_0");CHKERRQ(ierr);
466   ierr = PetscFECopyQuadrature(fe[0], feAux[1]);CHKERRQ(ierr);
467   ierr = PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &feAux[2]);CHKERRQ(ierr);
468   ierr = PetscObjectSetName((PetscObject) feAux[2], "flux_0");CHKERRQ(ierr);
469   ierr = PetscFECopyQuadrature(fe[0], feAux[2]);CHKERRQ(ierr);
470   /* Set discretization and boundary conditions for each mesh */
471   for (f = 0; f < Nf; ++f) {ierr = DMSetField(dm, f, NULL, (PetscObject) fe[f]);CHKERRQ(ierr);}
472   ierr = DMCreateDS(dm);CHKERRQ(ierr);
473   ierr = SetupProblem(dm, ctx);CHKERRQ(ierr);
474   while (cdm) {
475     ierr = SetupAuxDM(dm, NfAux, feAux, ctx);CHKERRQ(ierr);
476     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
477     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
478   }
479   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&fe[f]);CHKERRQ(ierr);}
480   for (f = 0; f < NfAux; ++f) {ierr = PetscFEDestroy(&feAux[f]);CHKERRQ(ierr);}
481   PetscFunctionReturn(0);
482 }
483 
484 int main(int argc, char **argv)
485 {
486   DM             dm;
487   TS             ts;
488   Vec            u, r;
489   AppCtx         ctx;
490   PetscReal      t       = 0.0;
491   PetscReal      L2error = 0.0;
492   PetscErrorCode ierr;
493   AppCtx        *ctxarr[5];
494 
495   ctxarr[0] = ctxarr[1] = ctxarr[2] = ctxarr[3] = ctxarr[4] = &ctx; /* each variable could have a different context */
496   s_ctx = &ctx;
497   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
498   ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr);
499   /* create mesh and problem */
500   ierr = CreateMesh(PETSC_COMM_WORLD, &ctx, &dm);CHKERRQ(ierr);
501   ierr = DMSetApplicationContext(dm, &ctx);CHKERRQ(ierr);
502   ierr = PetscMalloc1(5, &ctx.initialFuncs);CHKERRQ(ierr);
503   ierr = SetupDiscretization(dm, &ctx);CHKERRQ(ierr);
504   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
505   ierr = PetscObjectSetName((PetscObject) u, "u");CHKERRQ(ierr);
506   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
507   ierr = PetscObjectSetName((PetscObject) r, "r");CHKERRQ(ierr);
508   /* create TS */
509   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
510   ierr = TSSetDM(ts, dm);CHKERRQ(ierr);
511   ierr = TSSetApplicationContext(ts, &ctx);CHKERRQ(ierr);
512   ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx);CHKERRQ(ierr);
513   if (ctx.implicit) {
514     ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx);CHKERRQ(ierr);
515     ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx);CHKERRQ(ierr);
516   } else {
517     ierr = DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &ctx);CHKERRQ(ierr);
518   }
519   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
520   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
521   ierr = TSSetPostStep(ts, PostStep);CHKERRQ(ierr);
522   /* make solution & solve */
523   ierr = DMProjectFunction(dm, t, ctx.initialFuncs, (void **)ctxarr, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
524   ierr = TSSetSolution(ts,u);CHKERRQ(ierr);
525   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
526   ierr = PostStep(ts);CHKERRQ(ierr); /* print the initial state */
527   ierr = TSSolve(ts, u);CHKERRQ(ierr);
528   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
529   ierr = DMComputeL2Diff(dm, t, ctx.initialFuncs, (void **)ctxarr, u, &L2error);CHKERRQ(ierr);
530   if (L2error < 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
531   else                   {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", L2error);CHKERRQ(ierr);}
532   ierr = VecDestroy(&u);CHKERRQ(ierr);
533   ierr = VecDestroy(&r);CHKERRQ(ierr);
534   ierr = TSDestroy(&ts);CHKERRQ(ierr);
535   ierr = DMDestroy(&dm);CHKERRQ(ierr);
536   ierr = PetscFree(ctx.initialFuncs);CHKERRQ(ierr);
537   ierr = PetscFinalize();
538   return ierr;
539 }
540 
541 /*TEST
542 
543   test:
544     suffix: 0
545     args: -debug 1 -dm_refine 1 -dm_plex_simplex 0 -dm_plex_box_faces 3,3 -dm_plex_box_bd periodic,none -dm_plex_box_upper 2.0,6.283185307179586 \
546           -ts_max_steps 1 -ts_max_time 10. -ts_dt 1.0
547   test:
548     # Remapping with periodicity is broken
549     suffix: 1
550     args: -debug 1 -dm_plex_shape cylinder -dm_plex_dim 3 -dm_refine 1 -dm_refine_remap 0 -dm_plex_cylinder_bd periodic -dm_plex_boundary_label marker \
551            -ts_max_steps 1 -ts_max_time 10. -ts_dt 1.0
552 
553 TEST*/
554