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 PetscFunctionBeginUser; 168 options->debug = 1; 169 options->plotRef = PETSC_FALSE; 170 options->implicit = PETSC_FALSE; 171 options->mu = 0; 172 options->eta = 0; 173 options->beta = 1; 174 options->a = 1; 175 options->b = PETSC_PI; 176 options->Jop = 0; 177 options->m = 1; 178 options->eps = 1.e-6; 179 180 PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX"); 181 PetscCall(PetscOptionsInt("-debug", "The debugging level", "ex48.c", options->debug, &options->debug, NULL)); 182 PetscCall(PetscOptionsBool("-plot_ref", "Plot the reference fields", "ex48.c", options->plotRef, &options->plotRef, NULL)); 183 PetscCall(PetscOptionsReal("-mu", "mu", "ex48.c", options->mu, &options->mu, NULL)); 184 PetscCall(PetscOptionsReal("-eta", "eta", "ex48.c", options->eta, &options->eta, NULL)); 185 PetscCall(PetscOptionsReal("-beta", "beta", "ex48.c", options->beta, &options->beta, NULL)); 186 PetscCall(PetscOptionsReal("-Jop", "Jop", "ex48.c", options->Jop, &options->Jop, NULL)); 187 PetscCall(PetscOptionsReal("-m", "m", "ex48.c", options->m, &options->m, NULL)); 188 PetscCall(PetscOptionsReal("-eps", "eps", "ex48.c", options->eps, &options->eps, NULL)); 189 PetscCall(PetscOptionsBool("-implicit", "Use implicit time integrator", "ex48.c", options->implicit, &options->implicit, NULL)); 190 PetscOptionsEnd(); 191 options->ke = PetscSqrtScalar(options->Jop); 192 if (options->Jop==0.0) { 193 options->Jo = 1.0/PetscPowScalar(options->a,2); 194 } else { 195 options->Jo = options->Jop*PetscCosReal(options->ke*options->a)/(1.0-PetscCosReal(options->ke*options->a)); 196 } 197 options->ky = PETSC_PI*options->m/options->b; 198 if (PetscPowReal(options->ky, 2) < options->Jop) { 199 options->kx = PetscSqrtScalar(options->Jop-PetscPowScalar(options->ky,2)); 200 options->DeltaPrime = -2.0*options->kx*options->a*PetscCosReal(options->kx*options->a)/PetscSinReal(options->kx*options->a); 201 } else if (PetscPowReal(options->ky, 2) > options->Jop) { 202 options->kx = PetscSqrtScalar(PetscPowScalar(options->ky,2)-options->Jop); 203 options->DeltaPrime = -2.0*options->kx*options->a*PetscCoshReal(options->kx*options->a)/PetscSinhReal(options->kx*options->a); 204 } else { /*they're equal (or there's a NaN), lim(x*cot(x))_x->0=1*/ 205 options->kx = 0; 206 options->DeltaPrime = -2.0; 207 } 208 PetscCall(PetscPrintf(comm, "DeltaPrime=%g\n",(double)options->DeltaPrime)); 209 210 PetscFunctionReturn(0); 211 } 212 213 static void f_n(PetscInt dim, PetscInt Nf, PetscInt NfAux, 214 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 215 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 216 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 217 { 218 const PetscScalar *pn = &u[uOff[DENSITY]]; 219 *f0 = *pn; 220 } 221 222 static PetscErrorCode PostStep(TS ts) 223 { 224 DM dm; 225 AppCtx *ctx; 226 PetscInt stepi,num; 227 Vec X; 228 229 PetscFunctionBegin; 230 PetscCall(TSGetApplicationContext(ts, &ctx)); 231 if (ctx->debug<1) PetscFunctionReturn(0); 232 PetscCall(TSGetSolution(ts, &X)); 233 PetscCall(VecGetDM(X, &dm)); 234 PetscCall(TSGetStepNumber(ts, &stepi)); 235 PetscCall(DMGetOutputSequenceNumber(dm, &num, NULL)); 236 if (num < 0) PetscCall(DMSetOutputSequenceNumber(dm, 0, 0.0)); 237 PetscCall(PetscObjectSetName((PetscObject) X, "u")); 238 PetscCall(VecViewFromOptions(X, NULL, "-vec_view")); 239 /* print integrals */ 240 { 241 PetscDS prob; 242 DM plex; 243 PetscScalar den, tt[5]; 244 PetscCall(DMConvert(dm, DMPLEX, &plex)); 245 PetscCall(DMGetDS(plex, &prob)); 246 PetscCall(PetscDSSetObjective(prob, 0, &f_n)); 247 PetscCall(DMPlexComputeIntegralFEM(plex,X,tt,ctx)); 248 den = tt[0]; 249 PetscCall(DMDestroy(&plex)); 250 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "%" PetscInt_FMT ") total perturbed mass = %g\n", stepi, (double) PetscRealPart(den))); 251 } 252 PetscFunctionReturn(0); 253 } 254 255 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 256 { 257 PetscFunctionBeginUser; 258 PetscCall(DMCreate(comm, dm)); 259 PetscCall(DMSetType(*dm, DMPLEX)); 260 PetscCall(DMSetFromOptions(*dm)); 261 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 262 263 PetscCall(DMGetBoundingBox(*dm, ctx->lower, ctx->upper)); 264 ctx->a = (ctx->upper[0] - ctx->lower[0])/2.0; 265 ctx->b = (ctx->upper[1] - ctx->lower[1])/2.0; 266 PetscFunctionReturn(0); 267 } 268 269 static PetscErrorCode log_n_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 270 { 271 AppCtx *lctx = (AppCtx*)ctx; 272 u[0] = 2.*lctx->a + x[0]; 273 return 0; 274 } 275 276 static PetscErrorCode Omega_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 277 { 278 u[0] = 0.0; 279 return 0; 280 } 281 282 static PetscErrorCode psi_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 283 { 284 AppCtx *lctx = (AppCtx*)ctx; 285 /* This sets up a symmetrix By flux aroound the mid point in x, which represents a current density flux along z. The stability 286 is analytically known and reported in ProcessOptions. */ 287 if (lctx->ke!=0.0) { 288 u[0] = (PetscCosReal(lctx->ke*(x[0]-lctx->a))-PetscCosReal(lctx->ke*lctx->a))/(1.0-PetscCosReal(lctx->ke*lctx->a)); 289 } else { 290 u[0] = 1.0-PetscPowScalar((x[0]-lctx->a)/lctx->a,2); 291 } 292 return 0; 293 } 294 295 static PetscErrorCode initialSolution_n(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 296 { 297 u[0] = 0.0; 298 return 0; 299 } 300 301 static PetscErrorCode initialSolution_Omega(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 302 { 303 u[0] = 0.0; 304 return 0; 305 } 306 307 static PetscErrorCode initialSolution_psi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *a_ctx) 308 { 309 AppCtx *ctx = (AppCtx*)a_ctx; 310 PetscScalar r = ctx->eps*(PetscScalar) (rand()) / (PetscScalar) (RAND_MAX); 311 if (x[0] == ctx->lower[0] || x[0] == ctx->upper[0]) r = 0; 312 u[0] = r; 313 /* PetscPrintf(PETSC_COMM_WORLD, "rand psi %lf\n",u[0]); */ 314 return 0; 315 } 316 317 static PetscErrorCode initialSolution_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 318 { 319 u[0] = 0.0; 320 return 0; 321 } 322 323 static PetscErrorCode initialSolution_jz(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 324 { 325 u[0] = 0.0; 326 return 0; 327 } 328 329 static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) 330 { 331 PetscDS ds; 332 DMLabel label; 333 const PetscInt id = 1; 334 335 PetscFunctionBeginUser; 336 PetscCall(DMGetLabel(dm, "marker", &label)); 337 PetscCall(DMGetDS(dm, &ds)); 338 PetscCall(PetscDSSetResidual(ds, 0, f0_n, f1_n)); 339 PetscCall(PetscDSSetResidual(ds, 1, f0_Omega, f1_Omega)); 340 PetscCall(PetscDSSetResidual(ds, 2, f0_psi, f1_psi)); 341 PetscCall(PetscDSSetResidual(ds, 3, f0_phi, f1_phi)); 342 PetscCall(PetscDSSetResidual(ds, 4, f0_jz, f1_jz)); 343 ctx->initialFuncs[0] = initialSolution_n; 344 ctx->initialFuncs[1] = initialSolution_Omega; 345 ctx->initialFuncs[2] = initialSolution_psi; 346 ctx->initialFuncs[3] = initialSolution_phi; 347 ctx->initialFuncs[4] = initialSolution_jz; 348 for (PetscInt f = 0; f < 5; ++f) { 349 PetscCall(PetscDSSetImplicit(ds, f, ctx->implicit)); 350 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, f, 0, NULL, (void (*)(void)) ctx->initialFuncs[f], NULL, ctx, NULL)); 351 } 352 PetscCall(PetscDSSetContext(ds, 0, ctx)); 353 PetscFunctionReturn(0); 354 } 355 356 static PetscErrorCode SetupEquilibriumFields(DM dm, DM dmAux, AppCtx *ctx) 357 { 358 PetscErrorCode (*eqFuncs[3])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar [], void *) = {log_n_0, Omega_0, psi_0}; 359 Vec eq; 360 AppCtx *ctxarr[3]; 361 362 ctxarr[0] = ctxarr[1] = ctxarr[2] = ctx; /* each variable could have a different context */ 363 PetscFunctionBegin; 364 PetscCall(DMCreateLocalVector(dmAux, &eq)); 365 PetscCall(DMProjectFunctionLocal(dmAux, 0.0, eqFuncs, (void **)ctxarr, INSERT_ALL_VALUES, eq)); 366 PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, eq)); 367 if (ctx->plotRef) { /* plot reference functions */ 368 PetscViewer viewer = NULL; 369 PetscBool isHDF5,isVTK; 370 char buf[256]; 371 Vec global; 372 PetscInt dim; 373 374 PetscCall(DMGetDimension(dm, &dim)); 375 PetscCall(DMCreateGlobalVector(dmAux,&global)); 376 PetscCall(VecSet(global,.0)); /* BCs! */ 377 PetscCall(DMLocalToGlobalBegin(dmAux,eq,INSERT_VALUES,global)); 378 PetscCall(DMLocalToGlobalEnd(dmAux,eq,INSERT_VALUES,global)); 379 PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)dmAux),&viewer)); 380 #ifdef PETSC_HAVE_HDF5 381 PetscCall(PetscViewerSetType(viewer,PETSCVIEWERHDF5)); 382 #else 383 PetscCall(PetscViewerSetType(viewer,PETSCVIEWERVTK)); 384 #endif 385 PetscCall(PetscViewerSetFromOptions(viewer)); 386 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5)); 387 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK)); 388 if (isHDF5) { 389 PetscCall(PetscSNPrintf(buf, 256, "uEquilibrium-%" PetscInt_FMT "D.h5", dim)); 390 } else if (isVTK) { 391 PetscCall(PetscSNPrintf(buf, 256, "uEquilibrium-%" PetscInt_FMT "D.vtu", dim)); 392 PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU)); 393 } 394 PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE)); 395 PetscCall(PetscViewerFileSetName(viewer,buf)); 396 if (isHDF5) PetscCall(DMView(dmAux,viewer)); 397 /* view equilibrium fields, this will overwrite fine grids with coarse grids! */ 398 PetscCall(PetscObjectSetName((PetscObject) global, "u0")); 399 PetscCall(VecView(global,viewer)); 400 PetscCall(PetscViewerDestroy(&viewer)); 401 PetscCall(VecDestroy(&global)); 402 } 403 PetscCall(VecDestroy(&eq)); 404 PetscFunctionReturn(0); 405 } 406 407 static PetscErrorCode SetupAuxDM(DM dm, PetscInt NfAux, PetscFE feAux[], AppCtx *user) 408 { 409 DM dmAux, coordDM; 410 PetscInt f; 411 412 PetscFunctionBegin; 413 /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */ 414 PetscCall(DMGetCoordinateDM(dm, &coordDM)); 415 if (!feAux) PetscFunctionReturn(0); 416 PetscCall(DMClone(dm, &dmAux)); 417 PetscCall(DMSetCoordinateDM(dmAux, coordDM)); 418 for (f = 0; f < NfAux; ++f) PetscCall(DMSetField(dmAux, f, NULL, (PetscObject) feAux[f])); 419 PetscCall(DMCreateDS(dmAux)); 420 PetscCall(SetupEquilibriumFields(dm, dmAux, user)); 421 PetscCall(DMDestroy(&dmAux)); 422 PetscFunctionReturn(0); 423 } 424 425 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx) 426 { 427 DM cdm = dm; 428 PetscFE fe[5], feAux[3]; 429 PetscInt dim, Nf = 5, NfAux = 3, f; 430 PetscBool simplex; 431 MPI_Comm comm; 432 433 PetscFunctionBeginUser; 434 /* Create finite element */ 435 PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 436 PetscCall(DMGetDimension(dm, &dim)); 437 PetscCall(DMPlexIsSimplex(dm, &simplex)); 438 PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[0])); 439 PetscCall(PetscObjectSetName((PetscObject) fe[0], "density")); 440 PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[1])); 441 PetscCall(PetscObjectSetName((PetscObject) fe[1], "vorticity")); 442 PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 443 PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[2])); 444 PetscCall(PetscObjectSetName((PetscObject) fe[2], "flux")); 445 PetscCall(PetscFECopyQuadrature(fe[0], fe[2])); 446 PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[3])); 447 PetscCall(PetscObjectSetName((PetscObject) fe[3], "potential")); 448 PetscCall(PetscFECopyQuadrature(fe[0], fe[3])); 449 PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe[4])); 450 PetscCall(PetscObjectSetName((PetscObject) fe[4], "current")); 451 PetscCall(PetscFECopyQuadrature(fe[0], fe[4])); 452 453 PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &feAux[0])); 454 PetscCall(PetscObjectSetName((PetscObject) feAux[0], "n_0")); 455 PetscCall(PetscFECopyQuadrature(fe[0], feAux[0])); 456 PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &feAux[1])); 457 PetscCall(PetscObjectSetName((PetscObject) feAux[1], "vorticity_0")); 458 PetscCall(PetscFECopyQuadrature(fe[0], feAux[1])); 459 PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &feAux[2])); 460 PetscCall(PetscObjectSetName((PetscObject) feAux[2], "flux_0")); 461 PetscCall(PetscFECopyQuadrature(fe[0], feAux[2])); 462 /* Set discretization and boundary conditions for each mesh */ 463 for (f = 0; f < Nf; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject) fe[f])); 464 PetscCall(DMCreateDS(dm)); 465 PetscCall(SetupProblem(dm, ctx)); 466 while (cdm) { 467 PetscCall(SetupAuxDM(dm, NfAux, feAux, ctx)); 468 PetscCall(DMCopyDisc(dm, cdm)); 469 PetscCall(DMGetCoarseDM(cdm, &cdm)); 470 } 471 for (f = 0; f < Nf; ++f) PetscCall(PetscFEDestroy(&fe[f])); 472 for (f = 0; f < NfAux; ++f) PetscCall(PetscFEDestroy(&feAux[f])); 473 PetscFunctionReturn(0); 474 } 475 476 int main(int argc, char **argv) 477 { 478 DM dm; 479 TS ts; 480 Vec u, r; 481 AppCtx ctx; 482 PetscReal t = 0.0; 483 PetscReal L2error = 0.0; 484 AppCtx *ctxarr[5]; 485 486 ctxarr[0] = ctxarr[1] = ctxarr[2] = ctxarr[3] = ctxarr[4] = &ctx; /* each variable could have a different context */ 487 s_ctx = &ctx; 488 PetscFunctionBeginUser; 489 PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 490 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 491 /* create mesh and problem */ 492 PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm)); 493 PetscCall(DMSetApplicationContext(dm, &ctx)); 494 PetscCall(PetscMalloc1(5, &ctx.initialFuncs)); 495 PetscCall(SetupDiscretization(dm, &ctx)); 496 PetscCall(DMCreateGlobalVector(dm, &u)); 497 PetscCall(PetscObjectSetName((PetscObject) u, "u")); 498 PetscCall(VecDuplicate(u, &r)); 499 PetscCall(PetscObjectSetName((PetscObject) r, "r")); 500 /* create TS */ 501 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 502 PetscCall(TSSetDM(ts, dm)); 503 PetscCall(TSSetApplicationContext(ts, &ctx)); 504 PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx)); 505 if (ctx.implicit) { 506 PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx)); 507 PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx)); 508 } else { 509 PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &ctx)); 510 } 511 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 512 PetscCall(TSSetFromOptions(ts)); 513 PetscCall(TSSetPostStep(ts, PostStep)); 514 /* make solution & solve */ 515 PetscCall(DMProjectFunction(dm, t, ctx.initialFuncs, (void **)ctxarr, INSERT_ALL_VALUES, u)); 516 PetscCall(TSSetSolution(ts,u)); 517 PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 518 PetscCall(PostStep(ts)); /* print the initial state */ 519 PetscCall(TSSolve(ts, u)); 520 PetscCall(TSGetTime(ts, &t)); 521 PetscCall(DMComputeL2Diff(dm, t, ctx.initialFuncs, (void **)ctxarr, u, &L2error)); 522 if (L2error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n")); 523 else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", (double)L2error)); 524 PetscCall(VecDestroy(&u)); 525 PetscCall(VecDestroy(&r)); 526 PetscCall(TSDestroy(&ts)); 527 PetscCall(DMDestroy(&dm)); 528 PetscCall(PetscFree(ctx.initialFuncs)); 529 PetscCall(PetscFinalize()); 530 return 0; 531 } 532 533 /*TEST 534 535 test: 536 suffix: 0 537 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 \ 538 -ts_max_steps 1 -ts_max_time 10. -ts_dt 1.0 539 test: 540 # Remapping with periodicity is broken 541 suffix: 1 542 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 \ 543 -ts_max_steps 1 -ts_max_time 10. -ts_dt 1.0 544 545 TEST*/ 546