xref: /petsc/src/ts/tutorials/ex48.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Evolution of magnetic islands.\n\
2*c4762a1bSJed Brown 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*c4762a1bSJed Brown 
4*c4762a1bSJed Brown /*F
5*c4762a1bSJed Brown 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*c4762a1bSJed Brown \begin{equation}
7*c4762a1bSJed Brown   \begin{aligned}
8*c4762a1bSJed Brown     \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*c4762a1bSJed Brown   \partial_t \tilde\Omega   &= \left\{ \tilde\Omega, \tilde\phi \right\} + \beta \left\{ j_z, \tilde\psi \right\} + \mu \nabla^2_\perp \tilde\Omega \\
10*c4762a1bSJed Brown   \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*c4762a1bSJed Brown   \nabla^2_\perp\tilde\phi        &= \tilde\Omega \\
12*c4762a1bSJed Brown   j_z  &= -\nabla^2_\perp  \left(\tilde\psi + \psi_0  \right)\\
13*c4762a1bSJed Brown   \end{aligned}
14*c4762a1bSJed Brown \end{equation}
15*c4762a1bSJed Brown F*/
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown #include <petscdmplex.h>
18*c4762a1bSJed Brown #include <petscts.h>
19*c4762a1bSJed Brown #include <petscds.h>
20*c4762a1bSJed Brown #include <assert.h>
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown typedef struct {
23*c4762a1bSJed Brown   PetscInt       debug;             /* The debugging level */
24*c4762a1bSJed Brown   PetscBool      plotRef;           /* Plot the reference fields */
25*c4762a1bSJed Brown   /* Domain and mesh definition */
26*c4762a1bSJed Brown   PetscInt       dim;               /* The topological mesh dimension */
27*c4762a1bSJed Brown   char           filename[2048];    /* The optional ExodusII file */
28*c4762a1bSJed Brown   PetscBool      cell_simplex;           /* Simplicial mesh */
29*c4762a1bSJed Brown   DMBoundaryType boundary_types[3];
30*c4762a1bSJed Brown   PetscInt       cells[3];
31*c4762a1bSJed Brown   PetscInt       refine;
32*c4762a1bSJed Brown   /* geometry  */
33*c4762a1bSJed Brown   PetscReal      domain_lo[3], domain_hi[3];
34*c4762a1bSJed Brown   DMBoundaryType periodicity[3];              /* The domain periodicity */
35*c4762a1bSJed Brown   PetscReal      b0[3]; /* not used */
36*c4762a1bSJed Brown   /* Problem definition */
37*c4762a1bSJed Brown   PetscErrorCode (**initialFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
38*c4762a1bSJed Brown   PetscReal      mu, eta, beta;
39*c4762a1bSJed Brown   PetscReal      a,b,Jo,Jop,m,ke,kx,ky,DeltaPrime,eps;
40*c4762a1bSJed Brown   /* solver */
41*c4762a1bSJed Brown   PetscBool      implicit;
42*c4762a1bSJed Brown } AppCtx;
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown static AppCtx *s_ctx;
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown static PetscScalar poissonBracket(PetscInt dim, const PetscScalar df[], const PetscScalar dg[])
47*c4762a1bSJed Brown {
48*c4762a1bSJed Brown   PetscScalar ret = df[0]*dg[1] - df[1]*dg[0];
49*c4762a1bSJed Brown   return ret;
50*c4762a1bSJed Brown }
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown enum field_idx {DENSITY,OMEGA,PSI,PHI,JZ};
53*c4762a1bSJed Brown 
54*c4762a1bSJed Brown static void f0_n(PetscInt dim, PetscInt Nf, PetscInt NfAux,
55*c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
56*c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
57*c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
58*c4762a1bSJed Brown {
59*c4762a1bSJed Brown   const PetscScalar *pnDer   = &u_x[uOff_x[DENSITY]];
60*c4762a1bSJed Brown   const PetscScalar *ppsiDer = &u_x[uOff_x[PSI]];
61*c4762a1bSJed Brown   const PetscScalar *pphiDer = &u_x[uOff_x[PHI]];
62*c4762a1bSJed Brown   const PetscScalar *jzDer   = &u_x[uOff_x[JZ]];
63*c4762a1bSJed Brown   const PetscScalar *logRefDenDer = &a_x[aOff_x[DENSITY]];
64*c4762a1bSJed Brown   f0[0] += - poissonBracket(dim,pnDer, pphiDer) - s_ctx->beta*poissonBracket(dim,jzDer, ppsiDer) - poissonBracket(dim,logRefDenDer, pphiDer);
65*c4762a1bSJed Brown   if (u_t) f0[0] += u_t[DENSITY];
66*c4762a1bSJed Brown }
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown static void f1_n(PetscInt dim, PetscInt Nf, PetscInt NfAux,
69*c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
70*c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
71*c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
72*c4762a1bSJed Brown {
73*c4762a1bSJed Brown   const PetscScalar *pnDer = &u_x[uOff_x[DENSITY]];
74*c4762a1bSJed Brown   PetscInt           d;
75*c4762a1bSJed Brown 
76*c4762a1bSJed Brown   for (d = 0; d < dim-1; ++d) f1[d] = -s_ctx->mu*pnDer[d];
77*c4762a1bSJed Brown }
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown static void f0_Omega(PetscInt dim, PetscInt Nf, PetscInt NfAux,
80*c4762a1bSJed Brown                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
81*c4762a1bSJed Brown                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
82*c4762a1bSJed Brown                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
83*c4762a1bSJed Brown {
84*c4762a1bSJed Brown   const PetscScalar *pOmegaDer = &u_x[uOff_x[OMEGA]];
85*c4762a1bSJed Brown   const PetscScalar *ppsiDer   = &u_x[uOff_x[PSI]];
86*c4762a1bSJed Brown   const PetscScalar *pphiDer   = &u_x[uOff_x[PHI]];
87*c4762a1bSJed Brown   const PetscScalar *jzDer     = &u_x[uOff_x[JZ]];
88*c4762a1bSJed Brown 
89*c4762a1bSJed Brown   f0[0] += - poissonBracket(dim,pOmegaDer, pphiDer) - s_ctx->beta*poissonBracket(dim,jzDer, ppsiDer);
90*c4762a1bSJed Brown   if (u_t) f0[0] += u_t[OMEGA];
91*c4762a1bSJed Brown }
92*c4762a1bSJed Brown 
93*c4762a1bSJed Brown static void f1_Omega(PetscInt dim, PetscInt Nf, PetscInt NfAux,
94*c4762a1bSJed Brown                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
95*c4762a1bSJed Brown                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
96*c4762a1bSJed Brown                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
97*c4762a1bSJed Brown {
98*c4762a1bSJed Brown   const PetscScalar *pOmegaDer = &u_x[uOff_x[OMEGA]];
99*c4762a1bSJed Brown   PetscInt           d;
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown   for (d = 0; d < dim-1; ++d) f1[d] = -s_ctx->mu*pOmegaDer[d];
102*c4762a1bSJed Brown }
103*c4762a1bSJed Brown 
104*c4762a1bSJed Brown static void f0_psi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
105*c4762a1bSJed Brown                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
106*c4762a1bSJed Brown                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
107*c4762a1bSJed Brown                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
108*c4762a1bSJed Brown {
109*c4762a1bSJed Brown   const PetscScalar *pnDer     = &u_x[uOff_x[DENSITY]];
110*c4762a1bSJed Brown   const PetscScalar *ppsiDer   = &u_x[uOff_x[PSI]];
111*c4762a1bSJed Brown   const PetscScalar *pphiDer   = &u_x[uOff_x[PHI]];
112*c4762a1bSJed Brown   const PetscScalar *refPsiDer = &a_x[aOff_x[PSI]];
113*c4762a1bSJed Brown   const PetscScalar *logRefDenDer= &a_x[aOff_x[DENSITY]];
114*c4762a1bSJed Brown   PetscScalar       psiDer[3];
115*c4762a1bSJed Brown   PetscScalar       phi_n_Der[3];
116*c4762a1bSJed Brown   PetscInt          d;
117*c4762a1bSJed Brown   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 */
118*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) {
119*c4762a1bSJed Brown     psiDer[d]    = refPsiDer[d] + ppsiDer[d];
120*c4762a1bSJed Brown     phi_n_Der[d] = pphiDer[d]   - pnDer[d];
121*c4762a1bSJed Brown   }
122*c4762a1bSJed Brown   f0[0] = - poissonBracket(dim,psiDer, phi_n_Der) + poissonBracket(dim,logRefDenDer, ppsiDer);
123*c4762a1bSJed Brown   if (u_t) f0[0] += u_t[PSI];
124*c4762a1bSJed Brown }
125*c4762a1bSJed Brown 
126*c4762a1bSJed Brown static void f1_psi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
127*c4762a1bSJed Brown                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
128*c4762a1bSJed Brown                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
129*c4762a1bSJed Brown                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
130*c4762a1bSJed Brown {
131*c4762a1bSJed Brown   const PetscScalar *ppsi = &u_x[uOff_x[PSI]];
132*c4762a1bSJed Brown   PetscInt           d;
133*c4762a1bSJed Brown 
134*c4762a1bSJed Brown   for (d = 0; d < dim-1; ++d) f1[d] = -(s_ctx->eta/s_ctx->beta)*ppsi[d];
135*c4762a1bSJed Brown }
136*c4762a1bSJed Brown 
137*c4762a1bSJed Brown static void f0_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
138*c4762a1bSJed Brown                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
139*c4762a1bSJed Brown                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
140*c4762a1bSJed Brown                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
141*c4762a1bSJed Brown {
142*c4762a1bSJed Brown   f0[0] = -u[uOff[OMEGA]];
143*c4762a1bSJed Brown }
144*c4762a1bSJed Brown 
145*c4762a1bSJed Brown static void f1_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
146*c4762a1bSJed Brown                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
147*c4762a1bSJed Brown                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
148*c4762a1bSJed Brown                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
149*c4762a1bSJed Brown {
150*c4762a1bSJed Brown   const PetscScalar *pphi = &u_x[uOff_x[PHI]];
151*c4762a1bSJed Brown   PetscInt           d;
152*c4762a1bSJed Brown 
153*c4762a1bSJed Brown   for (d = 0; d < dim-1; ++d) f1[d] = pphi[d];
154*c4762a1bSJed Brown }
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown static void f0_jz(PetscInt dim, PetscInt Nf, PetscInt NfAux,
157*c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
158*c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
159*c4762a1bSJed Brown                   PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
160*c4762a1bSJed Brown {
161*c4762a1bSJed Brown   f0[0] = u[uOff[JZ]];
162*c4762a1bSJed Brown }
163*c4762a1bSJed Brown 
164*c4762a1bSJed Brown static void f1_jz(PetscInt dim, PetscInt Nf, PetscInt NfAux,
165*c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
166*c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
167*c4762a1bSJed Brown                   PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
168*c4762a1bSJed Brown {
169*c4762a1bSJed Brown   const PetscScalar *ppsi = &u_x[uOff_x[PSI]];
170*c4762a1bSJed Brown   const PetscScalar *refPsiDer = &a_x[aOff_x[PSI]]; /* aOff_x[PSI] == 2*PSI */
171*c4762a1bSJed Brown   PetscInt           d;
172*c4762a1bSJed Brown 
173*c4762a1bSJed Brown   for (d = 0; d < dim-1; ++d) f1[d] = ppsi[d] + refPsiDer[d];
174*c4762a1bSJed Brown }
175*c4762a1bSJed Brown 
176*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
177*c4762a1bSJed Brown {
178*c4762a1bSJed Brown   PetscBool      flg;
179*c4762a1bSJed Brown   PetscErrorCode ierr;
180*c4762a1bSJed Brown   PetscInt       ii, bd;
181*c4762a1bSJed Brown   PetscFunctionBeginUser;
182*c4762a1bSJed Brown   options->debug               = 1;
183*c4762a1bSJed Brown   options->plotRef             = PETSC_FALSE;
184*c4762a1bSJed Brown   options->dim                 = 2;
185*c4762a1bSJed Brown   options->filename[0]         = '\0';
186*c4762a1bSJed Brown   options->cell_simplex        = PETSC_FALSE;
187*c4762a1bSJed Brown   options->implicit            = PETSC_FALSE;
188*c4762a1bSJed Brown   options->refine              = 2;
189*c4762a1bSJed Brown   options->domain_lo[0]  = 0.0;
190*c4762a1bSJed Brown   options->domain_lo[1]  = 0.0;
191*c4762a1bSJed Brown   options->domain_lo[2]  = 0.0;
192*c4762a1bSJed Brown   options->domain_hi[0]  = 2.0;
193*c4762a1bSJed Brown   options->domain_hi[1]  = 2.0*PETSC_PI;
194*c4762a1bSJed Brown   options->domain_hi[2]  = 2.0;
195*c4762a1bSJed Brown   options->periodicity[0]    = DM_BOUNDARY_NONE;
196*c4762a1bSJed Brown   options->periodicity[1]    = DM_BOUNDARY_NONE;
197*c4762a1bSJed Brown   options->periodicity[2]    = DM_BOUNDARY_NONE;
198*c4762a1bSJed Brown   options->mu   = 0;
199*c4762a1bSJed Brown   options->eta  = 0;
200*c4762a1bSJed Brown   options->beta = 1;
201*c4762a1bSJed Brown   options->a = 1;
202*c4762a1bSJed Brown   options->b = PETSC_PI;
203*c4762a1bSJed Brown   options->Jop = 0;
204*c4762a1bSJed Brown   options->m = 1;
205*c4762a1bSJed Brown   options->eps = 1.e-6;
206*c4762a1bSJed Brown 
207*c4762a1bSJed Brown   for (ii = 0; ii < options->dim; ++ii) options->cells[ii] = 4;
208*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
209*c4762a1bSJed Brown   ierr = PetscOptionsInt("-debug", "The debugging level", "ex48.c", options->debug, &options->debug, NULL);CHKERRQ(ierr);
210*c4762a1bSJed Brown   ierr = PetscOptionsBool("-plot_ref", "Plot the reference fields", "ex48.c", options->plotRef, &options->plotRef, NULL);CHKERRQ(ierr);
211*c4762a1bSJed Brown   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex48.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
212*c4762a1bSJed Brown   if (options->dim < 2 || options->dim > 3) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Dim %D must be 2 or 3",options->dim);CHKERRQ(ierr);
213*c4762a1bSJed Brown   ierr = PetscOptionsInt("-dm_refine", "Hack to get refinement level for cylinder", "ex48.c", options->refine, &options->refine, NULL);CHKERRQ(ierr);
214*c4762a1bSJed Brown   ierr = PetscOptionsReal("-mu", "mu", "ex48.c", options->mu, &options->mu, NULL);CHKERRQ(ierr);
215*c4762a1bSJed Brown   ierr = PetscOptionsReal("-eta", "eta", "ex48.c", options->eta, &options->eta, NULL);CHKERRQ(ierr);
216*c4762a1bSJed Brown   ierr = PetscOptionsReal("-beta", "beta", "ex48.c", options->beta, &options->beta, NULL);CHKERRQ(ierr);
217*c4762a1bSJed Brown   ierr = PetscOptionsReal("-Jop", "Jop", "ex48.c", options->Jop, &options->Jop, NULL);CHKERRQ(ierr);
218*c4762a1bSJed Brown   ierr = PetscOptionsReal("-m", "m", "ex48.c", options->m, &options->m, NULL);CHKERRQ(ierr);
219*c4762a1bSJed Brown   ierr = PetscOptionsReal("-eps", "eps", "ex48.c", options->eps, &options->eps, NULL);CHKERRQ(ierr);
220*c4762a1bSJed Brown   ierr = PetscOptionsString("-f", "Exodus.II filename to read", "ex48.c", options->filename, options->filename, sizeof(options->filename), &flg);CHKERRQ(ierr);
221*c4762a1bSJed Brown   ierr = PetscOptionsBool("-cell_simplex", "Simplicial (true) or tensor (false) mesh", "ex48.c", options->cell_simplex, &options->cell_simplex, NULL);CHKERRQ(ierr);
222*c4762a1bSJed Brown   ierr = PetscOptionsBool("-implicit", "Use implicit time integrator", "ex48.c", options->implicit, &options->implicit, NULL);CHKERRQ(ierr);
223*c4762a1bSJed Brown   ii = options->dim;
224*c4762a1bSJed Brown   ierr = PetscOptionsRealArray("-domain_hi", "Domain size", "ex48.c", options->domain_hi, &ii, NULL);CHKERRQ(ierr);
225*c4762a1bSJed Brown   ii = options->dim;
226*c4762a1bSJed Brown   ierr = PetscOptionsRealArray("-domain_lo", "Domain size", "ex48.c", options->domain_lo, &ii, NULL);CHKERRQ(ierr);
227*c4762a1bSJed Brown   ii = options->dim;
228*c4762a1bSJed Brown   bd = options->periodicity[0];
229*c4762a1bSJed Brown   ierr = PetscOptionsEList("-x_periodicity", "The x-boundary periodicity", "ex48.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[0]], &bd, NULL);CHKERRQ(ierr);
230*c4762a1bSJed Brown   options->periodicity[0] = (DMBoundaryType) bd;
231*c4762a1bSJed Brown   bd = options->periodicity[1];
232*c4762a1bSJed Brown   ierr = PetscOptionsEList("-y_periodicity", "The y-boundary periodicity", "ex48.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[1]], &bd, NULL);CHKERRQ(ierr);
233*c4762a1bSJed Brown   options->periodicity[1] = (DMBoundaryType) bd;
234*c4762a1bSJed Brown   bd = options->periodicity[2];
235*c4762a1bSJed Brown   ierr = PetscOptionsEList("-z_periodicity", "The z-boundary periodicity", "ex48.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[2]], &bd, NULL);CHKERRQ(ierr);
236*c4762a1bSJed Brown   options->periodicity[2] = (DMBoundaryType) bd;
237*c4762a1bSJed Brown   ii = options->dim;
238*c4762a1bSJed Brown   ierr = PetscOptionsIntArray("-cells", "Number of cells in each dimension", "ex48.c", options->cells, &ii, NULL);CHKERRQ(ierr);
239*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
240*c4762a1bSJed Brown   options->a = (options->domain_hi[0]-options->domain_lo[0])/2.0;
241*c4762a1bSJed Brown   options->b = (options->domain_hi[1]-options->domain_lo[1])/2.0;
242*c4762a1bSJed Brown   for (ii = 0; ii < options->dim; ++ii) {
243*c4762a1bSJed Brown     if (options->domain_hi[ii] <= options->domain_lo[ii]) SETERRQ3(comm,PETSC_ERR_ARG_WRONG,"Domain %D lo=%g hi=%g",ii,options->domain_lo[ii],options->domain_hi[ii]);
244*c4762a1bSJed Brown   }
245*c4762a1bSJed Brown   options->ke = PetscSqrtScalar(options->Jop);
246*c4762a1bSJed Brown   if (options->Jop==0.0) {
247*c4762a1bSJed Brown     options->Jo = 1.0/PetscPowScalar(options->a,2);
248*c4762a1bSJed Brown   } else {
249*c4762a1bSJed Brown     options->Jo = options->Jop*PetscCosReal(options->ke*options->a)/(1.0-PetscCosReal(options->ke*options->a));
250*c4762a1bSJed Brown   }
251*c4762a1bSJed Brown   options->ky = PETSC_PI*options->m/options->b;
252*c4762a1bSJed Brown   if (PetscPowReal(options->ky, 2) < options->Jop) {
253*c4762a1bSJed Brown     options->kx = PetscSqrtScalar(options->Jop-PetscPowScalar(options->ky,2));
254*c4762a1bSJed Brown     options->DeltaPrime = -2.0*options->kx*options->a*PetscCosReal(options->kx*options->a)/PetscSinReal(options->kx*options->a);
255*c4762a1bSJed Brown   } else if (PetscPowReal(options->ky, 2) > options->Jop) {
256*c4762a1bSJed Brown     options->kx = PetscSqrtScalar(PetscPowScalar(options->ky,2)-options->Jop);
257*c4762a1bSJed Brown     options->DeltaPrime = -2.0*options->kx*options->a*PetscCoshReal(options->kx*options->a)/PetscSinhReal(options->kx*options->a);
258*c4762a1bSJed Brown   } else { /*they're equal (or there's a NaN), lim(x*cot(x))_x->0=1*/
259*c4762a1bSJed Brown     options->kx = 0;
260*c4762a1bSJed Brown     options->DeltaPrime = -2.0;
261*c4762a1bSJed Brown   }
262*c4762a1bSJed Brown   ierr = PetscPrintf(comm, "DeltaPrime=%g\n",options->DeltaPrime);CHKERRQ(ierr);
263*c4762a1bSJed Brown 
264*c4762a1bSJed Brown   PetscFunctionReturn(0);
265*c4762a1bSJed Brown }
266*c4762a1bSJed Brown 
267*c4762a1bSJed Brown static void f_n(PetscInt dim, PetscInt Nf, PetscInt NfAux,
268*c4762a1bSJed Brown                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
269*c4762a1bSJed Brown                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
270*c4762a1bSJed Brown                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
271*c4762a1bSJed Brown {
272*c4762a1bSJed Brown   const PetscScalar *pn = &u[uOff[DENSITY]];
273*c4762a1bSJed Brown   *f0 = *pn;
274*c4762a1bSJed Brown }
275*c4762a1bSJed Brown 
276*c4762a1bSJed Brown static PetscErrorCode PostStep(TS ts)
277*c4762a1bSJed Brown {
278*c4762a1bSJed Brown   PetscErrorCode    ierr;
279*c4762a1bSJed Brown   DM                dm;
280*c4762a1bSJed Brown   AppCtx            *ctx;
281*c4762a1bSJed Brown   PetscInt          stepi,num;
282*c4762a1bSJed Brown   Vec               X;
283*c4762a1bSJed Brown   PetscFunctionBegin;
284*c4762a1bSJed Brown   ierr = TSGetApplicationContext(ts, &ctx);CHKERRQ(ierr); assert(ctx);
285*c4762a1bSJed Brown   if (ctx->debug<1) PetscFunctionReturn(0);
286*c4762a1bSJed Brown   ierr = TSGetSolution(ts, &X);CHKERRQ(ierr);
287*c4762a1bSJed Brown   ierr = VecGetDM(X, &dm);CHKERRQ(ierr);
288*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts, &stepi);CHKERRQ(ierr);
289*c4762a1bSJed Brown   ierr = DMGetOutputSequenceNumber(dm, &num, NULL);CHKERRQ(ierr);
290*c4762a1bSJed Brown   if (num < 0) {ierr = DMSetOutputSequenceNumber(dm, 0, 0.0);CHKERRQ(ierr);}
291*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) X, "u");CHKERRQ(ierr);
292*c4762a1bSJed Brown   ierr = VecViewFromOptions(X, NULL, "-vec_view");CHKERRQ(ierr);
293*c4762a1bSJed Brown   /* print integrals */
294*c4762a1bSJed Brown   {
295*c4762a1bSJed Brown     PetscDS          prob;
296*c4762a1bSJed Brown     DM               plex;
297*c4762a1bSJed Brown     PetscScalar den, tt[5];
298*c4762a1bSJed Brown     ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
299*c4762a1bSJed Brown     ierr = DMGetDS(plex, &prob);CHKERRQ(ierr);
300*c4762a1bSJed Brown     ierr = PetscDSSetObjective(prob, 0, &f_n);CHKERRQ(ierr);
301*c4762a1bSJed Brown     ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
302*c4762a1bSJed Brown     den = tt[0];
303*c4762a1bSJed Brown     ierr = DMDestroy(&plex);CHKERRQ(ierr);
304*c4762a1bSJed Brown     PetscPrintf(PetscObjectComm((PetscObject)dm), "%D) total perturbed mass = %g\n", stepi, (double) PetscRealPart(den));CHKERRQ(ierr);
305*c4762a1bSJed Brown   }
306*c4762a1bSJed Brown   PetscFunctionReturn(0);
307*c4762a1bSJed Brown }
308*c4762a1bSJed Brown 
309*c4762a1bSJed Brown static PetscErrorCode CreateBCLabel(DM dm, const char name[])
310*c4762a1bSJed Brown {
311*c4762a1bSJed Brown   DMLabel        label;
312*c4762a1bSJed Brown   PetscErrorCode ierr;
313*c4762a1bSJed Brown   PetscFunctionBeginUser;
314*c4762a1bSJed Brown   ierr = DMCreateLabel(dm, name);CHKERRQ(ierr);
315*c4762a1bSJed Brown   ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
316*c4762a1bSJed Brown   ierr = DMPlexMarkBoundaryFaces(dm, 1, label);CHKERRQ(ierr);
317*c4762a1bSJed Brown   ierr = DMPlexLabelComplete(dm, label);CHKERRQ(ierr);
318*c4762a1bSJed Brown   PetscFunctionReturn(0);
319*c4762a1bSJed Brown }
320*c4762a1bSJed Brown 
321*c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
322*c4762a1bSJed Brown {
323*c4762a1bSJed Brown   PetscInt       dim      = ctx->dim;
324*c4762a1bSJed Brown   const char    *filename = ctx->filename;
325*c4762a1bSJed Brown   size_t         len;
326*c4762a1bSJed Brown   PetscMPIInt    numProcs;
327*c4762a1bSJed Brown   PetscErrorCode ierr;
328*c4762a1bSJed Brown 
329*c4762a1bSJed Brown   PetscFunctionBeginUser;
330*c4762a1bSJed Brown   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
331*c4762a1bSJed Brown   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
332*c4762a1bSJed Brown   if (len) {
333*c4762a1bSJed Brown     ierr = DMPlexCreateFromFile(comm, filename, PETSC_TRUE, dm);CHKERRQ(ierr);
334*c4762a1bSJed Brown   } else {
335*c4762a1bSJed Brown     PetscInt        d;
336*c4762a1bSJed Brown 
337*c4762a1bSJed Brown     /* create DM */
338*c4762a1bSJed Brown     if (ctx->cell_simplex && dim == 3) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot mesh a cylinder with simplices");
339*c4762a1bSJed Brown     if (dim==2) {
340*c4762a1bSJed Brown       PetscInt refineRatio, totCells = 1;
341*c4762a1bSJed Brown       if (ctx->cell_simplex) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot mesh 2D with simplices");
342*c4762a1bSJed Brown       refineRatio = PetscMax((PetscInt) (PetscPowReal(numProcs, 1.0/dim) + 0.1) - 1, 1);
343*c4762a1bSJed Brown       for (d = 0; d < dim; ++d) {
344*c4762a1bSJed Brown         if (ctx->cells[d] < refineRatio) ctx->cells[d] = refineRatio;
345*c4762a1bSJed Brown         if (ctx->periodicity[d]==DM_BOUNDARY_PERIODIC && ctx->cells[d]*refineRatio <= 2) refineRatio = 2;
346*c4762a1bSJed Brown       }
347*c4762a1bSJed Brown       for (d = 0; d < dim; ++d) {
348*c4762a1bSJed Brown         ctx->cells[d] *= refineRatio;
349*c4762a1bSJed Brown         totCells *= ctx->cells[d];
350*c4762a1bSJed Brown       }
351*c4762a1bSJed Brown       if (totCells % numProcs) SETERRQ2(comm,PETSC_ERR_ARG_WRONG,"Total cells %D not divisible by processes %D", totCells, numProcs);
352*c4762a1bSJed Brown       ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, ctx->cells, ctx->domain_lo, ctx->domain_hi, ctx->periodicity, PETSC_TRUE, dm);CHKERRQ(ierr);
353*c4762a1bSJed Brown     } else {
354*c4762a1bSJed Brown       if (ctx->periodicity[0]==DM_BOUNDARY_PERIODIC || ctx->periodicity[1]==DM_BOUNDARY_PERIODIC) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot do periodic in x or y in a cylinder");
355*c4762a1bSJed Brown       /* we stole dm_refine so clear it */
356*c4762a1bSJed Brown       ierr = PetscOptionsClearValue(NULL,"-dm_refine");CHKERRQ(ierr);
357*c4762a1bSJed Brown       ierr = DMPlexCreateHexCylinderMesh(comm, ctx->refine, ctx->periodicity[2], dm);CHKERRQ(ierr);
358*c4762a1bSJed Brown     }
359*c4762a1bSJed Brown   }
360*c4762a1bSJed Brown   {
361*c4762a1bSJed Brown     DM distributedMesh = NULL;
362*c4762a1bSJed Brown     /* Distribute mesh over processes */
363*c4762a1bSJed Brown     ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr);
364*c4762a1bSJed Brown     if (distributedMesh) {
365*c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
366*c4762a1bSJed Brown       *dm  = distributedMesh;
367*c4762a1bSJed Brown     }
368*c4762a1bSJed Brown   }
369*c4762a1bSJed Brown   {
370*c4762a1bSJed Brown     PetscBool hasLabel;
371*c4762a1bSJed Brown     ierr = DMHasLabel(*dm, "marker", &hasLabel);CHKERRQ(ierr);
372*c4762a1bSJed Brown     if (!hasLabel) {ierr = CreateBCLabel(*dm, "marker");CHKERRQ(ierr);}
373*c4762a1bSJed Brown   }
374*c4762a1bSJed Brown   {
375*c4762a1bSJed Brown     char      convType[256];
376*c4762a1bSJed Brown     PetscBool flg;
377*c4762a1bSJed Brown     ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr);
378*c4762a1bSJed Brown     ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex48",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr);
379*c4762a1bSJed Brown     ierr = PetscOptionsEnd();CHKERRQ(ierr);
380*c4762a1bSJed Brown     if (flg) {
381*c4762a1bSJed Brown       DM dmConv;
382*c4762a1bSJed Brown       ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr);
383*c4762a1bSJed Brown       if (dmConv) {
384*c4762a1bSJed Brown         ierr = DMDestroy(dm);CHKERRQ(ierr);
385*c4762a1bSJed Brown         *dm  = dmConv;
386*c4762a1bSJed Brown       }
387*c4762a1bSJed Brown     }
388*c4762a1bSJed Brown   }
389*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
390*c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
391*c4762a1bSJed Brown   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); /* needed for periodic */
392*c4762a1bSJed Brown   PetscFunctionReturn(0);
393*c4762a1bSJed Brown }
394*c4762a1bSJed Brown 
395*c4762a1bSJed Brown static PetscErrorCode log_n_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
396*c4762a1bSJed Brown {
397*c4762a1bSJed Brown   AppCtx *lctx = (AppCtx*)ctx;
398*c4762a1bSJed Brown   assert(ctx);
399*c4762a1bSJed Brown   u[0] = (lctx->domain_hi-lctx->domain_lo)+x[0];
400*c4762a1bSJed Brown   return 0;
401*c4762a1bSJed Brown }
402*c4762a1bSJed Brown 
403*c4762a1bSJed Brown static PetscErrorCode Omega_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
404*c4762a1bSJed Brown {
405*c4762a1bSJed Brown   u[0] = 0.0;
406*c4762a1bSJed Brown   return 0;
407*c4762a1bSJed Brown }
408*c4762a1bSJed Brown 
409*c4762a1bSJed Brown static PetscErrorCode psi_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
410*c4762a1bSJed Brown {
411*c4762a1bSJed Brown   AppCtx *lctx = (AppCtx*)ctx;
412*c4762a1bSJed Brown   assert(ctx);
413*c4762a1bSJed Brown   /* This sets up a symmetrix By flux aroound the mid point in x, which represents a current density flux along z.  The stability
414*c4762a1bSJed Brown      is analytically known and reported in ProcessOptions. */
415*c4762a1bSJed Brown   if (lctx->ke!=0.0) {
416*c4762a1bSJed Brown     u[0] = (PetscCosReal(lctx->ke*(x[0]-lctx->a))-PetscCosReal(lctx->ke*lctx->a))/(1.0-PetscCosReal(lctx->ke*lctx->a));
417*c4762a1bSJed Brown   } else {
418*c4762a1bSJed Brown     u[0] = 1.0-PetscPowScalar((x[0]-lctx->a)/lctx->a,2);
419*c4762a1bSJed Brown   }
420*c4762a1bSJed Brown   return 0;
421*c4762a1bSJed Brown }
422*c4762a1bSJed Brown 
423*c4762a1bSJed Brown static PetscErrorCode initialSolution_n(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
424*c4762a1bSJed Brown {
425*c4762a1bSJed Brown   u[0] = 0.0;
426*c4762a1bSJed Brown   return 0;
427*c4762a1bSJed Brown }
428*c4762a1bSJed Brown 
429*c4762a1bSJed Brown static PetscErrorCode initialSolution_Omega(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
430*c4762a1bSJed Brown {
431*c4762a1bSJed Brown   u[0] = 0.0;
432*c4762a1bSJed Brown   return 0;
433*c4762a1bSJed Brown }
434*c4762a1bSJed Brown 
435*c4762a1bSJed Brown static PetscErrorCode initialSolution_psi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *a_ctx)
436*c4762a1bSJed Brown {
437*c4762a1bSJed Brown   AppCtx *ctx = (AppCtx*)a_ctx;
438*c4762a1bSJed Brown   PetscScalar r = ctx->eps*(PetscScalar) (rand()) / (PetscScalar) (RAND_MAX);
439*c4762a1bSJed Brown   assert(ctx);
440*c4762a1bSJed Brown   if (x[0] == ctx->domain_lo[0] || x[0] == ctx->domain_hi[0]) r = 0;
441*c4762a1bSJed Brown   u[0] = r;
442*c4762a1bSJed Brown   /* PetscPrintf(PETSC_COMM_WORLD, "rand psi %lf\n",u[0]); */
443*c4762a1bSJed Brown   return 0;
444*c4762a1bSJed Brown }
445*c4762a1bSJed Brown 
446*c4762a1bSJed Brown static PetscErrorCode initialSolution_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
447*c4762a1bSJed Brown {
448*c4762a1bSJed Brown   u[0] = 0.0;
449*c4762a1bSJed Brown   return 0;
450*c4762a1bSJed Brown }
451*c4762a1bSJed Brown 
452*c4762a1bSJed Brown static PetscErrorCode initialSolution_jz(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
453*c4762a1bSJed Brown {
454*c4762a1bSJed Brown   u[0] = 0.0;
455*c4762a1bSJed Brown   return 0;
456*c4762a1bSJed Brown }
457*c4762a1bSJed Brown 
458*c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
459*c4762a1bSJed Brown {
460*c4762a1bSJed Brown   PetscDS        prob;
461*c4762a1bSJed Brown   const PetscInt id = 1;
462*c4762a1bSJed Brown   PetscErrorCode ierr, f;
463*c4762a1bSJed Brown 
464*c4762a1bSJed Brown   PetscFunctionBeginUser;
465*c4762a1bSJed Brown   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
466*c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 0, f0_n,     f1_n);CHKERRQ(ierr);
467*c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 1, f0_Omega, f1_Omega);CHKERRQ(ierr);
468*c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 2, f0_psi,   f1_psi);CHKERRQ(ierr);
469*c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 3, f0_phi,   f1_phi);CHKERRQ(ierr);
470*c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 4, f0_jz,    f1_jz);CHKERRQ(ierr);
471*c4762a1bSJed Brown   ctx->initialFuncs[0] = initialSolution_n;
472*c4762a1bSJed Brown   ctx->initialFuncs[1] = initialSolution_Omega;
473*c4762a1bSJed Brown   ctx->initialFuncs[2] = initialSolution_psi;
474*c4762a1bSJed Brown   ctx->initialFuncs[3] = initialSolution_phi;
475*c4762a1bSJed Brown   ctx->initialFuncs[4] = initialSolution_jz;
476*c4762a1bSJed Brown   for (f = 0; f < 5; ++f) {
477*c4762a1bSJed Brown     ierr = PetscDSSetImplicit( prob, f, ctx->implicit);CHKERRQ(ierr);
478*c4762a1bSJed Brown     ierr = PetscDSAddBoundary( prob, DM_BC_ESSENTIAL, "wall", "marker", f, 0, NULL, (void (*)(void)) ctx->initialFuncs[f], 1, &id, ctx);CHKERRQ(ierr);
479*c4762a1bSJed Brown   }
480*c4762a1bSJed Brown   ierr = PetscDSSetContext(prob, 0, ctx);CHKERRQ(ierr);
481*c4762a1bSJed Brown   PetscFunctionReturn(0);
482*c4762a1bSJed Brown }
483*c4762a1bSJed Brown 
484*c4762a1bSJed Brown static PetscErrorCode SetupEquilibriumFields(DM dm, DM dmAux, AppCtx *ctx)
485*c4762a1bSJed Brown {
486*c4762a1bSJed Brown   PetscErrorCode (*eqFuncs[3])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar [], void *) = {log_n_0, Omega_0, psi_0};
487*c4762a1bSJed Brown   Vec            eq;
488*c4762a1bSJed Brown   PetscErrorCode ierr;
489*c4762a1bSJed Brown   AppCtx *ctxarr[3];
490*c4762a1bSJed Brown 
491*c4762a1bSJed Brown   ctxarr[0] = ctxarr[1] = ctxarr[2] = ctx; /* each variable could have a different context */
492*c4762a1bSJed Brown   PetscFunctionBegin;
493*c4762a1bSJed Brown   ierr = DMCreateLocalVector(dmAux, &eq);CHKERRQ(ierr);
494*c4762a1bSJed Brown   ierr = DMProjectFunctionLocal(dmAux, 0.0, eqFuncs, (void **)ctxarr, INSERT_ALL_VALUES, eq);CHKERRQ(ierr);
495*c4762a1bSJed Brown   ierr = PetscObjectCompose((PetscObject) dm, "A", (PetscObject) eq);CHKERRQ(ierr);
496*c4762a1bSJed Brown   if (ctx->plotRef) {  /* plot reference functions */
497*c4762a1bSJed Brown     PetscViewer       viewer = NULL;
498*c4762a1bSJed Brown     PetscBool         isHDF5,isVTK;
499*c4762a1bSJed Brown     char              buf[256];
500*c4762a1bSJed Brown     Vec               global;
501*c4762a1bSJed Brown     ierr = DMCreateGlobalVector(dmAux,&global);CHKERRQ(ierr);
502*c4762a1bSJed Brown     ierr = VecSet(global,.0);CHKERRQ(ierr); /* BCs! */
503*c4762a1bSJed Brown     ierr = DMLocalToGlobalBegin(dmAux,eq,INSERT_VALUES,global);CHKERRQ(ierr);
504*c4762a1bSJed Brown     ierr = DMLocalToGlobalEnd(dmAux,eq,INSERT_VALUES,global);CHKERRQ(ierr);
505*c4762a1bSJed Brown     ierr = PetscViewerCreate(PetscObjectComm((PetscObject)dmAux),&viewer);CHKERRQ(ierr);
506*c4762a1bSJed Brown #ifdef PETSC_HAVE_HDF5
507*c4762a1bSJed Brown     ierr = PetscViewerSetType(viewer,PETSCVIEWERHDF5);CHKERRQ(ierr);
508*c4762a1bSJed Brown #else
509*c4762a1bSJed Brown     ierr = PetscViewerSetType(viewer,PETSCVIEWERVTK);CHKERRQ(ierr);
510*c4762a1bSJed Brown #endif
511*c4762a1bSJed Brown     ierr = PetscViewerSetFromOptions(viewer);CHKERRQ(ierr);
512*c4762a1bSJed Brown     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5);CHKERRQ(ierr);
513*c4762a1bSJed Brown     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK);CHKERRQ(ierr);
514*c4762a1bSJed Brown     if (isHDF5) {
515*c4762a1bSJed Brown       ierr = PetscSNPrintf(buf, 256, "uEquilibrium-%dD.h5", ctx->dim);CHKERRQ(ierr);
516*c4762a1bSJed Brown     } else if (isVTK) {
517*c4762a1bSJed Brown       ierr = PetscSNPrintf(buf, 256, "uEquilibrium-%dD.vtu", ctx->dim);CHKERRQ(ierr);
518*c4762a1bSJed Brown       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);CHKERRQ(ierr);
519*c4762a1bSJed Brown     }
520*c4762a1bSJed Brown     ierr = PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
521*c4762a1bSJed Brown     ierr = PetscViewerFileSetName(viewer,buf);CHKERRQ(ierr);
522*c4762a1bSJed Brown     if (isHDF5) {ierr = DMView(dmAux,viewer);CHKERRQ(ierr);}
523*c4762a1bSJed Brown     /* view equilibrium fields, this will overwrite fine grids with coarse grids! */
524*c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) global, "u0");CHKERRQ(ierr);
525*c4762a1bSJed Brown     ierr = VecView(global,viewer);CHKERRQ(ierr);
526*c4762a1bSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
527*c4762a1bSJed Brown     ierr = VecDestroy(&global);CHKERRQ(ierr);
528*c4762a1bSJed Brown   }
529*c4762a1bSJed Brown   ierr = VecDestroy(&eq);CHKERRQ(ierr);
530*c4762a1bSJed Brown   PetscFunctionReturn(0);
531*c4762a1bSJed Brown }
532*c4762a1bSJed Brown 
533*c4762a1bSJed Brown static PetscErrorCode SetupAuxDM(DM dm, PetscInt NfAux, PetscFE feAux[], AppCtx *user)
534*c4762a1bSJed Brown {
535*c4762a1bSJed Brown   DM             dmAux, coordDM;
536*c4762a1bSJed Brown   PetscInt       f;
537*c4762a1bSJed Brown   PetscErrorCode ierr;
538*c4762a1bSJed Brown 
539*c4762a1bSJed Brown   PetscFunctionBegin;
540*c4762a1bSJed Brown   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
541*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
542*c4762a1bSJed Brown   if (!feAux) PetscFunctionReturn(0);
543*c4762a1bSJed Brown   ierr = DMClone(dm, &dmAux);CHKERRQ(ierr);
544*c4762a1bSJed Brown   ierr = PetscObjectCompose((PetscObject) dm, "dmAux", (PetscObject) dmAux);CHKERRQ(ierr);
545*c4762a1bSJed Brown   ierr = DMSetCoordinateDM(dmAux, coordDM);CHKERRQ(ierr);
546*c4762a1bSJed Brown   for (f = 0; f < NfAux; ++f) {ierr = DMSetField(dmAux, f, NULL, (PetscObject) feAux[f]);CHKERRQ(ierr);}
547*c4762a1bSJed Brown   ierr = DMCreateDS(dmAux);CHKERRQ(ierr);
548*c4762a1bSJed Brown   ierr = SetupEquilibriumFields(dm, dmAux, user);CHKERRQ(ierr);
549*c4762a1bSJed Brown   ierr = DMDestroy(&dmAux);CHKERRQ(ierr);
550*c4762a1bSJed Brown   PetscFunctionReturn(0);
551*c4762a1bSJed Brown }
552*c4762a1bSJed Brown 
553*c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
554*c4762a1bSJed Brown {
555*c4762a1bSJed Brown   DM              cdm = dm;
556*c4762a1bSJed Brown   const PetscInt  dim = ctx->dim;
557*c4762a1bSJed Brown   PetscFE         fe[5], feAux[3];
558*c4762a1bSJed Brown   PetscInt        Nf = 5, NfAux = 3, f;
559*c4762a1bSJed Brown   PetscBool       cell_simplex = ctx->cell_simplex;
560*c4762a1bSJed Brown   MPI_Comm        comm;
561*c4762a1bSJed Brown   PetscErrorCode  ierr;
562*c4762a1bSJed Brown 
563*c4762a1bSJed Brown   PetscFunctionBeginUser;
564*c4762a1bSJed Brown   /* Create finite element */
565*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
566*c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, cell_simplex, NULL, -1, &fe[0]);CHKERRQ(ierr);
567*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[0], "density");CHKERRQ(ierr);
568*c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, cell_simplex, NULL, -1, &fe[1]);CHKERRQ(ierr);
569*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[1], "vorticity");CHKERRQ(ierr);
570*c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
571*c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, cell_simplex, NULL, -1, &fe[2]);CHKERRQ(ierr);
572*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[2], "flux");CHKERRQ(ierr);
573*c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr);
574*c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, cell_simplex, NULL, -1, &fe[3]);CHKERRQ(ierr);
575*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[3], "potential");CHKERRQ(ierr);
576*c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], fe[3]);CHKERRQ(ierr);
577*c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, cell_simplex, NULL, -1, &fe[4]);CHKERRQ(ierr);
578*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[4], "current");CHKERRQ(ierr);
579*c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], fe[4]);CHKERRQ(ierr);
580*c4762a1bSJed Brown 
581*c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, cell_simplex, NULL, -1, &feAux[0]);CHKERRQ(ierr);
582*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) feAux[0], "n_0");CHKERRQ(ierr);
583*c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], feAux[0]);CHKERRQ(ierr);
584*c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, cell_simplex, NULL, -1, &feAux[1]);CHKERRQ(ierr);
585*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) feAux[1], "vorticity_0");CHKERRQ(ierr);
586*c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], feAux[1]);CHKERRQ(ierr);
587*c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, cell_simplex, NULL, -1, &feAux[2]);CHKERRQ(ierr);
588*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) feAux[2], "flux_0");CHKERRQ(ierr);
589*c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], feAux[2]);CHKERRQ(ierr);
590*c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
591*c4762a1bSJed Brown   for (f = 0; f < Nf; ++f) {ierr = DMSetField(dm, f, NULL, (PetscObject) fe[f]);CHKERRQ(ierr);}
592*c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
593*c4762a1bSJed Brown   ierr = SetupProblem(dm, ctx);CHKERRQ(ierr);
594*c4762a1bSJed Brown   while (cdm) {
595*c4762a1bSJed Brown     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
596*c4762a1bSJed Brown     ierr = SetupAuxDM(dm, NfAux, feAux, ctx);CHKERRQ(ierr);
597*c4762a1bSJed Brown     {
598*c4762a1bSJed Brown       PetscBool hasLabel;
599*c4762a1bSJed Brown 
600*c4762a1bSJed Brown       ierr = DMHasLabel(cdm, "marker", &hasLabel);CHKERRQ(ierr);
601*c4762a1bSJed Brown       if (!hasLabel) {ierr = CreateBCLabel(cdm, "marker");CHKERRQ(ierr);}
602*c4762a1bSJed Brown     }
603*c4762a1bSJed Brown     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
604*c4762a1bSJed Brown   }
605*c4762a1bSJed Brown   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&fe[f]);CHKERRQ(ierr);}
606*c4762a1bSJed Brown   for (f = 0; f < NfAux; ++f) {ierr = PetscFEDestroy(&feAux[f]);CHKERRQ(ierr);}
607*c4762a1bSJed Brown   PetscFunctionReturn(0);
608*c4762a1bSJed Brown }
609*c4762a1bSJed Brown 
610*c4762a1bSJed Brown int main(int argc, char **argv)
611*c4762a1bSJed Brown {
612*c4762a1bSJed Brown   DM             dm;
613*c4762a1bSJed Brown   TS             ts;
614*c4762a1bSJed Brown   Vec            u, r;
615*c4762a1bSJed Brown   AppCtx         ctx;
616*c4762a1bSJed Brown   PetscReal      t       = 0.0;
617*c4762a1bSJed Brown   PetscReal      L2error = 0.0;
618*c4762a1bSJed Brown   PetscErrorCode ierr;
619*c4762a1bSJed Brown   AppCtx        *ctxarr[5];
620*c4762a1bSJed Brown 
621*c4762a1bSJed Brown   ctxarr[0] = ctxarr[1] = ctxarr[2] = ctxarr[3] = ctxarr[4] = &ctx; /* each variable could have a different context */
622*c4762a1bSJed Brown   s_ctx = &ctx;
623*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
624*c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr);
625*c4762a1bSJed Brown   /* create mesh and problem */
626*c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &ctx, &dm);CHKERRQ(ierr);
627*c4762a1bSJed Brown   ierr = DMSetApplicationContext(dm, &ctx);CHKERRQ(ierr);
628*c4762a1bSJed Brown   ierr = PetscMalloc1(5, &ctx.initialFuncs);CHKERRQ(ierr);
629*c4762a1bSJed Brown   ierr = SetupDiscretization(dm, &ctx);CHKERRQ(ierr);
630*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
631*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) u, "u");CHKERRQ(ierr);
632*c4762a1bSJed Brown   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
633*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) r, "r");CHKERRQ(ierr);
634*c4762a1bSJed Brown   /* create TS */
635*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
636*c4762a1bSJed Brown   ierr = TSSetDM(ts, dm);CHKERRQ(ierr);
637*c4762a1bSJed Brown   ierr = TSSetApplicationContext(ts, &ctx);CHKERRQ(ierr);
638*c4762a1bSJed Brown   ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx);CHKERRQ(ierr);
639*c4762a1bSJed Brown   if (ctx.implicit) {
640*c4762a1bSJed Brown     ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx);CHKERRQ(ierr);
641*c4762a1bSJed Brown     ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx);CHKERRQ(ierr);
642*c4762a1bSJed Brown   } else {
643*c4762a1bSJed Brown     ierr = DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &ctx);CHKERRQ(ierr);
644*c4762a1bSJed Brown   }
645*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
646*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
647*c4762a1bSJed Brown   ierr = TSSetPostStep(ts, PostStep);CHKERRQ(ierr);
648*c4762a1bSJed Brown   /* make solution & solve */
649*c4762a1bSJed Brown   ierr = DMProjectFunction(dm, t, ctx.initialFuncs, (void **)ctxarr, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
650*c4762a1bSJed Brown   ierr = TSSetSolution(ts,u);CHKERRQ(ierr);
651*c4762a1bSJed Brown   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
652*c4762a1bSJed Brown   ierr = PostStep(ts);CHKERRQ(ierr); /* print the initial state */
653*c4762a1bSJed Brown   ierr = TSSolve(ts, u);CHKERRQ(ierr);
654*c4762a1bSJed Brown   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
655*c4762a1bSJed Brown   ierr = DMComputeL2Diff(dm, t, ctx.initialFuncs, (void **)ctxarr, u, &L2error);CHKERRQ(ierr);
656*c4762a1bSJed Brown   if (L2error < 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
657*c4762a1bSJed Brown   else                   {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", L2error);CHKERRQ(ierr);}
658*c4762a1bSJed Brown #if 0
659*c4762a1bSJed Brown   {
660*c4762a1bSJed Brown     PetscReal res = 0.0;
661*c4762a1bSJed Brown     /* Check discretization error */
662*c4762a1bSJed Brown     ierr = VecViewFromOptions(u, NULL, "-initial_guess_view");CHKERRQ(ierr);
663*c4762a1bSJed Brown     ierr = DMComputeL2Diff(dm, 0.0, ctx.exactFuncs, NULL, u, &error);CHKERRQ(ierr);
664*c4762a1bSJed Brown     if (error < 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
665*c4762a1bSJed Brown     else                 {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", error);CHKERRQ(ierr);}
666*c4762a1bSJed Brown     /* Check residual */
667*c4762a1bSJed Brown     ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
668*c4762a1bSJed Brown     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
669*c4762a1bSJed Brown     ierr = VecViewFromOptions(r, NULL, "-initial_residual_view");CHKERRQ(ierr);
670*c4762a1bSJed Brown     ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
671*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", res);CHKERRQ(ierr);
672*c4762a1bSJed Brown     /* Check Jacobian */
673*c4762a1bSJed Brown     {
674*c4762a1bSJed Brown       Mat A;
675*c4762a1bSJed Brown       Vec b;
676*c4762a1bSJed Brown 
677*c4762a1bSJed Brown       ierr = SNESGetJacobian(snes, &A, NULL, NULL, NULL);CHKERRQ(ierr);
678*c4762a1bSJed Brown       ierr = SNESComputeJacobian(snes, u, A, A);CHKERRQ(ierr);
679*c4762a1bSJed Brown       ierr = VecDuplicate(u, &b);CHKERRQ(ierr);
680*c4762a1bSJed Brown       ierr = VecSet(r, 0.0);CHKERRQ(ierr);
681*c4762a1bSJed Brown       ierr = SNESComputeFunction(snes, r, b);CHKERRQ(ierr);
682*c4762a1bSJed Brown       ierr = MatMult(A, u, r);CHKERRQ(ierr);
683*c4762a1bSJed Brown       ierr = VecAXPY(r, 1.0, b);CHKERRQ(ierr);
684*c4762a1bSJed Brown       ierr = VecDestroy(&b);CHKERRQ(ierr);
685*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");CHKERRQ(ierr);
686*c4762a1bSJed Brown       ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
687*c4762a1bSJed Brown       ierr = VecViewFromOptions(r, NULL, "-linear_residual_view");CHKERRQ(ierr);
688*c4762a1bSJed Brown       ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
689*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", res);CHKERRQ(ierr);
690*c4762a1bSJed Brown     }
691*c4762a1bSJed Brown   }
692*c4762a1bSJed Brown #endif
693*c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
694*c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
695*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
696*c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
697*c4762a1bSJed Brown   ierr = PetscFree(ctx.initialFuncs);CHKERRQ(ierr);
698*c4762a1bSJed Brown   ierr = PetscFinalize();
699*c4762a1bSJed Brown   return ierr;
700*c4762a1bSJed Brown }
701*c4762a1bSJed Brown 
702*c4762a1bSJed Brown /*TEST
703*c4762a1bSJed Brown 
704*c4762a1bSJed Brown   test:
705*c4762a1bSJed Brown     suffix: 0
706*c4762a1bSJed Brown     args: -debug 1 -dim 2 -dm_refine 1 -x_periodicity PERIODIC -ts_max_steps 1 -ts_max_time 10. -ts_dt 1.0
707*c4762a1bSJed Brown   test:
708*c4762a1bSJed Brown     suffix: 1
709*c4762a1bSJed Brown     args: -debug 1 -dim 3 -dm_refine 1 -z_periodicity PERIODIC -ts_max_steps 1 -ts_max_time 10. -ts_dt 1.0 -domain_lo -2,-1,-1 -domain_hi 2,1,1
710*c4762a1bSJed Brown 
711*c4762a1bSJed Brown TEST*/
712