1a4e36a32SMatthew G. Knepley static char help[] = "Time-dependent reactive low Mach Flow in 2d and 3d channels with finite elements.\n\ 2a4e36a32SMatthew G. Knepley We solve the reactive low Mach flow problem in a rectangular domain\n\ 3a4e36a32SMatthew G. Knepley using a parallel unstructured mesh (DMPLEX) to discretize the flow\n\ 4a4e36a32SMatthew G. Knepley and particles (DWSWARM) to discretize the chemical species.\n\n\n"; 5a4e36a32SMatthew G. Knepley 6a4e36a32SMatthew G. Knepley /*F 7a4e36a32SMatthew G. Knepley This low Mach flow is time-dependent isoviscous Navier-Stokes flow. We discretize using the 8a4e36a32SMatthew G. Knepley finite element method on an unstructured mesh. The weak form equations are 9a4e36a32SMatthew G. Knepley 10a4e36a32SMatthew G. Knepley \begin{align*} 11a4e36a32SMatthew G. Knepley < q, \nabla\cdot u > = 0 12a4e36a32SMatthew G. Knepley <v, du/dt> + <v, u \cdot \nabla u> + < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p > - < v, f > = 0 13a4e36a32SMatthew G. Knepley < w, u \cdot \nabla T > + < \nabla w, \alpha \nabla T > - < w, Q > = 0 14a4e36a32SMatthew G. Knepley \end{align*} 15a4e36a32SMatthew G. Knepley 16a4e36a32SMatthew G. Knepley where $\nu$ is the kinematic viscosity and $\alpha$ is thermal diffusivity. 17a4e36a32SMatthew G. Knepley 18a4e36a32SMatthew G. Knepley For visualization, use 19a4e36a32SMatthew G. Knepley 20a4e36a32SMatthew G. Knepley -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append 21a4e36a32SMatthew G. Knepley 22a4e36a32SMatthew G. Knepley The particles can be visualized using 23a4e36a32SMatthew G. Knepley 24a4e36a32SMatthew G. Knepley -part_dm_view draw -part_dm_view_swarm_radius 0.03 25a4e36a32SMatthew G. Knepley 26a4e36a32SMatthew G. Knepley F*/ 27a4e36a32SMatthew G. Knepley 28a4e36a32SMatthew G. Knepley #include <petscdmplex.h> 29a4e36a32SMatthew G. Knepley #include <petscdmswarm.h> 30a4e36a32SMatthew G. Knepley #include <petscts.h> 31a4e36a32SMatthew G. Knepley #include <petscds.h> 32a4e36a32SMatthew G. Knepley #include <petscbag.h> 33a4e36a32SMatthew G. Knepley 34*9371c9d4SSatish Balay typedef enum { 35*9371c9d4SSatish Balay SOL_TRIG_TRIG, 36*9371c9d4SSatish Balay NUM_SOL_TYPES 37*9371c9d4SSatish Balay } SolType; 38a4e36a32SMatthew G. Knepley const char *solTypes[NUM_SOL_TYPES + 1] = {"trig_trig", "unknown"}; 39a4e36a32SMatthew G. Knepley 40*9371c9d4SSatish Balay typedef enum { 41*9371c9d4SSatish Balay PART_LAYOUT_CELL, 42*9371c9d4SSatish Balay PART_LAYOUT_BOX, 43*9371c9d4SSatish Balay NUM_PART_LAYOUT_TYPES 44*9371c9d4SSatish Balay } PartLayoutType; 45a4e36a32SMatthew G. Knepley const char *partLayoutTypes[NUM_PART_LAYOUT_TYPES + 1] = {"cell", "box", "unknown"}; 46a4e36a32SMatthew G. Knepley 47a4e36a32SMatthew G. Knepley typedef struct { 48a4e36a32SMatthew G. Knepley PetscReal nu; /* Kinematic viscosity */ 49a4e36a32SMatthew G. Knepley PetscReal alpha; /* Thermal diffusivity */ 50a4e36a32SMatthew G. Knepley PetscReal T_in; /* Inlet temperature*/ 51a4e36a32SMatthew G. Knepley PetscReal omega; /* Rotation speed in MMS benchmark */ 52a4e36a32SMatthew G. Knepley } Parameter; 53a4e36a32SMatthew G. Knepley 54a4e36a32SMatthew G. Knepley typedef struct { 55a4e36a32SMatthew G. Knepley /* Problem definition */ 56a4e36a32SMatthew G. Knepley PetscBag bag; /* Holds problem parameters */ 57a4e36a32SMatthew G. Knepley SolType solType; /* MMS solution type */ 58a4e36a32SMatthew G. Knepley PartLayoutType partLayout; /* Type of particle distribution */ 59a4e36a32SMatthew G. Knepley PetscInt Npc; /* The initial number of particles per cell */ 60a4e36a32SMatthew G. Knepley PetscReal partLower[3]; /* Lower left corner of particle box */ 61a4e36a32SMatthew G. Knepley PetscReal partUpper[3]; /* Upper right corner of particle box */ 62a4e36a32SMatthew G. Knepley PetscInt Npb; /* The initial number of particles per box dimension */ 63a4e36a32SMatthew G. Knepley } AppCtx; 64a4e36a32SMatthew G. Knepley 65a4e36a32SMatthew G. Knepley typedef struct { 66a4e36a32SMatthew G. Knepley PetscReal ti; /* The time for ui, at the beginning of the advection solve */ 67a4e36a32SMatthew G. Knepley PetscReal tf; /* The time for uf, at the end of the advection solve */ 68a4e36a32SMatthew G. Knepley Vec ui; /* The PDE solution field at ti */ 69a4e36a32SMatthew G. Knepley Vec uf; /* The PDE solution field at tf */ 70a4e36a32SMatthew G. Knepley Vec x0; /* The initial particle positions at t = 0 */ 71a4e36a32SMatthew G. Knepley PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 72a4e36a32SMatthew G. Knepley AppCtx *ctx; /* Context for exact solution */ 73a4e36a32SMatthew G. Knepley } AdvCtx; 74a4e36a32SMatthew G. Knepley 75*9371c9d4SSatish Balay static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 76a4e36a32SMatthew G. Knepley PetscInt d; 77a4e36a32SMatthew G. Knepley for (d = 0; d < Nc; ++d) u[d] = 0.0; 78a4e36a32SMatthew G. Knepley return 0; 79a4e36a32SMatthew G. Knepley } 80a4e36a32SMatthew G. Knepley 81*9371c9d4SSatish Balay static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 82a4e36a32SMatthew G. Knepley PetscInt d; 83a4e36a32SMatthew G. Knepley for (d = 0; d < Nc; ++d) u[d] = 1.0; 84a4e36a32SMatthew G. Knepley return 0; 85a4e36a32SMatthew G. Knepley } 86a4e36a32SMatthew G. Knepley 87a4e36a32SMatthew G. Knepley /* 88a4e36a32SMatthew G. Knepley CASE: trigonometric-trigonometric 89a4e36a32SMatthew G. Knepley In 2D we use exact solution: 90a4e36a32SMatthew G. Knepley 91a4e36a32SMatthew G. Knepley x = r0 cos(w t + theta0) r0 = sqrt(x0^2 + y0^2) 92a4e36a32SMatthew G. Knepley y = r0 sin(w t + theta0) theta0 = arctan(y0/x0) 93a4e36a32SMatthew G. Knepley u = -w r0 sin(theta0) = -w y 94a4e36a32SMatthew G. Knepley v = w r0 cos(theta0) = w x 95a4e36a32SMatthew G. Knepley p = x + y - 1 96a4e36a32SMatthew G. Knepley T = t + x + y 97a4e36a32SMatthew G. Knepley f = <1, 1> 98a4e36a32SMatthew G. Knepley Q = 1 + w (x - y)/r 99a4e36a32SMatthew G. Knepley 100a4e36a32SMatthew G. Knepley so that 101a4e36a32SMatthew G. Knepley 102a4e36a32SMatthew G. Knepley \nabla \cdot u = 0 + 0 = 0 103a4e36a32SMatthew G. Knepley 104a4e36a32SMatthew G. Knepley f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p 105a4e36a32SMatthew G. Knepley = <0, 0> + u_i d_i u_j - \nu 0 + <1, 1> 106a4e36a32SMatthew G. Knepley = <1, 1> + w^2 <-y, x> . <<0, 1>, <-1, 0>> 107a4e36a32SMatthew G. Knepley = <1, 1> + w^2 <-x, -y> 108a4e36a32SMatthew G. Knepley = <1, 1> - w^2 <x, y> 109a4e36a32SMatthew G. Knepley 110a4e36a32SMatthew G. Knepley Q = dT/dt + u \cdot \nabla T - \alpha \Delta T 111a4e36a32SMatthew G. Knepley = 1 + <u, v> . <1, 1> - \alpha 0 112a4e36a32SMatthew G. Knepley = 1 + u + v 113a4e36a32SMatthew G. Knepley */ 114*9371c9d4SSatish Balay static PetscErrorCode trig_trig_x(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *x, void *ctx) { 115a4e36a32SMatthew G. Knepley const PetscReal x0 = X[0]; 116a4e36a32SMatthew G. Knepley const PetscReal y0 = X[1]; 117a4e36a32SMatthew G. Knepley const PetscReal R0 = PetscSqrtReal(x0 * x0 + y0 * y0); 118a4e36a32SMatthew G. Knepley const PetscReal theta0 = PetscAtan2Real(y0, x0); 119a4e36a32SMatthew G. Knepley Parameter *p = (Parameter *)ctx; 120a4e36a32SMatthew G. Knepley 121a4e36a32SMatthew G. Knepley x[0] = R0 * PetscCosReal(p->omega * time + theta0); 122a4e36a32SMatthew G. Knepley x[1] = R0 * PetscSinReal(p->omega * time + theta0); 123a4e36a32SMatthew G. Knepley return 0; 124a4e36a32SMatthew G. Knepley } 125*9371c9d4SSatish Balay static PetscErrorCode trig_trig_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) { 126a4e36a32SMatthew G. Knepley Parameter *p = (Parameter *)ctx; 127a4e36a32SMatthew G. Knepley 128a4e36a32SMatthew G. Knepley u[0] = -p->omega * X[1]; 129a4e36a32SMatthew G. Knepley u[1] = p->omega * X[0]; 130a4e36a32SMatthew G. Knepley return 0; 131a4e36a32SMatthew G. Knepley } 132*9371c9d4SSatish Balay static PetscErrorCode trig_trig_u_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) { 133a4e36a32SMatthew G. Knepley u[0] = 0.0; 134a4e36a32SMatthew G. Knepley u[1] = 0.0; 135a4e36a32SMatthew G. Knepley return 0; 136a4e36a32SMatthew G. Knepley } 137a4e36a32SMatthew G. Knepley 138*9371c9d4SSatish Balay static PetscErrorCode trig_trig_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) { 139a4e36a32SMatthew G. Knepley p[0] = X[0] + X[1] - 1.0; 140a4e36a32SMatthew G. Knepley return 0; 141a4e36a32SMatthew G. Knepley } 142a4e36a32SMatthew G. Knepley 143*9371c9d4SSatish Balay static PetscErrorCode trig_trig_T(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) { 144a4e36a32SMatthew G. Knepley T[0] = time + X[0] + X[1]; 145a4e36a32SMatthew G. Knepley return 0; 146a4e36a32SMatthew G. Knepley } 147*9371c9d4SSatish Balay static PetscErrorCode trig_trig_T_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) { 148a4e36a32SMatthew G. Knepley T[0] = 1.0; 149a4e36a32SMatthew G. Knepley return 0; 150a4e36a32SMatthew G. Knepley } 151a4e36a32SMatthew G. Knepley 152*9371c9d4SSatish Balay static void f0_trig_trig_v(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[]) { 153a4e36a32SMatthew G. Knepley const PetscReal omega = PetscRealPart(constants[3]); 154a4e36a32SMatthew G. Knepley PetscInt Nc = dim; 155a4e36a32SMatthew G. Knepley PetscInt c, d; 156a4e36a32SMatthew G. Knepley 157a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[d] = u_t[uOff[0] + d]; 158a4e36a32SMatthew G. Knepley 159a4e36a32SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 160a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[c] += u[d] * u_x[c * dim + d]; 161a4e36a32SMatthew G. Knepley } 162a4e36a32SMatthew G. Knepley f0[0] -= 1.0 - omega * omega * X[0]; 163a4e36a32SMatthew G. Knepley f0[1] -= 1.0 - omega * omega * X[1]; 164a4e36a32SMatthew G. Knepley } 165a4e36a32SMatthew G. Knepley 166*9371c9d4SSatish Balay static void f0_trig_trig_w(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[]) { 167a4e36a32SMatthew G. Knepley const PetscReal omega = PetscRealPart(constants[3]); 168a4e36a32SMatthew G. Knepley PetscInt d; 169a4e36a32SMatthew G. Knepley 170a4e36a32SMatthew G. Knepley for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0] + d] * u_x[uOff_x[2] + d]; 171a4e36a32SMatthew G. Knepley f0[0] += u_t[uOff[2]] - (1.0 + omega * (X[0] - X[1])); 172a4e36a32SMatthew G. Knepley } 173a4e36a32SMatthew G. Knepley 174*9371c9d4SSatish Balay static void f0_q(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[]) { 175a4e36a32SMatthew G. Knepley PetscInt d; 176a4e36a32SMatthew G. Knepley for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d * dim + d]; 177a4e36a32SMatthew G. Knepley } 178a4e36a32SMatthew G. Knepley 179a4e36a32SMatthew G. Knepley /*f1_v = \nu[grad(u) + grad(u)^T] - pI */ 180*9371c9d4SSatish Balay static void f1_v(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[]) { 181a4e36a32SMatthew G. Knepley const PetscReal nu = PetscRealPart(constants[0]); 182a4e36a32SMatthew G. Knepley const PetscInt Nc = dim; 183a4e36a32SMatthew G. Knepley PetscInt c, d; 184a4e36a32SMatthew G. Knepley 185a4e36a32SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 186*9371c9d4SSatish Balay for (d = 0; d < dim; ++d) { f1[c * dim + d] = nu * (u_x[c * dim + d] + u_x[d * dim + c]); } 187a4e36a32SMatthew G. Knepley f1[c * dim + c] -= u[uOff[1]]; 188a4e36a32SMatthew G. Knepley } 189a4e36a32SMatthew G. Knepley } 190a4e36a32SMatthew G. Knepley 191*9371c9d4SSatish Balay static void f1_w(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[]) { 192a4e36a32SMatthew G. Knepley const PetscReal alpha = PetscRealPart(constants[1]); 193a4e36a32SMatthew G. Knepley PetscInt d; 194a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = alpha * u_x[uOff_x[2] + d]; 195a4e36a32SMatthew G. Knepley } 196a4e36a32SMatthew G. Knepley 197a4e36a32SMatthew G. Knepley /* Jacobians */ 198*9371c9d4SSatish Balay static void g1_qu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) { 199a4e36a32SMatthew G. Knepley PetscInt d; 200a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; 201a4e36a32SMatthew G. Knepley } 202a4e36a32SMatthew G. Knepley 203*9371c9d4SSatish Balay static void g0_vu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) { 204a4e36a32SMatthew G. Knepley PetscInt c, d; 205a4e36a32SMatthew G. Knepley const PetscInt Nc = dim; 206a4e36a32SMatthew G. Knepley 207a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g0[d * dim + d] = u_tShift; 208a4e36a32SMatthew G. Knepley 209a4e36a32SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 210*9371c9d4SSatish Balay for (d = 0; d < dim; ++d) { g0[c * Nc + d] += u_x[c * Nc + d]; } 211a4e36a32SMatthew G. Knepley } 212a4e36a32SMatthew G. Knepley } 213a4e36a32SMatthew G. Knepley 214*9371c9d4SSatish Balay static void g1_vu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) { 215a4e36a32SMatthew G. Knepley PetscInt NcI = dim; 216a4e36a32SMatthew G. Knepley PetscInt NcJ = dim; 217a4e36a32SMatthew G. Knepley PetscInt c, d, e; 218a4e36a32SMatthew G. Knepley 219a4e36a32SMatthew G. Knepley for (c = 0; c < NcI; ++c) { 220a4e36a32SMatthew G. Knepley for (d = 0; d < NcJ; ++d) { 221a4e36a32SMatthew G. Knepley for (e = 0; e < dim; ++e) { 222*9371c9d4SSatish Balay if (c == d) { g1[(c * NcJ + d) * dim + e] += u[e]; } 223a4e36a32SMatthew G. Knepley } 224a4e36a32SMatthew G. Knepley } 225a4e36a32SMatthew G. Knepley } 226a4e36a32SMatthew G. Knepley } 227a4e36a32SMatthew G. Knepley 228*9371c9d4SSatish Balay static void g2_vp(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) { 229a4e36a32SMatthew G. Knepley PetscInt d; 230a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0; 231a4e36a32SMatthew G. Knepley } 232a4e36a32SMatthew G. Knepley 233*9371c9d4SSatish Balay static void g3_vu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) { 234a4e36a32SMatthew G. Knepley const PetscReal nu = PetscRealPart(constants[0]); 235a4e36a32SMatthew G. Knepley const PetscInt Nc = dim; 236a4e36a32SMatthew G. Knepley PetscInt c, d; 237a4e36a32SMatthew G. Knepley 238a4e36a32SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 239a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) { 240a4e36a32SMatthew G. Knepley g3[((c * Nc + c) * dim + d) * dim + d] += nu; // gradU 241a4e36a32SMatthew G. Knepley g3[((c * Nc + d) * dim + d) * dim + c] += nu; // gradU transpose 242a4e36a32SMatthew G. Knepley } 243a4e36a32SMatthew G. Knepley } 244a4e36a32SMatthew G. Knepley } 245a4e36a32SMatthew G. Knepley 246*9371c9d4SSatish Balay static void g0_wT(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) { 247a4e36a32SMatthew G. Knepley PetscInt d; 248a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g0[d] = u_tShift; 249a4e36a32SMatthew G. Knepley } 250a4e36a32SMatthew G. Knepley 251*9371c9d4SSatish Balay static void g0_wu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) { 252a4e36a32SMatthew G. Knepley PetscInt d; 253a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2] + d]; 254a4e36a32SMatthew G. Knepley } 255a4e36a32SMatthew G. Knepley 256*9371c9d4SSatish Balay static void g1_wT(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) { 257a4e36a32SMatthew G. Knepley PetscInt d; 258a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g1[d] = u[uOff[0] + d]; 259a4e36a32SMatthew G. Knepley } 260a4e36a32SMatthew G. Knepley 261*9371c9d4SSatish Balay static void g3_wT(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) { 262a4e36a32SMatthew G. Knepley const PetscReal alpha = PetscRealPart(constants[1]); 263a4e36a32SMatthew G. Knepley PetscInt d; 264a4e36a32SMatthew G. Knepley 265a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d * dim + d] = alpha; 266a4e36a32SMatthew G. Knepley } 267a4e36a32SMatthew G. Knepley 268*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 269a4e36a32SMatthew G. Knepley PetscInt sol, pl, n; 270a4e36a32SMatthew G. Knepley 271a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 272a4e36a32SMatthew G. Knepley options->solType = SOL_TRIG_TRIG; 273a4e36a32SMatthew G. Knepley options->partLayout = PART_LAYOUT_CELL; 274a4e36a32SMatthew G. Knepley options->Npc = 1; 275a4e36a32SMatthew G. Knepley options->Npb = 1; 276a4e36a32SMatthew G. Knepley options->partLower[0] = options->partLower[1] = options->partLower[2] = 0.; 277a4e36a32SMatthew G. Knepley options->partUpper[0] = options->partUpper[1] = options->partUpper[2] = 1.; 278d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX"); 279a4e36a32SMatthew G. Knepley sol = options->solType; 2809566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-sol_type", "The solution type", "ex77.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL)); 281a4e36a32SMatthew G. Knepley options->solType = (SolType)sol; 282a4e36a32SMatthew G. Knepley pl = options->partLayout; 2839566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-part_layout_type", "The particle layout type", "ex77.c", partLayoutTypes, NUM_PART_LAYOUT_TYPES, partLayoutTypes[options->partLayout], &pl, NULL)); 284a4e36a32SMatthew G. Knepley options->partLayout = (PartLayoutType)pl; 2859566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-Npc", "The initial number of particles per cell", "ex77.c", options->Npc, &options->Npc, NULL)); 286a4e36a32SMatthew G. Knepley n = 3; 2879566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-part_lower", "The lower left corner of the particle box", "ex77.c", options->partLower, &n, NULL)); 288a4e36a32SMatthew G. Knepley n = 3; 2899566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-part_upper", "The upper right corner of the particle box", "ex77.c", options->partUpper, &n, NULL)); 2909566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-Npb", "The initial number of particles per box dimension", "ex77.c", options->Npb, &options->Npb, NULL)); 291d0609cedSBarry Smith PetscOptionsEnd(); 292a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 293a4e36a32SMatthew G. Knepley } 294a4e36a32SMatthew G. Knepley 295*9371c9d4SSatish Balay static PetscErrorCode SetupParameters(AppCtx *user) { 296a4e36a32SMatthew G. Knepley PetscBag bag; 297a4e36a32SMatthew G. Knepley Parameter *p; 298a4e36a32SMatthew G. Knepley 299a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 300a4e36a32SMatthew G. Knepley /* setup PETSc parameter bag */ 3019566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)&p)); 3029566063dSJacob Faibussowitsch PetscCall(PetscBagSetName(user->bag, "par", "Low Mach flow parameters")); 303a4e36a32SMatthew G. Knepley bag = user->bag; 3049566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity")); 3059566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity")); 3069566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->T_in, 1.0, "T_in", "Inlet temperature")); 3079566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->omega, 1.0, "omega", "Rotation speed in MMS benchmark")); 308a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 309a4e36a32SMatthew G. Knepley } 310a4e36a32SMatthew G. Knepley 311*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 312a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 3139566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 3149566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 3159566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 3169566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 317a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 318a4e36a32SMatthew G. Knepley } 319a4e36a32SMatthew G. Knepley 320*9371c9d4SSatish Balay static PetscErrorCode SetupProblem(DM dm, AppCtx *user) { 321a4e36a32SMatthew G. Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 322a4e36a32SMatthew G. Knepley PetscErrorCode (*exactFuncs_t[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 323a4e36a32SMatthew G. Knepley PetscDS prob; 32445480ffeSMatthew G. Knepley DMLabel label; 325a4e36a32SMatthew G. Knepley Parameter *ctx; 326a4e36a32SMatthew G. Knepley PetscInt id; 327a4e36a32SMatthew G. Knepley 328a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 3299566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 3309566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 331a4e36a32SMatthew G. Knepley switch (user->solType) { 332a4e36a32SMatthew G. Knepley case SOL_TRIG_TRIG: 3339566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_trig_trig_v, f1_v)); 3349566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 2, f0_trig_trig_w, f1_w)); 335a4e36a32SMatthew G. Knepley 336a4e36a32SMatthew G. Knepley exactFuncs[0] = trig_trig_u; 337a4e36a32SMatthew G. Knepley exactFuncs[1] = trig_trig_p; 338a4e36a32SMatthew G. Knepley exactFuncs[2] = trig_trig_T; 339a4e36a32SMatthew G. Knepley exactFuncs_t[0] = trig_trig_u_t; 340a4e36a32SMatthew G. Knepley exactFuncs_t[1] = NULL; 341a4e36a32SMatthew G. Knepley exactFuncs_t[2] = trig_trig_T_t; 342a4e36a32SMatthew G. Knepley break; 34363a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 344a4e36a32SMatthew G. Knepley } 345a4e36a32SMatthew G. Knepley 3469566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 1, f0_q, NULL)); 347a4e36a32SMatthew G. Knepley 3489566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu)); 3499566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL)); 3509566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL)); 3519566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL)); 3529566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 2, 2, g0_wT, g1_wT, NULL, g3_wT)); 353a4e36a32SMatthew G. Knepley /* Setup constants */ 354a4e36a32SMatthew G. Knepley { 355a4e36a32SMatthew G. Knepley Parameter *param; 356a4e36a32SMatthew G. Knepley PetscScalar constants[4]; 357a4e36a32SMatthew G. Knepley 3589566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 359a4e36a32SMatthew G. Knepley 360a4e36a32SMatthew G. Knepley constants[0] = param->nu; 361a4e36a32SMatthew G. Knepley constants[1] = param->alpha; 362a4e36a32SMatthew G. Knepley constants[2] = param->T_in; 363a4e36a32SMatthew G. Knepley constants[3] = param->omega; 3649566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 4, constants)); 365a4e36a32SMatthew G. Knepley } 366a4e36a32SMatthew G. Knepley /* Setup Boundary Conditions */ 3679566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)&ctx)); 368a4e36a32SMatthew G. Knepley id = 3; 3699566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], (void (*)(void))exactFuncs_t[0], ctx, NULL)); 370a4e36a32SMatthew G. Knepley id = 1; 3719566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], (void (*)(void))exactFuncs_t[0], ctx, NULL)); 372a4e36a32SMatthew G. Knepley id = 2; 3739566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], (void (*)(void))exactFuncs_t[0], ctx, NULL)); 374a4e36a32SMatthew G. Knepley id = 4; 3759566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], (void (*)(void))exactFuncs_t[0], ctx, NULL)); 376a4e36a32SMatthew G. Knepley id = 3; 3779566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void))exactFuncs[2], (void (*)(void))exactFuncs_t[2], ctx, NULL)); 378a4e36a32SMatthew G. Knepley id = 1; 3799566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void))exactFuncs[2], (void (*)(void))exactFuncs_t[2], ctx, NULL)); 380a4e36a32SMatthew G. Knepley id = 2; 3819566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void))exactFuncs[2], (void (*)(void))exactFuncs_t[2], ctx, NULL)); 382a4e36a32SMatthew G. Knepley id = 4; 3839566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void))exactFuncs[2], (void (*)(void))exactFuncs_t[2], ctx, NULL)); 384a4e36a32SMatthew G. Knepley 385a4e36a32SMatthew G. Knepley /*setup exact solution.*/ 3869566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx)); 3879566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx)); 3889566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx)); 3899566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(prob, 0, exactFuncs_t[0], ctx)); 3909566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(prob, 1, exactFuncs_t[1], ctx)); 3919566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(prob, 2, exactFuncs_t[2], ctx)); 392a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 393a4e36a32SMatthew G. Knepley } 394a4e36a32SMatthew G. Knepley 395a4e36a32SMatthew G. Knepley /* x_t = v 396a4e36a32SMatthew G. Knepley 397a4e36a32SMatthew G. Knepley Note that here we use the velocity field at t_{n+1} to advect the particles from 398a4e36a32SMatthew G. Knepley t_n to t_{n+1}. If we use both of these fields, we could use Crank-Nicholson or 399a4e36a32SMatthew G. Knepley the method of characteristics. 400a4e36a32SMatthew G. Knepley */ 401*9371c9d4SSatish Balay static PetscErrorCode FreeStreaming(TS ts, PetscReal t, Vec X, Vec F, void *ctx) { 402a4e36a32SMatthew G. Knepley AdvCtx *adv = (AdvCtx *)ctx; 403a4e36a32SMatthew G. Knepley Vec u = adv->ui; 404a4e36a32SMatthew G. Knepley DM sdm, dm, vdm; 405a4e36a32SMatthew G. Knepley Vec vel, locvel, pvel; 406a4e36a32SMatthew G. Knepley IS vis; 407a4e36a32SMatthew G. Knepley DMInterpolationInfo ictx; 408a4e36a32SMatthew G. Knepley const PetscScalar *coords, *v; 409a4e36a32SMatthew G. Knepley PetscScalar *f; 410a4e36a32SMatthew G. Knepley PetscInt vf[1] = {0}; 411a4e36a32SMatthew G. Knepley PetscInt dim, Np; 412a4e36a32SMatthew G. Knepley 413a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 4149566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &sdm)); 4159566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sdm, &dm)); 4169566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(sdm, &pvel)); 4179566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(sdm, &Np)); 4189566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 419a4e36a32SMatthew G. Knepley /* Get local velocity */ 4209566063dSJacob Faibussowitsch PetscCall(DMCreateSubDM(dm, 1, vf, &vis, &vdm)); 4219566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(u, vis, &vel)); 4229566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(vdm, &locvel)); 4239566063dSJacob Faibussowitsch PetscCall(DMPlexInsertBoundaryValues(vdm, PETSC_TRUE, locvel, adv->ti, NULL, NULL, NULL)); 4249566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(vdm, vel, INSERT_VALUES, locvel)); 4259566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(vdm, vel, INSERT_VALUES, locvel)); 4269566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(u, vis, &vel)); 4279566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vis)); 428a4e36a32SMatthew G. Knepley /* Interpolate velocity */ 4299566063dSJacob Faibussowitsch PetscCall(DMInterpolationCreate(PETSC_COMM_SELF, &ictx)); 4309566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetDim(ictx, dim)); 4319566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetDof(ictx, dim)); 4329566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &coords)); 4339566063dSJacob Faibussowitsch PetscCall(DMInterpolationAddPoints(ictx, Np, (PetscReal *)coords)); 4349566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &coords)); 43552aa1562SMatthew G. Knepley /* Particles that lie outside the domain should be dropped, 43652aa1562SMatthew G. Knepley whereas particles that move to another partition should trigger a migration */ 4379566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetUp(ictx, vdm, PETSC_FALSE, PETSC_TRUE)); 4389566063dSJacob Faibussowitsch PetscCall(VecSet(pvel, 0.)); 4399566063dSJacob Faibussowitsch PetscCall(DMInterpolationEvaluate(ictx, vdm, locvel, pvel)); 4409566063dSJacob Faibussowitsch PetscCall(DMInterpolationDestroy(&ictx)); 4419566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(vdm, &locvel)); 4429566063dSJacob Faibussowitsch PetscCall(DMDestroy(&vdm)); 443a4e36a32SMatthew G. Knepley 4449566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 4459566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pvel, &v)); 4469566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(f, v, Np * dim)); 4479566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pvel, &v)); 4489566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 4499566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(sdm, &pvel)); 450a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 451a4e36a32SMatthew G. Knepley } 452a4e36a32SMatthew G. Knepley 453*9371c9d4SSatish Balay static PetscErrorCode SetInitialParticleConditions(TS ts, Vec u) { 454a4e36a32SMatthew G. Knepley AppCtx *user; 455a4e36a32SMatthew G. Knepley void *ctx; 456a4e36a32SMatthew G. Knepley DM dm; 457a4e36a32SMatthew G. Knepley PetscScalar *coords; 458a4e36a32SMatthew G. Knepley PetscReal x[3], dx[3]; 459a4e36a32SMatthew G. Knepley PetscInt n[3]; 4605f80ce2aSJacob Faibussowitsch PetscInt dim, d, i, j, k; 461a4e36a32SMatthew G. Knepley 4627510d9b0SBarry Smith PetscFunctionBeginUser; 4639566063dSJacob Faibussowitsch PetscCall(TSGetApplicationContext(ts, &ctx)); 464a4e36a32SMatthew G. Knepley user = ((AdvCtx *)ctx)->ctx; 4659566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 4669566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 467a4e36a32SMatthew G. Knepley switch (user->partLayout) { 468*9371c9d4SSatish Balay case PART_LAYOUT_CELL: PetscCall(DMSwarmSetPointCoordinatesRandom(dm, user->Npc)); break; 469a4e36a32SMatthew G. Knepley case PART_LAYOUT_BOX: 470a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) { 471a4e36a32SMatthew G. Knepley n[d] = user->Npb; 472a4e36a32SMatthew G. Knepley dx[d] = (user->partUpper[d] - user->partLower[d]) / PetscMax(1, n[d] - 1); 473a4e36a32SMatthew G. Knepley } 4749566063dSJacob Faibussowitsch PetscCall(VecGetArray(u, &coords)); 475a4e36a32SMatthew G. Knepley switch (dim) { 476a4e36a32SMatthew G. Knepley case 2: 477a4e36a32SMatthew G. Knepley x[0] = user->partLower[0]; 478a4e36a32SMatthew G. Knepley for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 479a4e36a32SMatthew G. Knepley x[1] = user->partLower[1]; 480a4e36a32SMatthew G. Knepley for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 481a4e36a32SMatthew G. Knepley const PetscInt p = j * n[0] + i; 482a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) coords[p * dim + d] = x[d]; 483a4e36a32SMatthew G. Knepley } 484a4e36a32SMatthew G. Knepley } 485a4e36a32SMatthew G. Knepley break; 486a4e36a32SMatthew G. Knepley case 3: 487a4e36a32SMatthew G. Knepley x[0] = user->partLower[0]; 488a4e36a32SMatthew G. Knepley for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 489a4e36a32SMatthew G. Knepley x[1] = user->partLower[1]; 490a4e36a32SMatthew G. Knepley for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 491a4e36a32SMatthew G. Knepley x[2] = user->partLower[2]; 492a4e36a32SMatthew G. Knepley for (k = 0; k < n[2]; ++k, x[2] += dx[2]) { 493a4e36a32SMatthew G. Knepley const PetscInt p = (k * n[1] + j) * n[0] + i; 494a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) coords[p * dim + d] = x[d]; 495a4e36a32SMatthew G. Knepley } 496a4e36a32SMatthew G. Knepley } 497a4e36a32SMatthew G. Knepley } 498a4e36a32SMatthew G. Knepley break; 49963a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Do not support particle layout in dimension %" PetscInt_FMT, dim); 500a4e36a32SMatthew G. Knepley } 5019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(u, &coords)); 502a4e36a32SMatthew G. Knepley break; 50398921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]); 504a4e36a32SMatthew G. Knepley } 505a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 506a4e36a32SMatthew G. Knepley } 507a4e36a32SMatthew G. Knepley 508*9371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, DM sdm, AppCtx *user) { 509a4e36a32SMatthew G. Knepley DM cdm = dm; 510a4e36a32SMatthew G. Knepley PetscFE fe[3]; 511a4e36a32SMatthew G. Knepley Parameter *param; 512a4e36a32SMatthew G. Knepley PetscInt *cellid, n[3]; 513a4e36a32SMatthew G. Knepley PetscReal x[3], dx[3]; 514a4e36a32SMatthew G. Knepley PetscScalar *coords; 515a4e36a32SMatthew G. Knepley DMPolytopeType ct; 516a4e36a32SMatthew G. Knepley PetscInt dim, d, cStart, cEnd, c, Np, p, i, j, k; 517a4e36a32SMatthew G. Knepley PetscBool simplex; 518a4e36a32SMatthew G. Knepley MPI_Comm comm; 519a4e36a32SMatthew G. Knepley 520a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 5219566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 5229566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 5239566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 524a4e36a32SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 525a4e36a32SMatthew G. Knepley /* Create finite element */ 5269566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 5279566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0])); 5289566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity")); 529a4e36a32SMatthew G. Knepley 5309566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1])); 5319566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 5329566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure")); 533a4e36a32SMatthew G. Knepley 5349566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2])); 5359566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[2])); 5369566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[2], "temperature")); 537a4e36a32SMatthew G. Knepley 538a4e36a32SMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */ 5399566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0])); 5409566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1])); 5419566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 2, NULL, (PetscObject)fe[2])); 5429566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 5439566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user)); 5449566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 545a4e36a32SMatthew G. Knepley while (cdm) { 5469566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 5479566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 548a4e36a32SMatthew G. Knepley } 5499566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0])); 5509566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1])); 5519566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[2])); 552a4e36a32SMatthew G. Knepley 553a4e36a32SMatthew G. Knepley { 554a4e36a32SMatthew G. Knepley PetscObject pressure; 555a4e36a32SMatthew G. Knepley MatNullSpace nullspacePres; 556a4e36a32SMatthew G. Knepley 5579566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 1, NULL, &pressure)); 5589566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres)); 5599566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject)nullspacePres)); 5609566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nullspacePres)); 561a4e36a32SMatthew G. Knepley } 562a4e36a32SMatthew G. Knepley 563a4e36a32SMatthew G. Knepley /* Setup particle information */ 5649566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(sdm, DMSWARM_PIC)); 5659566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(sdm, "mass", 1, PETSC_REAL)); 5669566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(sdm)); 567a4e36a32SMatthew G. Knepley switch (user->partLayout) { 568a4e36a32SMatthew G. Knepley case PART_LAYOUT_CELL: 5699566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(sdm, (cEnd - cStart) * user->Npc, 0)); 5709566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(sdm)); 5719566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 572a4e36a32SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 573a4e36a32SMatthew G. Knepley for (p = 0; p < user->Npc; ++p) { 574a4e36a32SMatthew G. Knepley const PetscInt n = c * user->Npc + p; 575a4e36a32SMatthew G. Knepley 576a4e36a32SMatthew G. Knepley cellid[n] = c; 577a4e36a32SMatthew G. Knepley } 578a4e36a32SMatthew G. Knepley } 5799566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 5809566063dSJacob Faibussowitsch PetscCall(DMSwarmSetPointCoordinatesRandom(sdm, user->Npc)); 581a4e36a32SMatthew G. Knepley break; 582a4e36a32SMatthew G. Knepley case PART_LAYOUT_BOX: 583a4e36a32SMatthew G. Knepley Np = 1; 584a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) { 585a4e36a32SMatthew G. Knepley n[d] = user->Npb; 586a4e36a32SMatthew G. Knepley dx[d] = (user->partUpper[d] - user->partLower[d]) / PetscMax(1, n[d] - 1); 587a4e36a32SMatthew G. Knepley Np *= n[d]; 588a4e36a32SMatthew G. Knepley } 5899566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(sdm, Np, 0)); 5909566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(sdm)); 5919566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 592a4e36a32SMatthew G. Knepley switch (dim) { 593a4e36a32SMatthew G. Knepley case 2: 594a4e36a32SMatthew G. Knepley x[0] = user->partLower[0]; 595a4e36a32SMatthew G. Knepley for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 596a4e36a32SMatthew G. Knepley x[1] = user->partLower[1]; 597a4e36a32SMatthew G. Knepley for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 598a4e36a32SMatthew G. Knepley const PetscInt p = j * n[0] + i; 599a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) coords[p * dim + d] = x[d]; 600a4e36a32SMatthew G. Knepley } 601a4e36a32SMatthew G. Knepley } 602a4e36a32SMatthew G. Knepley break; 603a4e36a32SMatthew G. Knepley case 3: 604a4e36a32SMatthew G. Knepley x[0] = user->partLower[0]; 605a4e36a32SMatthew G. Knepley for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 606a4e36a32SMatthew G. Knepley x[1] = user->partLower[1]; 607a4e36a32SMatthew G. Knepley for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 608a4e36a32SMatthew G. Knepley x[2] = user->partLower[2]; 609a4e36a32SMatthew G. Knepley for (k = 0; k < n[2]; ++k, x[2] += dx[2]) { 610a4e36a32SMatthew G. Knepley const PetscInt p = (k * n[1] + j) * n[0] + i; 611a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) coords[p * dim + d] = x[d]; 612a4e36a32SMatthew G. Knepley } 613a4e36a32SMatthew G. Knepley } 614a4e36a32SMatthew G. Knepley } 615a4e36a32SMatthew G. Knepley break; 61663a3b9bcSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_SUP, "Do not support particle layout in dimension %" PetscInt_FMT, dim); 617a4e36a32SMatthew G. Knepley } 6189566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 6199566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 620a4e36a32SMatthew G. Knepley for (p = 0; p < Np; ++p) cellid[p] = 0; 6219566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 6229566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate(sdm, PETSC_TRUE)); 623a4e36a32SMatthew G. Knepley break; 62498921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]); 625a4e36a32SMatthew G. Knepley } 6269566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sdm, "Particles")); 6279566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(sdm, NULL, "-dm_view")); 628a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 629a4e36a32SMatthew G. Knepley } 630a4e36a32SMatthew G. Knepley 631*9371c9d4SSatish Balay static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) { 632a4e36a32SMatthew G. Knepley Vec vec; 633a4e36a32SMatthew G. Knepley PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 634a4e36a32SMatthew G. Knepley 635a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 63663a3b9bcSJacob Faibussowitsch PetscCheck(ofield == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %" PetscInt_FMT, ofield); 637a4e36a32SMatthew G. Knepley funcs[nfield] = constant; 6389566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &vec)); 6399566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec)); 6409566063dSJacob Faibussowitsch PetscCall(VecNormalize(vec, NULL)); 6419566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)vec, "Pressure Null Space")); 6429566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(vec, NULL, "-pressure_nullspace_view")); 6439566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullSpace)); 6449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec)); 645a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 646a4e36a32SMatthew G. Knepley } 647a4e36a32SMatthew G. Knepley 648*9371c9d4SSatish Balay static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u) { 649a4e36a32SMatthew G. Knepley DM dm; 650a4e36a32SMatthew G. Knepley MatNullSpace nullsp; 651a4e36a32SMatthew G. Knepley 6527510d9b0SBarry Smith PetscFunctionBeginUser; 6539566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 6549566063dSJacob Faibussowitsch PetscCall(CreatePressureNullSpace(dm, 1, 1, &nullsp)); 6559566063dSJacob Faibussowitsch PetscCall(MatNullSpaceRemove(nullsp, u)); 6569566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nullsp)); 657a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 658a4e36a32SMatthew G. Knepley } 659a4e36a32SMatthew G. Knepley 660a4e36a32SMatthew G. Knepley /* Make the discrete pressure discretely divergence free */ 661*9371c9d4SSatish Balay static PetscErrorCode RemoveDiscretePressureNullspace(TS ts) { 662a4e36a32SMatthew G. Knepley Vec u; 663a4e36a32SMatthew G. Knepley 6647510d9b0SBarry Smith PetscFunctionBeginUser; 6659566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &u)); 6669566063dSJacob Faibussowitsch PetscCall(RemoveDiscretePressureNullspace_Private(ts, u)); 667a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 668a4e36a32SMatthew G. Knepley } 669a4e36a32SMatthew G. Knepley 670*9371c9d4SSatish Balay static PetscErrorCode SetInitialConditions(TS ts, Vec u) { 671a4e36a32SMatthew G. Knepley DM dm; 672a4e36a32SMatthew G. Knepley PetscReal t; 673a4e36a32SMatthew G. Knepley 6747510d9b0SBarry Smith PetscFunctionBeginUser; 6759566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 6769566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 6779566063dSJacob Faibussowitsch PetscCall(DMComputeExactSolution(dm, t, u, NULL)); 6789566063dSJacob Faibussowitsch PetscCall(RemoveDiscretePressureNullspace_Private(ts, u)); 679a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 680a4e36a32SMatthew G. Knepley } 681a4e36a32SMatthew G. Knepley 682*9371c9d4SSatish Balay static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) { 683a4e36a32SMatthew G. Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 684a4e36a32SMatthew G. Knepley void *ctxs[3]; 685a4e36a32SMatthew G. Knepley DM dm; 686a4e36a32SMatthew G. Knepley PetscDS ds; 687a4e36a32SMatthew G. Knepley Vec v; 688a4e36a32SMatthew G. Knepley PetscReal ferrors[3]; 689b22b7e2eSMatthew G. Knepley PetscInt tl, l, f; 690a4e36a32SMatthew G. Knepley 691a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 6929566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 6939566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 694a4e36a32SMatthew G. Knepley 6959566063dSJacob Faibussowitsch for (f = 0; f < 3; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f])); 6969566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors)); 6979566063dSJacob Faibussowitsch PetscCall(PetscObjectGetTabLevel((PetscObject)ts, &tl)); 6989566063dSJacob Faibussowitsch for (l = 0; l < tl; ++l) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t")); 6999566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [%2.3g, %2.3g, %2.3g]\n", (int)step, (double)crtime, (double)ferrors[0], (double)ferrors[1], (double)ferrors[2])); 700a4e36a32SMatthew G. Knepley 7019566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &u)); 7029566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "Numerical Solution")); 7039566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 7049566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &u)); 705a4e36a32SMatthew G. Knepley 7069566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &v)); 7079566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, v)); 7089566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution")); 7099566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view")); 7109566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &v)); 711a4e36a32SMatthew G. Knepley 712a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 713a4e36a32SMatthew G. Knepley } 714a4e36a32SMatthew G. Knepley 715d2ae0b49SMatthew G. Knepley /* Note that adv->x0 will not be correct after migration */ 716*9371c9d4SSatish Balay static PetscErrorCode ComputeParticleError(TS ts, Vec u, Vec e) { 717a4e36a32SMatthew G. Knepley AdvCtx *adv; 718a4e36a32SMatthew G. Knepley DM sdm; 719a4e36a32SMatthew G. Knepley Parameter *param; 720a4e36a32SMatthew G. Knepley const PetscScalar *xp0, *xp; 721a4e36a32SMatthew G. Knepley PetscScalar *ep; 722a4e36a32SMatthew G. Knepley PetscReal time; 723a4e36a32SMatthew G. Knepley PetscInt dim, Np, p; 724a4e36a32SMatthew G. Knepley MPI_Comm comm; 725a4e36a32SMatthew G. Knepley 726a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 7279566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &time)); 7289566063dSJacob Faibussowitsch PetscCall(TSGetApplicationContext(ts, &adv)); 7299566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(adv->ctx->bag, (void **)¶m)); 7309566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 7319566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &sdm)); 7329566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sdm, &dim)); 7339566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(sdm, &Np)); 7349566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(adv->x0, &xp0)); 7359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(u, &xp)); 7369566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(e, &ep)); 737a4e36a32SMatthew G. Knepley for (p = 0; p < Np; ++p) { 738a4e36a32SMatthew G. Knepley PetscScalar x[3]; 739a4e36a32SMatthew G. Knepley PetscReal x0[3]; 740a4e36a32SMatthew G. Knepley PetscInt d; 741a4e36a32SMatthew G. Knepley 742a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p * dim + d]); 7439566063dSJacob Faibussowitsch PetscCall(adv->exact(dim, time, x0, 1, x, param)); 744a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) ep[p * dim + d] += x[d] - xp[p * dim + d]; 745a4e36a32SMatthew G. Knepley } 7469566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(adv->x0, &xp0)); 7479566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(u, &xp)); 7489566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(e, &ep)); 749a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 750a4e36a32SMatthew G. Knepley } 751a4e36a32SMatthew G. Knepley 752*9371c9d4SSatish Balay static PetscErrorCode MonitorParticleError(TS ts, PetscInt step, PetscReal time, Vec u, void *ctx) { 753a4e36a32SMatthew G. Knepley AdvCtx *adv = (AdvCtx *)ctx; 754a4e36a32SMatthew G. Knepley DM sdm; 755a4e36a32SMatthew G. Knepley Parameter *param; 756a4e36a32SMatthew G. Knepley const PetscScalar *xp0, *xp; 757a4e36a32SMatthew G. Knepley PetscReal error = 0.0; 758b22b7e2eSMatthew G. Knepley PetscInt dim, tl, l, Np, p; 759a4e36a32SMatthew G. Knepley MPI_Comm comm; 760a4e36a32SMatthew G. Knepley 761a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 7629566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(adv->ctx->bag, (void **)¶m)); 7639566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 7649566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &sdm)); 7659566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sdm, &dim)); 7669566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(sdm, &Np)); 7679566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(adv->x0, &xp0)); 7689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(u, &xp)); 769a4e36a32SMatthew G. Knepley for (p = 0; p < Np; ++p) { 770a4e36a32SMatthew G. Knepley PetscScalar x[3]; 771a4e36a32SMatthew G. Knepley PetscReal x0[3]; 772a4e36a32SMatthew G. Knepley PetscReal perror = 0.0; 773a4e36a32SMatthew G. Knepley PetscInt d; 774a4e36a32SMatthew G. Knepley 775a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p * dim + d]); 7769566063dSJacob Faibussowitsch PetscCall(adv->exact(dim, time, x0, 1, x, param)); 777a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) perror += PetscSqr(PetscRealPart(x[d] - xp[p * dim + d])); 778a4e36a32SMatthew G. Knepley error += perror; 779a4e36a32SMatthew G. Knepley } 7809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(adv->x0, &xp0)); 7819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(u, &xp)); 7829566063dSJacob Faibussowitsch PetscCall(PetscObjectGetTabLevel((PetscObject)ts, &tl)); 7839566063dSJacob Faibussowitsch for (l = 0; l < tl; ++l) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t")); 7849566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Timestep: %04d time = %-8.4g \t L_2 Particle Error: [%2.3g]\n", (int)step, (double)time, (double)error)); 785a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 786a4e36a32SMatthew G. Knepley } 787a4e36a32SMatthew G. Knepley 788*9371c9d4SSatish Balay static PetscErrorCode AdvectParticles(TS ts) { 789a4e36a32SMatthew G. Knepley TS sts; 790a4e36a32SMatthew G. Knepley DM sdm; 791ad95fc09SMatthew G. Knepley Vec coordinates; 792a4e36a32SMatthew G. Knepley AdvCtx *adv; 793a4e36a32SMatthew G. Knepley PetscReal time; 794d2ae0b49SMatthew G. Knepley PetscBool lreset, reset; 795d2ae0b49SMatthew G. Knepley PetscInt dim, n, N, newn, newN; 796a4e36a32SMatthew G. Knepley 797a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 7989566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)ts, "_SwarmTS", (PetscObject *)&sts)); 7999566063dSJacob Faibussowitsch PetscCall(TSGetDM(sts, &sdm)); 8009566063dSJacob Faibussowitsch PetscCall(TSGetRHSFunction(sts, NULL, NULL, (void **)&adv)); 8019566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sdm, &dim)); 8029566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(sdm, &N)); 8039566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(sdm, &n)); 8049566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sdm, DMSwarmPICField_coor, &coordinates)); 8059566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &time)); 8069566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(sts, time)); 807a4e36a32SMatthew G. Knepley adv->tf = time; 8089566063dSJacob Faibussowitsch PetscCall(TSSolve(sts, coordinates)); 8099566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sdm, DMSwarmPICField_coor, &coordinates)); 8109566063dSJacob Faibussowitsch PetscCall(VecCopy(adv->uf, adv->ui)); 811a4e36a32SMatthew G. Knepley adv->ti = adv->tf; 812a4e36a32SMatthew G. Knepley 8139566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate(sdm, PETSC_TRUE)); 8149566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(sdm, &newN)); 8159566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(sdm, &newn)); 816d2ae0b49SMatthew G. Knepley lreset = (n != newn || N != newN) ? PETSC_TRUE : PETSC_FALSE; 8179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&lreset, &reset, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sts))); 818ad95fc09SMatthew G. Knepley if (reset) { 8199566063dSJacob Faibussowitsch PetscCall(TSReset(sts)); 8209566063dSJacob Faibussowitsch PetscCall(DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor)); 821ad95fc09SMatthew G. Knepley } 8229566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(sdm, NULL, "-dm_view")); 823a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 824a4e36a32SMatthew G. Knepley } 825a4e36a32SMatthew G. Knepley 826*9371c9d4SSatish Balay int main(int argc, char **argv) { 827a4e36a32SMatthew G. Knepley DM dm, sdm; 828a4e36a32SMatthew G. Knepley TS ts, sts; 829d2ae0b49SMatthew G. Knepley Vec u, xtmp; 830a4e36a32SMatthew G. Knepley AppCtx user; 831a4e36a32SMatthew G. Knepley AdvCtx adv; 832a4e36a32SMatthew G. Knepley PetscReal t; 833a4e36a32SMatthew G. Knepley PetscInt dim; 834a4e36a32SMatthew G. Knepley 835327415f7SBarry Smith PetscFunctionBeginUser; 8369566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 8379566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 8389566063dSJacob Faibussowitsch PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag)); 8399566063dSJacob Faibussowitsch PetscCall(SetupParameters(&user)); 8409566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 8419566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 8429566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm)); 8439566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &user)); 844a4e36a32SMatthew G. Knepley /* Discretize chemical species */ 8459566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &sdm)); 8469566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)sdm, "part_")); 8479566063dSJacob Faibussowitsch PetscCall(DMSetType(sdm, DMSWARM)); 8489566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 8499566063dSJacob Faibussowitsch PetscCall(DMSetDimension(sdm, dim)); 8509566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(sdm, dm)); 851a4e36a32SMatthew G. Knepley /* Setup problem */ 8529566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, sdm, &user)); 8539566063dSJacob Faibussowitsch PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 854a4e36a32SMatthew G. Knepley 8559566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 8569566063dSJacob Faibussowitsch PetscCall(DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace)); 857a4e36a32SMatthew G. Knepley 8589566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 8599566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 8609566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 8619566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 8629566063dSJacob Faibussowitsch PetscCall(TSSetPreStep(ts, RemoveDiscretePressureNullspace)); 8639566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, MonitorError, &user, NULL)); 8649566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 865a4e36a32SMatthew G. Knepley 8669566063dSJacob Faibussowitsch PetscCall(TSSetComputeInitialCondition(ts, SetInitialConditions)); /* Must come after SetFromOptions() */ 8679566063dSJacob Faibussowitsch PetscCall(SetInitialConditions(ts, u)); 8689566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 8699566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dm, 0, t)); 8709566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u)); 871a4e36a32SMatthew G. Knepley 872a4e36a32SMatthew G. Knepley /* Setup particle position integrator */ 8739566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &sts)); 8749566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)sts, "part_")); 8759566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)sts, (PetscObject)ts, 1)); 8769566063dSJacob Faibussowitsch PetscCall(TSSetDM(sts, sdm)); 8779566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(sts, TS_NONLINEAR)); 8789566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(sts, TS_EXACTFINALTIME_MATCHSTEP)); 8799566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(sts, MonitorParticleError, &adv, NULL)); 8809566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(sts)); 8819566063dSJacob Faibussowitsch PetscCall(TSSetApplicationContext(sts, &adv)); 8829566063dSJacob Faibussowitsch PetscCall(TSSetComputeExactError(sts, ComputeParticleError)); 8839566063dSJacob Faibussowitsch PetscCall(TSSetComputeInitialCondition(sts, SetInitialParticleConditions)); 884a4e36a32SMatthew G. Knepley adv.ti = t; 885a4e36a32SMatthew G. Knepley adv.uf = u; 8869566063dSJacob Faibussowitsch PetscCall(VecDuplicate(adv.uf, &adv.ui)); 8879566063dSJacob Faibussowitsch PetscCall(VecCopy(u, adv.ui)); 8889566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(sts, NULL, FreeStreaming, &adv)); 8899566063dSJacob Faibussowitsch PetscCall(TSSetPostStep(ts, AdvectParticles)); 8909566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)ts, "_SwarmTS", (PetscObject)sts)); 8919566063dSJacob Faibussowitsch PetscCall(DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor)); 8929566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(sdm, &adv.x0)); 8939566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp)); 8949566063dSJacob Faibussowitsch PetscCall(VecCopy(xtmp, adv.x0)); 8959566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp)); 896a4e36a32SMatthew G. Knepley switch (user.solType) { 897a4e36a32SMatthew G. Knepley case SOL_TRIG_TRIG: adv.exact = trig_trig_x; break; 89863a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)sdm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", solTypes[PetscMin(user.solType, NUM_SOL_TYPES)], user.solType); 899a4e36a32SMatthew G. Knepley } 900a4e36a32SMatthew G. Knepley adv.ctx = &user; 901a4e36a32SMatthew G. Knepley 9029566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 9039566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u)); 9049566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "Numerical Solution")); 905a4e36a32SMatthew G. Knepley 9069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 9079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&adv.x0)); 9089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&adv.ui)); 9099566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 9109566063dSJacob Faibussowitsch PetscCall(DMDestroy(&sdm)); 9119566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 9129566063dSJacob Faibussowitsch PetscCall(TSDestroy(&sts)); 9139566063dSJacob Faibussowitsch PetscCall(PetscBagDestroy(&user.bag)); 9149566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 915b122ec5aSJacob Faibussowitsch return 0; 916a4e36a32SMatthew G. Knepley } 917a4e36a32SMatthew G. Knepley 918a4e36a32SMatthew G. Knepley /*TEST 919a4e36a32SMatthew G. Knepley 920a4e36a32SMatthew G. Knepley # Swarm does not work with complex 921a4e36a32SMatthew G. Knepley test: 922a4e36a32SMatthew G. Knepley suffix: 2d_tri_p2_p1_p1_tconvp 923a4e36a32SMatthew G. Knepley requires: triangle !single !complex 924a4e36a32SMatthew G. Knepley args: -dm_plex_separate_marker -sol_type trig_trig -dm_refine 2 \ 925a4e36a32SMatthew G. Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 926a4e36a32SMatthew G. Knepley -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 -ts_monitor_cancel \ 927a4e36a32SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 928a4e36a32SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 929a4e36a32SMatthew G. Knepley -fieldsplit_0_pc_type lu \ 930a4e36a32SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \ 931a4e36a32SMatthew G. Knepley -omega 0.5 -part_layout_type box -part_lower 0.25,0.25 -part_upper 0.75,0.75 -Npb 5 \ 932a4e36a32SMatthew G. Knepley -part_ts_max_steps 2 -part_ts_dt 0.05 -part_ts_convergence_estimate -convest_num_refine 1 -part_ts_monitor_cancel 933b22b7e2eSMatthew G. Knepley test: 934b22b7e2eSMatthew G. Knepley suffix: 2d_tri_p2_p1_p1_exit 935b22b7e2eSMatthew G. Knepley requires: triangle !single !complex 936b22b7e2eSMatthew G. Knepley args: -dm_plex_separate_marker -sol_type trig_trig -dm_refine 1 \ 937b22b7e2eSMatthew G. Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 938b22b7e2eSMatthew G. Knepley -dmts_check .001 -ts_max_steps 10 -ts_dt 0.1 \ 939b22b7e2eSMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 940b22b7e2eSMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 941b22b7e2eSMatthew G. Knepley -fieldsplit_0_pc_type lu \ 942b22b7e2eSMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \ 943b22b7e2eSMatthew G. Knepley -omega 0.5 -part_layout_type box -part_lower 0.25,0.25 -part_upper 0.75,0.75 -Npb 5 \ 944b22b7e2eSMatthew G. Knepley -part_ts_max_steps 20 -part_ts_dt 0.05 945a4e36a32SMatthew G. Knepley 946a4e36a32SMatthew G. Knepley TEST*/ 947