1c4762a1bSJed Brown static char help[] = "Evolution of magnetic islands.\n\ 2c4762a1bSJed 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"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*F 5c4762a1bSJed 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$. 6c4762a1bSJed Brown \begin{equation} 7c4762a1bSJed Brown \begin{aligned} 8c4762a1bSJed 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 \\ 9c4762a1bSJed Brown \partial_t \tilde\Omega &= \left\{ \tilde\Omega, \tilde\phi \right\} + \beta \left\{ j_z, \tilde\psi \right\} + \mu \nabla^2_\perp \tilde\Omega \\ 10c4762a1bSJed 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 \\ 11c4762a1bSJed Brown \nabla^2_\perp\tilde\phi &= \tilde\Omega \\ 12c4762a1bSJed Brown j_z &= -\nabla^2_\perp \left(\tilde\psi + \psi_0 \right)\\ 13c4762a1bSJed Brown \end{aligned} 14c4762a1bSJed Brown \end{equation} 15c4762a1bSJed Brown F*/ 16c4762a1bSJed Brown 17c4762a1bSJed Brown #include <petscdmplex.h> 18c4762a1bSJed Brown #include <petscts.h> 19c4762a1bSJed Brown #include <petscds.h> 20c4762a1bSJed Brown 21c4762a1bSJed Brown typedef struct { 22c4762a1bSJed Brown PetscInt debug; /* The debugging level */ 23c4762a1bSJed Brown PetscBool plotRef; /* Plot the reference fields */ 2430602db0SMatthew G. Knepley PetscReal lower[3], upper[3]; 25c4762a1bSJed Brown /* Problem definition */ 26c4762a1bSJed Brown PetscErrorCode (**initialFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 27c4762a1bSJed Brown PetscReal mu, eta, beta; 28c4762a1bSJed Brown PetscReal a, b, Jo, Jop, m, ke, kx, ky, DeltaPrime, eps; 29c4762a1bSJed Brown /* solver */ 30c4762a1bSJed Brown PetscBool implicit; 31c4762a1bSJed Brown } AppCtx; 32c4762a1bSJed Brown 33c4762a1bSJed Brown static AppCtx *s_ctx; 34c4762a1bSJed Brown 35d71ae5a4SJacob Faibussowitsch static PetscScalar poissonBracket(PetscInt dim, const PetscScalar df[], const PetscScalar dg[]) 36d71ae5a4SJacob Faibussowitsch { 37c4762a1bSJed Brown PetscScalar ret = df[0] * dg[1] - df[1] * dg[0]; 38c4762a1bSJed Brown return ret; 39c4762a1bSJed Brown } 40c4762a1bSJed Brown 419371c9d4SSatish Balay enum field_idx { 429371c9d4SSatish Balay DENSITY, 439371c9d4SSatish Balay OMEGA, 449371c9d4SSatish Balay PSI, 459371c9d4SSatish Balay PHI, 469371c9d4SSatish Balay JZ 479371c9d4SSatish Balay }; 48c4762a1bSJed Brown 49d71ae5a4SJacob Faibussowitsch static void f0_n(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 50d71ae5a4SJacob Faibussowitsch { 51c4762a1bSJed Brown const PetscScalar *pnDer = &u_x[uOff_x[DENSITY]]; 52c4762a1bSJed Brown const PetscScalar *ppsiDer = &u_x[uOff_x[PSI]]; 53c4762a1bSJed Brown const PetscScalar *pphiDer = &u_x[uOff_x[PHI]]; 54c4762a1bSJed Brown const PetscScalar *jzDer = &u_x[uOff_x[JZ]]; 55c4762a1bSJed Brown const PetscScalar *logRefDenDer = &a_x[aOff_x[DENSITY]]; 56c4762a1bSJed Brown f0[0] += -poissonBracket(dim, pnDer, pphiDer) - s_ctx->beta * poissonBracket(dim, jzDer, ppsiDer) - poissonBracket(dim, logRefDenDer, pphiDer); 57c4762a1bSJed Brown if (u_t) f0[0] += u_t[DENSITY]; 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60d71ae5a4SJacob Faibussowitsch static void f1_n(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 61d71ae5a4SJacob Faibussowitsch { 62c4762a1bSJed Brown const PetscScalar *pnDer = &u_x[uOff_x[DENSITY]]; 63c4762a1bSJed Brown PetscInt d; 64c4762a1bSJed Brown 65c4762a1bSJed Brown for (d = 0; d < dim - 1; ++d) f1[d] = -s_ctx->mu * pnDer[d]; 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68d71ae5a4SJacob Faibussowitsch static void f0_Omega(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 69d71ae5a4SJacob Faibussowitsch { 70c4762a1bSJed Brown const PetscScalar *pOmegaDer = &u_x[uOff_x[OMEGA]]; 71c4762a1bSJed Brown const PetscScalar *ppsiDer = &u_x[uOff_x[PSI]]; 72c4762a1bSJed Brown const PetscScalar *pphiDer = &u_x[uOff_x[PHI]]; 73c4762a1bSJed Brown const PetscScalar *jzDer = &u_x[uOff_x[JZ]]; 74c4762a1bSJed Brown 75c4762a1bSJed Brown f0[0] += -poissonBracket(dim, pOmegaDer, pphiDer) - s_ctx->beta * poissonBracket(dim, jzDer, ppsiDer); 76c4762a1bSJed Brown if (u_t) f0[0] += u_t[OMEGA]; 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79d71ae5a4SJacob Faibussowitsch static void f1_Omega(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 80d71ae5a4SJacob Faibussowitsch { 81c4762a1bSJed Brown const PetscScalar *pOmegaDer = &u_x[uOff_x[OMEGA]]; 82c4762a1bSJed Brown PetscInt d; 83c4762a1bSJed Brown 84c4762a1bSJed Brown for (d = 0; d < dim - 1; ++d) f1[d] = -s_ctx->mu * pOmegaDer[d]; 85c4762a1bSJed Brown } 86c4762a1bSJed Brown 87d71ae5a4SJacob Faibussowitsch static void f0_psi(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 88d71ae5a4SJacob Faibussowitsch { 89c4762a1bSJed Brown const PetscScalar *pnDer = &u_x[uOff_x[DENSITY]]; 90c4762a1bSJed Brown const PetscScalar *ppsiDer = &u_x[uOff_x[PSI]]; 91c4762a1bSJed Brown const PetscScalar *pphiDer = &u_x[uOff_x[PHI]]; 92c4762a1bSJed Brown const PetscScalar *refPsiDer = &a_x[aOff_x[PSI]]; 93c4762a1bSJed Brown const PetscScalar *logRefDenDer = &a_x[aOff_x[DENSITY]]; 94c4762a1bSJed Brown PetscScalar psiDer[3]; 95c4762a1bSJed Brown PetscScalar phi_n_Der[3]; 96c4762a1bSJed Brown PetscInt d; 979371c9d4SSatish Balay if (dim < 2) { 989371c9d4SSatish Balay MPI_Abort(MPI_COMM_WORLD, 1); 999371c9d4SSatish Balay return; 1009371c9d4SSatish Balay } /* this is needed so that the clang static analyzer does not generate a warning about variables used by not set */ 101c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 102c4762a1bSJed Brown psiDer[d] = refPsiDer[d] + ppsiDer[d]; 103c4762a1bSJed Brown phi_n_Der[d] = pphiDer[d] - pnDer[d]; 104c4762a1bSJed Brown } 105c4762a1bSJed Brown f0[0] = -poissonBracket(dim, psiDer, phi_n_Der) + poissonBracket(dim, logRefDenDer, ppsiDer); 106c4762a1bSJed Brown if (u_t) f0[0] += u_t[PSI]; 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109d71ae5a4SJacob Faibussowitsch static void f1_psi(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 110d71ae5a4SJacob Faibussowitsch { 111c4762a1bSJed Brown const PetscScalar *ppsi = &u_x[uOff_x[PSI]]; 112c4762a1bSJed Brown PetscInt d; 113c4762a1bSJed Brown 114c4762a1bSJed Brown for (d = 0; d < dim - 1; ++d) f1[d] = -(s_ctx->eta / s_ctx->beta) * ppsi[d]; 115c4762a1bSJed Brown } 116c4762a1bSJed Brown 117d71ae5a4SJacob Faibussowitsch static void f0_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 118d71ae5a4SJacob Faibussowitsch { 119c4762a1bSJed Brown f0[0] = -u[uOff[OMEGA]]; 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122d71ae5a4SJacob Faibussowitsch static void f1_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 123d71ae5a4SJacob Faibussowitsch { 124c4762a1bSJed Brown const PetscScalar *pphi = &u_x[uOff_x[PHI]]; 125c4762a1bSJed Brown PetscInt d; 126c4762a1bSJed Brown 127c4762a1bSJed Brown for (d = 0; d < dim - 1; ++d) f1[d] = pphi[d]; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130d71ae5a4SJacob Faibussowitsch static void f0_jz(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 131d71ae5a4SJacob Faibussowitsch { 132c4762a1bSJed Brown f0[0] = u[uOff[JZ]]; 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135d71ae5a4SJacob Faibussowitsch static void f1_jz(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 136d71ae5a4SJacob Faibussowitsch { 137c4762a1bSJed Brown const PetscScalar *ppsi = &u_x[uOff_x[PSI]]; 138c4762a1bSJed Brown const PetscScalar *refPsiDer = &a_x[aOff_x[PSI]]; /* aOff_x[PSI] == 2*PSI */ 139c4762a1bSJed Brown PetscInt d; 140c4762a1bSJed Brown 141c4762a1bSJed Brown for (d = 0; d < dim - 1; ++d) f1[d] = ppsi[d] + refPsiDer[d]; 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 145d71ae5a4SJacob Faibussowitsch { 146c4762a1bSJed Brown PetscFunctionBeginUser; 147c4762a1bSJed Brown options->debug = 1; 148c4762a1bSJed Brown options->plotRef = PETSC_FALSE; 149c4762a1bSJed Brown options->implicit = PETSC_FALSE; 150c4762a1bSJed Brown options->mu = 0; 151c4762a1bSJed Brown options->eta = 0; 152c4762a1bSJed Brown options->beta = 1; 153c4762a1bSJed Brown options->a = 1; 154c4762a1bSJed Brown options->b = PETSC_PI; 155c4762a1bSJed Brown options->Jop = 0; 156c4762a1bSJed Brown options->m = 1; 157c4762a1bSJed Brown options->eps = 1.e-6; 158c4762a1bSJed Brown 159d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX"); 1609566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-debug", "The debugging level", "ex48.c", options->debug, &options->debug, NULL)); 1619566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-plot_ref", "Plot the reference fields", "ex48.c", options->plotRef, &options->plotRef, NULL)); 1629566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mu", "mu", "ex48.c", options->mu, &options->mu, NULL)); 1639566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-eta", "eta", "ex48.c", options->eta, &options->eta, NULL)); 1649566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-beta", "beta", "ex48.c", options->beta, &options->beta, NULL)); 1659566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-Jop", "Jop", "ex48.c", options->Jop, &options->Jop, NULL)); 1669566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-m", "m", "ex48.c", options->m, &options->m, NULL)); 1679566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-eps", "eps", "ex48.c", options->eps, &options->eps, NULL)); 1689566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-implicit", "Use implicit time integrator", "ex48.c", options->implicit, &options->implicit, NULL)); 169d0609cedSBarry Smith PetscOptionsEnd(); 170c4762a1bSJed Brown options->ke = PetscSqrtScalar(options->Jop); 171c4762a1bSJed Brown if (options->Jop == 0.0) { 172c4762a1bSJed Brown options->Jo = 1.0 / PetscPowScalar(options->a, 2); 173c4762a1bSJed Brown } else { 174c4762a1bSJed Brown options->Jo = options->Jop * PetscCosReal(options->ke * options->a) / (1.0 - PetscCosReal(options->ke * options->a)); 175c4762a1bSJed Brown } 176c4762a1bSJed Brown options->ky = PETSC_PI * options->m / options->b; 177c4762a1bSJed Brown if (PetscPowReal(options->ky, 2) < options->Jop) { 178c4762a1bSJed Brown options->kx = PetscSqrtScalar(options->Jop - PetscPowScalar(options->ky, 2)); 179c4762a1bSJed Brown options->DeltaPrime = -2.0 * options->kx * options->a * PetscCosReal(options->kx * options->a) / PetscSinReal(options->kx * options->a); 180c4762a1bSJed Brown } else if (PetscPowReal(options->ky, 2) > options->Jop) { 181c4762a1bSJed Brown options->kx = PetscSqrtScalar(PetscPowScalar(options->ky, 2) - options->Jop); 182c4762a1bSJed Brown options->DeltaPrime = -2.0 * options->kx * options->a * PetscCoshReal(options->kx * options->a) / PetscSinhReal(options->kx * options->a); 183c4762a1bSJed Brown } else { /*they're equal (or there's a NaN), lim(x*cot(x))_x->0=1*/ 184c4762a1bSJed Brown options->kx = 0; 185c4762a1bSJed Brown options->DeltaPrime = -2.0; 186c4762a1bSJed Brown } 18763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "DeltaPrime=%g\n", (double)options->DeltaPrime)); 188c4762a1bSJed Brown 189*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 190c4762a1bSJed Brown } 191c4762a1bSJed Brown 192d71ae5a4SJacob Faibussowitsch static void f_n(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 193d71ae5a4SJacob Faibussowitsch { 194c4762a1bSJed Brown const PetscScalar *pn = &u[uOff[DENSITY]]; 195c4762a1bSJed Brown *f0 = *pn; 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 198d71ae5a4SJacob Faibussowitsch static PetscErrorCode PostStep(TS ts) 199d71ae5a4SJacob Faibussowitsch { 200c4762a1bSJed Brown DM dm; 201c4762a1bSJed Brown AppCtx *ctx; 202c4762a1bSJed Brown PetscInt stepi, num; 203c4762a1bSJed Brown Vec X; 2045f80ce2aSJacob Faibussowitsch 2057510d9b0SBarry Smith PetscFunctionBeginUser; 2069566063dSJacob Faibussowitsch PetscCall(TSGetApplicationContext(ts, &ctx)); 207*3ba16761SJacob Faibussowitsch if (ctx->debug < 1) PetscFunctionReturn(PETSC_SUCCESS); 2089566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &X)); 2099566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dm)); 2109566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &stepi)); 2119566063dSJacob Faibussowitsch PetscCall(DMGetOutputSequenceNumber(dm, &num, NULL)); 2129566063dSJacob Faibussowitsch if (num < 0) PetscCall(DMSetOutputSequenceNumber(dm, 0, 0.0)); 2139566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)X, "u")); 2149566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X, NULL, "-vec_view")); 215c4762a1bSJed Brown /* print integrals */ 216c4762a1bSJed Brown { 217c4762a1bSJed Brown PetscDS prob; 218c4762a1bSJed Brown DM plex; 219c4762a1bSJed Brown PetscScalar den, tt[5]; 2209566063dSJacob Faibussowitsch PetscCall(DMConvert(dm, DMPLEX, &plex)); 2219566063dSJacob Faibussowitsch PetscCall(DMGetDS(plex, &prob)); 2229566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f_n)); 2239566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plex, X, tt, ctx)); 224c4762a1bSJed Brown den = tt[0]; 2259566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 22663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "%" PetscInt_FMT ") total perturbed mass = %g\n", stepi, (double)PetscRealPart(den))); 227c4762a1bSJed Brown } 228*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 229c4762a1bSJed Brown } 230c4762a1bSJed Brown 231d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 232d71ae5a4SJacob Faibussowitsch { 233c4762a1bSJed Brown PetscFunctionBeginUser; 2349566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 2359566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 2369566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2379566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 23830602db0SMatthew G. Knepley 2399566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(*dm, ctx->lower, ctx->upper)); 24030602db0SMatthew G. Knepley ctx->a = (ctx->upper[0] - ctx->lower[0]) / 2.0; 24130602db0SMatthew G. Knepley ctx->b = (ctx->upper[1] - ctx->lower[1]) / 2.0; 242*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown 245d71ae5a4SJacob Faibussowitsch static PetscErrorCode log_n_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 246d71ae5a4SJacob Faibussowitsch { 247c4762a1bSJed Brown AppCtx *lctx = (AppCtx *)ctx; 24830602db0SMatthew G. Knepley u[0] = 2. * lctx->a + x[0]; 249*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252d71ae5a4SJacob Faibussowitsch static PetscErrorCode Omega_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 253d71ae5a4SJacob Faibussowitsch { 254c4762a1bSJed Brown u[0] = 0.0; 255*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258d71ae5a4SJacob Faibussowitsch static PetscErrorCode psi_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 259d71ae5a4SJacob Faibussowitsch { 260c4762a1bSJed Brown AppCtx *lctx = (AppCtx *)ctx; 261c4762a1bSJed Brown /* This sets up a symmetrix By flux aroound the mid point in x, which represents a current density flux along z. The stability 262c4762a1bSJed Brown is analytically known and reported in ProcessOptions. */ 263c4762a1bSJed Brown if (lctx->ke != 0.0) { 264c4762a1bSJed Brown u[0] = (PetscCosReal(lctx->ke * (x[0] - lctx->a)) - PetscCosReal(lctx->ke * lctx->a)) / (1.0 - PetscCosReal(lctx->ke * lctx->a)); 265c4762a1bSJed Brown } else { 266c4762a1bSJed Brown u[0] = 1.0 - PetscPowScalar((x[0] - lctx->a) / lctx->a, 2); 267c4762a1bSJed Brown } 268*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 269c4762a1bSJed Brown } 270c4762a1bSJed Brown 271d71ae5a4SJacob Faibussowitsch static PetscErrorCode initialSolution_n(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 272d71ae5a4SJacob Faibussowitsch { 273c4762a1bSJed Brown u[0] = 0.0; 274*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 275c4762a1bSJed Brown } 276c4762a1bSJed Brown 277d71ae5a4SJacob Faibussowitsch static PetscErrorCode initialSolution_Omega(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 278d71ae5a4SJacob Faibussowitsch { 279c4762a1bSJed Brown u[0] = 0.0; 280*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 281c4762a1bSJed Brown } 282c4762a1bSJed Brown 283d71ae5a4SJacob Faibussowitsch static PetscErrorCode initialSolution_psi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *a_ctx) 284d71ae5a4SJacob Faibussowitsch { 285c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)a_ctx; 286c4762a1bSJed Brown PetscScalar r = ctx->eps * (PetscScalar)(rand()) / (PetscScalar)(RAND_MAX); 28730602db0SMatthew G. Knepley if (x[0] == ctx->lower[0] || x[0] == ctx->upper[0]) r = 0; 288c4762a1bSJed Brown u[0] = r; 289c4762a1bSJed Brown /* PetscPrintf(PETSC_COMM_WORLD, "rand psi %lf\n",u[0]); */ 290*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 293d71ae5a4SJacob Faibussowitsch static PetscErrorCode initialSolution_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 294d71ae5a4SJacob Faibussowitsch { 295c4762a1bSJed Brown u[0] = 0.0; 296*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 297c4762a1bSJed Brown } 298c4762a1bSJed Brown 299d71ae5a4SJacob Faibussowitsch static PetscErrorCode initialSolution_jz(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 300d71ae5a4SJacob Faibussowitsch { 301c4762a1bSJed Brown u[0] = 0.0; 302*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 303c4762a1bSJed Brown } 304c4762a1bSJed Brown 305d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) 306d71ae5a4SJacob Faibussowitsch { 30745480ffeSMatthew G. Knepley PetscDS ds; 30845480ffeSMatthew G. Knepley DMLabel label; 309c4762a1bSJed Brown const PetscInt id = 1; 310c4762a1bSJed Brown 311c4762a1bSJed Brown PetscFunctionBeginUser; 3129566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 3139566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3149566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_n, f1_n)); 3159566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_Omega, f1_Omega)); 3169566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 2, f0_psi, f1_psi)); 3179566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 3, f0_phi, f1_phi)); 3189566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 4, f0_jz, f1_jz)); 319c4762a1bSJed Brown ctx->initialFuncs[0] = initialSolution_n; 320c4762a1bSJed Brown ctx->initialFuncs[1] = initialSolution_Omega; 321c4762a1bSJed Brown ctx->initialFuncs[2] = initialSolution_psi; 322c4762a1bSJed Brown ctx->initialFuncs[3] = initialSolution_phi; 323c4762a1bSJed Brown ctx->initialFuncs[4] = initialSolution_jz; 3245f80ce2aSJacob Faibussowitsch for (PetscInt f = 0; f < 5; ++f) { 3259566063dSJacob Faibussowitsch PetscCall(PetscDSSetImplicit(ds, f, ctx->implicit)); 3269566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, f, 0, NULL, (void (*)(void))ctx->initialFuncs[f], NULL, ctx, NULL)); 327c4762a1bSJed Brown } 3289566063dSJacob Faibussowitsch PetscCall(PetscDSSetContext(ds, 0, ctx)); 329*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 330c4762a1bSJed Brown } 331c4762a1bSJed Brown 332d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupEquilibriumFields(DM dm, DM dmAux, AppCtx *ctx) 333d71ae5a4SJacob Faibussowitsch { 334c4762a1bSJed Brown PetscErrorCode (*eqFuncs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *) = {log_n_0, Omega_0, psi_0}; 335c4762a1bSJed Brown Vec eq; 336c4762a1bSJed Brown AppCtx *ctxarr[3]; 337c4762a1bSJed Brown 338c4762a1bSJed Brown ctxarr[0] = ctxarr[1] = ctxarr[2] = ctx; /* each variable could have a different context */ 3397510d9b0SBarry Smith PetscFunctionBeginUser; 3409566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dmAux, &eq)); 3419566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dmAux, 0.0, eqFuncs, (void **)ctxarr, INSERT_ALL_VALUES, eq)); 3429566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, eq)); 343c4762a1bSJed Brown if (ctx->plotRef) { /* plot reference functions */ 344c4762a1bSJed Brown PetscViewer viewer = NULL; 345c4762a1bSJed Brown PetscBool isHDF5, isVTK; 346c4762a1bSJed Brown char buf[256]; 347c4762a1bSJed Brown Vec global; 34830602db0SMatthew G. Knepley PetscInt dim; 34930602db0SMatthew G. Knepley 3509566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 3519566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmAux, &global)); 3529566063dSJacob Faibussowitsch PetscCall(VecSet(global, .0)); /* BCs! */ 3539566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dmAux, eq, INSERT_VALUES, global)); 3549566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dmAux, eq, INSERT_VALUES, global)); 3559566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)dmAux), &viewer)); 356c4762a1bSJed Brown #ifdef PETSC_HAVE_HDF5 3579566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERHDF5)); 358c4762a1bSJed Brown #else 3599566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERVTK)); 360c4762a1bSJed Brown #endif 3619566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions(viewer)); 3629566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &isHDF5)); 3639566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isVTK)); 364c4762a1bSJed Brown if (isHDF5) { 36563a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, 256, "uEquilibrium-%" PetscInt_FMT "D.h5", dim)); 366c4762a1bSJed Brown } else if (isVTK) { 36763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, 256, "uEquilibrium-%" PetscInt_FMT "D.vtu", dim)); 3689566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_VTK_VTU)); 369c4762a1bSJed Brown } 3709566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_WRITE)); 3719566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer, buf)); 3729566063dSJacob Faibussowitsch if (isHDF5) PetscCall(DMView(dmAux, viewer)); 373c4762a1bSJed Brown /* view equilibrium fields, this will overwrite fine grids with coarse grids! */ 3749566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)global, "u0")); 3759566063dSJacob Faibussowitsch PetscCall(VecView(global, viewer)); 3769566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 3779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 378c4762a1bSJed Brown } 3799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&eq)); 380*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 381c4762a1bSJed Brown } 382c4762a1bSJed Brown 383d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupAuxDM(DM dm, PetscInt NfAux, PetscFE feAux[], AppCtx *user) 384d71ae5a4SJacob Faibussowitsch { 385c4762a1bSJed Brown DM dmAux, coordDM; 386c4762a1bSJed Brown PetscInt f; 387c4762a1bSJed Brown 3887510d9b0SBarry Smith PetscFunctionBeginUser; 389c4762a1bSJed Brown /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */ 3909566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &coordDM)); 391*3ba16761SJacob Faibussowitsch if (!feAux) PetscFunctionReturn(PETSC_SUCCESS); 3929566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmAux)); 3939566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDM(dmAux, coordDM)); 3949566063dSJacob Faibussowitsch for (f = 0; f < NfAux; ++f) PetscCall(DMSetField(dmAux, f, NULL, (PetscObject)feAux[f])); 3959566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmAux)); 3969566063dSJacob Faibussowitsch PetscCall(SetupEquilibriumFields(dm, dmAux, user)); 3979566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmAux)); 398*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 399c4762a1bSJed Brown } 400c4762a1bSJed Brown 401d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx) 402d71ae5a4SJacob Faibussowitsch { 403c4762a1bSJed Brown DM cdm = dm; 404c4762a1bSJed Brown PetscFE fe[5], feAux[3]; 40530602db0SMatthew G. Knepley PetscInt dim, Nf = 5, NfAux = 3, f; 40630602db0SMatthew G. Knepley PetscBool simplex; 407c4762a1bSJed Brown MPI_Comm comm; 408c4762a1bSJed Brown 409c4762a1bSJed Brown PetscFunctionBeginUser; 410c4762a1bSJed Brown /* Create finite element */ 4119566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 4129566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4139566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 4149566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[0])); 4159566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "density")); 4169566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[1])); 4179566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "vorticity")); 4189566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 4199566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[2])); 4209566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[2], "flux")); 4219566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[2])); 4229566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[3])); 4239566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[3], "potential")); 4249566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[3])); 4259566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[4])); 4269566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[4], "current")); 4279566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[4])); 428c4762a1bSJed Brown 4299566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &feAux[0])); 4309566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)feAux[0], "n_0")); 4319566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], feAux[0])); 4329566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &feAux[1])); 4339566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)feAux[1], "vorticity_0")); 4349566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], feAux[1])); 4359566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &feAux[2])); 4369566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)feAux[2], "flux_0")); 4379566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], feAux[2])); 438c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 4399566063dSJacob Faibussowitsch for (f = 0; f < Nf; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe[f])); 4409566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 4419566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, ctx)); 442c4762a1bSJed Brown while (cdm) { 4439566063dSJacob Faibussowitsch PetscCall(SetupAuxDM(dm, NfAux, feAux, ctx)); 4449566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 4459566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 446c4762a1bSJed Brown } 4479566063dSJacob Faibussowitsch for (f = 0; f < Nf; ++f) PetscCall(PetscFEDestroy(&fe[f])); 4489566063dSJacob Faibussowitsch for (f = 0; f < NfAux; ++f) PetscCall(PetscFEDestroy(&feAux[f])); 449*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 450c4762a1bSJed Brown } 451c4762a1bSJed Brown 452d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 453d71ae5a4SJacob Faibussowitsch { 454c4762a1bSJed Brown DM dm; 455c4762a1bSJed Brown TS ts; 456c4762a1bSJed Brown Vec u, r; 457c4762a1bSJed Brown AppCtx ctx; 458c4762a1bSJed Brown PetscReal t = 0.0; 459c4762a1bSJed Brown PetscReal L2error = 0.0; 460c4762a1bSJed Brown AppCtx *ctxarr[5]; 461c4762a1bSJed Brown 462c4762a1bSJed Brown ctxarr[0] = ctxarr[1] = ctxarr[2] = ctxarr[3] = ctxarr[4] = &ctx; /* each variable could have a different context */ 463c4762a1bSJed Brown s_ctx = &ctx; 464327415f7SBarry Smith PetscFunctionBeginUser; 4659566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 4669566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 467c4762a1bSJed Brown /* create mesh and problem */ 4689566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm)); 4699566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &ctx)); 4709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5, &ctx.initialFuncs)); 4719566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &ctx)); 4729566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 4739566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "u")); 4749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 4759566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r, "r")); 476c4762a1bSJed Brown /* create TS */ 4779566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 4789566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm)); 4799566063dSJacob Faibussowitsch PetscCall(TSSetApplicationContext(ts, &ctx)); 4809566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx)); 481c4762a1bSJed Brown if (ctx.implicit) { 4829566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx)); 4839566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx)); 484c4762a1bSJed Brown } else { 4859566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &ctx)); 486c4762a1bSJed Brown } 4879566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 4889566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 4899566063dSJacob Faibussowitsch PetscCall(TSSetPostStep(ts, PostStep)); 490c4762a1bSJed Brown /* make solution & solve */ 4919566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, t, ctx.initialFuncs, (void **)ctxarr, INSERT_ALL_VALUES, u)); 4929566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, u)); 4939566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 4949566063dSJacob Faibussowitsch PetscCall(PostStep(ts)); /* print the initial state */ 4959566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 4969566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 4979566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, t, ctx.initialFuncs, (void **)ctxarr, u, &L2error)); 4989566063dSJacob Faibussowitsch if (L2error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n")); 49963a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", (double)L2error)); 5009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 5019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 5029566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 5039566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5049566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx.initialFuncs)); 5059566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 506b122ec5aSJacob Faibussowitsch return 0; 507c4762a1bSJed Brown } 508c4762a1bSJed Brown 509c4762a1bSJed Brown /*TEST 510c4762a1bSJed Brown 511c4762a1bSJed Brown test: 512c4762a1bSJed Brown suffix: 0 51330602db0SMatthew G. Knepley 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 \ 51430602db0SMatthew G. Knepley -ts_max_steps 1 -ts_max_time 10. -ts_dt 1.0 515c4762a1bSJed Brown test: 51661a622f3SMatthew G. Knepley # Remapping with periodicity is broken 517c4762a1bSJed Brown suffix: 1 51861a622f3SMatthew G. Knepley 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 \ 51930602db0SMatthew G. Knepley -ts_max_steps 1 -ts_max_time 10. -ts_dt 1.0 520c4762a1bSJed Brown 521c4762a1bSJed Brown TEST*/ 522