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 34a4e36a32SMatthew G. Knepley typedef enum {SOL_TRIG_TRIG, NUM_SOL_TYPES} SolType; 35a4e36a32SMatthew G. Knepley const char *solTypes[NUM_SOL_TYPES+1] = {"trig_trig", "unknown"}; 36a4e36a32SMatthew G. Knepley 37a4e36a32SMatthew G. Knepley typedef enum {PART_LAYOUT_CELL, PART_LAYOUT_BOX, NUM_PART_LAYOUT_TYPES} PartLayoutType; 38a4e36a32SMatthew G. Knepley const char *partLayoutTypes[NUM_PART_LAYOUT_TYPES+1] = {"cell", "box", "unknown"}; 39a4e36a32SMatthew G. Knepley 40a4e36a32SMatthew G. Knepley typedef struct { 41a4e36a32SMatthew G. Knepley PetscReal nu; /* Kinematic viscosity */ 42a4e36a32SMatthew G. Knepley PetscReal alpha; /* Thermal diffusivity */ 43a4e36a32SMatthew G. Knepley PetscReal T_in; /* Inlet temperature*/ 44a4e36a32SMatthew G. Knepley PetscReal omega; /* Rotation speed in MMS benchmark */ 45a4e36a32SMatthew G. Knepley } Parameter; 46a4e36a32SMatthew G. Knepley 47a4e36a32SMatthew G. Knepley typedef struct { 48a4e36a32SMatthew G. Knepley /* Problem definition */ 49a4e36a32SMatthew G. Knepley PetscBag bag; /* Holds problem parameters */ 50a4e36a32SMatthew G. Knepley SolType solType; /* MMS solution type */ 51a4e36a32SMatthew G. Knepley PartLayoutType partLayout; /* Type of particle distribution */ 52a4e36a32SMatthew G. Knepley PetscInt Npc; /* The initial number of particles per cell */ 53a4e36a32SMatthew G. Knepley PetscReal partLower[3]; /* Lower left corner of particle box */ 54a4e36a32SMatthew G. Knepley PetscReal partUpper[3]; /* Upper right corner of particle box */ 55a4e36a32SMatthew G. Knepley PetscInt Npb; /* The initial number of particles per box dimension */ 56a4e36a32SMatthew G. Knepley } AppCtx; 57a4e36a32SMatthew G. Knepley 58a4e36a32SMatthew G. Knepley typedef struct { 59a4e36a32SMatthew G. Knepley PetscReal ti; /* The time for ui, at the beginning of the advection solve */ 60a4e36a32SMatthew G. Knepley PetscReal tf; /* The time for uf, at the end of the advection solve */ 61a4e36a32SMatthew G. Knepley Vec ui; /* The PDE solution field at ti */ 62a4e36a32SMatthew G. Knepley Vec uf; /* The PDE solution field at tf */ 63a4e36a32SMatthew G. Knepley Vec x0; /* The initial particle positions at t = 0 */ 64a4e36a32SMatthew G. Knepley PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 65a4e36a32SMatthew G. Knepley AppCtx *ctx; /* Context for exact solution */ 66a4e36a32SMatthew G. Knepley } AdvCtx; 67a4e36a32SMatthew G. Knepley 68a4e36a32SMatthew G. Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 69a4e36a32SMatthew G. Knepley { 70a4e36a32SMatthew G. Knepley PetscInt d; 71a4e36a32SMatthew G. Knepley for (d = 0; d < Nc; ++d) u[d] = 0.0; 72a4e36a32SMatthew G. Knepley return 0; 73a4e36a32SMatthew G. Knepley } 74a4e36a32SMatthew G. Knepley 75a4e36a32SMatthew G. Knepley static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 76a4e36a32SMatthew G. Knepley { 77a4e36a32SMatthew G. Knepley PetscInt d; 78a4e36a32SMatthew G. Knepley for (d = 0; d < Nc; ++d) u[d] = 1.0; 79a4e36a32SMatthew G. Knepley return 0; 80a4e36a32SMatthew G. Knepley } 81a4e36a32SMatthew G. Knepley 82a4e36a32SMatthew G. Knepley /* 83a4e36a32SMatthew G. Knepley CASE: trigonometric-trigonometric 84a4e36a32SMatthew G. Knepley In 2D we use exact solution: 85a4e36a32SMatthew G. Knepley 86a4e36a32SMatthew G. Knepley x = r0 cos(w t + theta0) r0 = sqrt(x0^2 + y0^2) 87a4e36a32SMatthew G. Knepley y = r0 sin(w t + theta0) theta0 = arctan(y0/x0) 88a4e36a32SMatthew G. Knepley u = -w r0 sin(theta0) = -w y 89a4e36a32SMatthew G. Knepley v = w r0 cos(theta0) = w x 90a4e36a32SMatthew G. Knepley p = x + y - 1 91a4e36a32SMatthew G. Knepley T = t + x + y 92a4e36a32SMatthew G. Knepley f = <1, 1> 93a4e36a32SMatthew G. Knepley Q = 1 + w (x - y)/r 94a4e36a32SMatthew G. Knepley 95a4e36a32SMatthew G. Knepley so that 96a4e36a32SMatthew G. Knepley 97a4e36a32SMatthew G. Knepley \nabla \cdot u = 0 + 0 = 0 98a4e36a32SMatthew G. Knepley 99a4e36a32SMatthew G. Knepley f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p 100a4e36a32SMatthew G. Knepley = <0, 0> + u_i d_i u_j - \nu 0 + <1, 1> 101a4e36a32SMatthew G. Knepley = <1, 1> + w^2 <-y, x> . <<0, 1>, <-1, 0>> 102a4e36a32SMatthew G. Knepley = <1, 1> + w^2 <-x, -y> 103a4e36a32SMatthew G. Knepley = <1, 1> - w^2 <x, y> 104a4e36a32SMatthew G. Knepley 105a4e36a32SMatthew G. Knepley Q = dT/dt + u \cdot \nabla T - \alpha \Delta T 106a4e36a32SMatthew G. Knepley = 1 + <u, v> . <1, 1> - \alpha 0 107a4e36a32SMatthew G. Knepley = 1 + u + v 108a4e36a32SMatthew G. Knepley */ 109a4e36a32SMatthew G. Knepley static PetscErrorCode trig_trig_x(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *x, void *ctx) 110a4e36a32SMatthew G. Knepley { 111a4e36a32SMatthew G. Knepley const PetscReal x0 = X[0]; 112a4e36a32SMatthew G. Knepley const PetscReal y0 = X[1]; 113a4e36a32SMatthew G. Knepley const PetscReal R0 = PetscSqrtReal(x0*x0 + y0*y0); 114a4e36a32SMatthew G. Knepley const PetscReal theta0 = PetscAtan2Real(y0, x0); 115a4e36a32SMatthew G. Knepley Parameter *p = (Parameter *) ctx; 116a4e36a32SMatthew G. Knepley 117a4e36a32SMatthew G. Knepley x[0] = R0*PetscCosReal(p->omega*time + theta0); 118a4e36a32SMatthew G. Knepley x[1] = R0*PetscSinReal(p->omega*time + theta0); 119a4e36a32SMatthew G. Knepley return 0; 120a4e36a32SMatthew G. Knepley } 121a4e36a32SMatthew G. Knepley static PetscErrorCode trig_trig_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 122a4e36a32SMatthew G. Knepley { 123a4e36a32SMatthew G. Knepley Parameter *p = (Parameter *) ctx; 124a4e36a32SMatthew G. Knepley 125a4e36a32SMatthew G. Knepley u[0] = -p->omega*X[1]; 126a4e36a32SMatthew G. Knepley u[1] = p->omega*X[0]; 127a4e36a32SMatthew G. Knepley return 0; 128a4e36a32SMatthew G. Knepley } 129a4e36a32SMatthew G. Knepley static PetscErrorCode trig_trig_u_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 130a4e36a32SMatthew G. Knepley { 131a4e36a32SMatthew G. Knepley u[0] = 0.0; 132a4e36a32SMatthew G. Knepley u[1] = 0.0; 133a4e36a32SMatthew G. Knepley return 0; 134a4e36a32SMatthew G. Knepley } 135a4e36a32SMatthew G. Knepley 136a4e36a32SMatthew G. Knepley static PetscErrorCode trig_trig_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 137a4e36a32SMatthew G. Knepley { 138a4e36a32SMatthew G. Knepley p[0] = X[0] + X[1] - 1.0; 139a4e36a32SMatthew G. Knepley return 0; 140a4e36a32SMatthew G. Knepley } 141a4e36a32SMatthew G. Knepley 142a4e36a32SMatthew G. Knepley static PetscErrorCode trig_trig_T(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 143a4e36a32SMatthew G. Knepley { 144a4e36a32SMatthew G. Knepley T[0] = time + X[0] + X[1]; 145a4e36a32SMatthew G. Knepley return 0; 146a4e36a32SMatthew G. Knepley } 147a4e36a32SMatthew G. Knepley static PetscErrorCode trig_trig_T_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 148a4e36a32SMatthew G. Knepley { 149a4e36a32SMatthew G. Knepley T[0] = 1.0; 150a4e36a32SMatthew G. Knepley return 0; 151a4e36a32SMatthew G. Knepley } 152a4e36a32SMatthew G. Knepley 153a4e36a32SMatthew G. Knepley static void f0_trig_trig_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 154a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 155a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 156a4e36a32SMatthew G. Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 157a4e36a32SMatthew G. Knepley { 158a4e36a32SMatthew G. Knepley const PetscReal omega = PetscRealPart(constants[3]); 159a4e36a32SMatthew G. Knepley PetscInt Nc = dim; 160a4e36a32SMatthew G. Knepley PetscInt c, d; 161a4e36a32SMatthew G. Knepley 162a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[d] = u_t[uOff[0]+d]; 163a4e36a32SMatthew G. Knepley 164a4e36a32SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 165a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 166a4e36a32SMatthew G. Knepley } 167a4e36a32SMatthew G. Knepley f0[0] -= 1.0 - omega*omega*X[0]; 168a4e36a32SMatthew G. Knepley f0[1] -= 1.0 - omega*omega*X[1]; 169a4e36a32SMatthew G. Knepley } 170a4e36a32SMatthew G. Knepley 171a4e36a32SMatthew G. Knepley static void f0_trig_trig_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 172a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 173a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 174a4e36a32SMatthew G. Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 175a4e36a32SMatthew G. Knepley { 176a4e36a32SMatthew G. Knepley const PetscReal omega = PetscRealPart(constants[3]); 177a4e36a32SMatthew G. Knepley PetscInt d; 178a4e36a32SMatthew G. Knepley 179a4e36a32SMatthew G. Knepley for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 180a4e36a32SMatthew G. Knepley f0[0] += u_t[uOff[2]] - (1.0 + omega*(X[0] - X[1])); 181a4e36a32SMatthew G. Knepley } 182a4e36a32SMatthew G. Knepley 183a4e36a32SMatthew G. Knepley static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 184a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 185a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 186a4e36a32SMatthew G. Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 187a4e36a32SMatthew G. Knepley { 188a4e36a32SMatthew G. Knepley PetscInt d; 189a4e36a32SMatthew G. Knepley for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; 190a4e36a32SMatthew G. Knepley } 191a4e36a32SMatthew G. Knepley 192a4e36a32SMatthew G. Knepley /*f1_v = \nu[grad(u) + grad(u)^T] - pI */ 193a4e36a32SMatthew G. Knepley static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 194a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 195a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 196a4e36a32SMatthew G. Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 197a4e36a32SMatthew G. Knepley { 198a4e36a32SMatthew G. Knepley const PetscReal nu = PetscRealPart(constants[0]); 199a4e36a32SMatthew G. Knepley const PetscInt Nc = dim; 200a4e36a32SMatthew G. Knepley PetscInt c, d; 201a4e36a32SMatthew G. Knepley 202a4e36a32SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 203a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) { 204a4e36a32SMatthew G. Knepley f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]); 205a4e36a32SMatthew G. Knepley } 206a4e36a32SMatthew G. Knepley f1[c*dim+c] -= u[uOff[1]]; 207a4e36a32SMatthew G. Knepley } 208a4e36a32SMatthew G. Knepley } 209a4e36a32SMatthew G. Knepley 210a4e36a32SMatthew G. Knepley static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 211a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 212a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 213a4e36a32SMatthew G. Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 214a4e36a32SMatthew G. Knepley { 215a4e36a32SMatthew G. Knepley const PetscReal alpha = PetscRealPart(constants[1]); 216a4e36a32SMatthew G. Knepley PetscInt d; 217a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d]; 218a4e36a32SMatthew G. Knepley } 219a4e36a32SMatthew G. Knepley 220a4e36a32SMatthew G. Knepley /*Jacobians*/ 221a4e36a32SMatthew G. Knepley static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 222a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 223a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 224a4e36a32SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 225a4e36a32SMatthew G. Knepley { 226a4e36a32SMatthew G. Knepley PetscInt d; 227a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; 228a4e36a32SMatthew G. Knepley } 229a4e36a32SMatthew G. Knepley 230a4e36a32SMatthew G. Knepley static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 231a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 232a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 233a4e36a32SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 234a4e36a32SMatthew G. Knepley { 235a4e36a32SMatthew G. Knepley PetscInt c, d; 236a4e36a32SMatthew G. Knepley const PetscInt Nc = dim; 237a4e36a32SMatthew G. Knepley 238a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g0[d*dim+d] = u_tShift; 239a4e36a32SMatthew G. Knepley 240a4e36a32SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 241a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) { 242a4e36a32SMatthew G. Knepley g0[c*Nc+d] += u_x[c*Nc+d]; 243a4e36a32SMatthew G. Knepley } 244a4e36a32SMatthew G. Knepley } 245a4e36a32SMatthew G. Knepley } 246a4e36a32SMatthew G. Knepley 247a4e36a32SMatthew G. Knepley static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 248a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 249a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 250a4e36a32SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 251a4e36a32SMatthew G. Knepley { 252a4e36a32SMatthew G. Knepley PetscInt NcI = dim; 253a4e36a32SMatthew G. Knepley PetscInt NcJ = dim; 254a4e36a32SMatthew G. Knepley PetscInt c, d, e; 255a4e36a32SMatthew G. Knepley 256a4e36a32SMatthew G. Knepley for (c = 0; c < NcI; ++c) { 257a4e36a32SMatthew G. Knepley for (d = 0; d < NcJ; ++d) { 258a4e36a32SMatthew G. Knepley for (e = 0; e < dim; ++e) { 259a4e36a32SMatthew G. Knepley if (c == d) { 260a4e36a32SMatthew G. Knepley g1[(c*NcJ+d)*dim+e] += u[e]; 261a4e36a32SMatthew G. Knepley } 262a4e36a32SMatthew G. Knepley } 263a4e36a32SMatthew G. Knepley } 264a4e36a32SMatthew G. Knepley } 265a4e36a32SMatthew G. Knepley } 266a4e36a32SMatthew G. Knepley 267a4e36a32SMatthew G. Knepley static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 268a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 269a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 270a4e36a32SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 271a4e36a32SMatthew G. Knepley { 272a4e36a32SMatthew G. Knepley PetscInt d; 273a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; 274a4e36a32SMatthew G. Knepley } 275a4e36a32SMatthew G. Knepley 276a4e36a32SMatthew G. Knepley static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 277a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 278a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 279a4e36a32SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 280a4e36a32SMatthew G. Knepley { 281a4e36a32SMatthew G. Knepley const PetscReal nu = PetscRealPart(constants[0]); 282a4e36a32SMatthew G. Knepley const PetscInt Nc = dim; 283a4e36a32SMatthew G. Knepley PetscInt c, d; 284a4e36a32SMatthew G. Knepley 285a4e36a32SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 286a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) { 287a4e36a32SMatthew G. Knepley g3[((c*Nc+c)*dim+d)*dim+d] += nu; // gradU 288a4e36a32SMatthew G. Knepley g3[((c*Nc+d)*dim+d)*dim+c] += nu; // gradU transpose 289a4e36a32SMatthew G. Knepley } 290a4e36a32SMatthew G. Knepley } 291a4e36a32SMatthew G. Knepley } 292a4e36a32SMatthew G. Knepley 293a4e36a32SMatthew G. Knepley static void g0_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 294a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 295a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 296a4e36a32SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 297a4e36a32SMatthew G. Knepley { 298a4e36a32SMatthew G. Knepley PetscInt d; 299a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g0[d] = u_tShift; 300a4e36a32SMatthew G. Knepley } 301a4e36a32SMatthew G. Knepley 302a4e36a32SMatthew G. Knepley static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 303a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 304a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 305a4e36a32SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 306a4e36a32SMatthew G. Knepley { 307a4e36a32SMatthew G. Knepley PetscInt d; 308a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d]; 309a4e36a32SMatthew G. Knepley } 310a4e36a32SMatthew G. Knepley 311a4e36a32SMatthew G. Knepley static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 312a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 313a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 314a4e36a32SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 315a4e36a32SMatthew G. Knepley { 316a4e36a32SMatthew G. Knepley PetscInt d; 317a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d]; 318a4e36a32SMatthew G. Knepley } 319a4e36a32SMatthew G. Knepley 320a4e36a32SMatthew G. Knepley static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 321a4e36a32SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 322a4e36a32SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 323a4e36a32SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 324a4e36a32SMatthew G. Knepley { 325a4e36a32SMatthew G. Knepley const PetscReal alpha = PetscRealPart(constants[1]); 326a4e36a32SMatthew G. Knepley PetscInt d; 327a4e36a32SMatthew G. Knepley 328a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha; 329a4e36a32SMatthew G. Knepley } 330a4e36a32SMatthew G. Knepley 331a4e36a32SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 332a4e36a32SMatthew G. Knepley { 333a4e36a32SMatthew G. Knepley PetscInt sol, pl, n; 334a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 335a4e36a32SMatthew G. Knepley 336a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 337a4e36a32SMatthew G. Knepley options->solType = SOL_TRIG_TRIG; 338a4e36a32SMatthew G. Knepley options->partLayout = PART_LAYOUT_CELL; 339a4e36a32SMatthew G. Knepley options->Npc = 1; 340a4e36a32SMatthew G. Knepley options->Npb = 1; 341a4e36a32SMatthew G. Knepley 342a4e36a32SMatthew G. Knepley options->partLower[0] = options->partLower[1] = options->partLower[2] = 0.; 343a4e36a32SMatthew G. Knepley options->partUpper[0] = options->partUpper[1] = options->partUpper[2] = 1.; 344a4e36a32SMatthew G. Knepley ierr = PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX");CHKERRQ(ierr); 345a4e36a32SMatthew G. Knepley sol = options->solType; 346a4e36a32SMatthew G. Knepley ierr = PetscOptionsEList("-sol_type", "The solution type", "ex77.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr); 347a4e36a32SMatthew G. Knepley options->solType = (SolType) sol; 348a4e36a32SMatthew G. Knepley pl = options->partLayout; 349a4e36a32SMatthew G. Knepley ierr = PetscOptionsEList("-part_layout_type", "The particle layout type", "ex77.c", partLayoutTypes, NUM_PART_LAYOUT_TYPES, partLayoutTypes[options->partLayout], &pl, NULL);CHKERRQ(ierr); 350a4e36a32SMatthew G. Knepley options->partLayout = (PartLayoutType) pl; 351a4e36a32SMatthew G. Knepley ierr = PetscOptionsInt("-Npc", "The initial number of particles per cell", "ex77.c", options->Npc, &options->Npc, NULL);CHKERRQ(ierr); 352a4e36a32SMatthew G. Knepley n = 3; 353a4e36a32SMatthew G. Knepley ierr = PetscOptionsRealArray("-part_lower", "The lower left corner of the particle box", "ex77.c", options->partLower, &n, NULL);CHKERRQ(ierr); 354a4e36a32SMatthew G. Knepley n = 3; 355a4e36a32SMatthew G. Knepley ierr = PetscOptionsRealArray("-part_upper", "The upper right corner of the particle box", "ex77.c", options->partUpper, &n, NULL);CHKERRQ(ierr); 356a4e36a32SMatthew G. Knepley ierr = PetscOptionsInt("-Npb", "The initial number of particles per box dimension", "ex77.c", options->Npb, &options->Npb, NULL);CHKERRQ(ierr); 357*1e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 358a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 359a4e36a32SMatthew G. Knepley } 360a4e36a32SMatthew G. Knepley 361a4e36a32SMatthew G. Knepley static PetscErrorCode SetupParameters(AppCtx *user) 362a4e36a32SMatthew G. Knepley { 363a4e36a32SMatthew G. Knepley PetscBag bag; 364a4e36a32SMatthew G. Knepley Parameter *p; 365a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 366a4e36a32SMatthew G. Knepley 367a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 368a4e36a32SMatthew G. Knepley /* setup PETSc parameter bag */ 369a4e36a32SMatthew G. Knepley ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr); 370a4e36a32SMatthew G. Knepley ierr = PetscBagSetName(user->bag, "par", "Low Mach flow parameters");CHKERRQ(ierr); 371a4e36a32SMatthew G. Knepley bag = user->bag; 372a4e36a32SMatthew G. Knepley ierr = PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity");CHKERRQ(ierr); 373a4e36a32SMatthew G. Knepley ierr = PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity");CHKERRQ(ierr); 374a4e36a32SMatthew G. Knepley ierr = PetscBagRegisterReal(bag, &p->T_in, 1.0, "T_in", "Inlet temperature");CHKERRQ(ierr); 375a4e36a32SMatthew G. Knepley ierr = PetscBagRegisterReal(bag, &p->omega, 1.0, "omega", "Rotation speed in MMS benchmark");CHKERRQ(ierr); 376a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 377a4e36a32SMatthew G. Knepley } 378a4e36a32SMatthew G. Knepley 379a4e36a32SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 380a4e36a32SMatthew G. Knepley { 381a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 382a4e36a32SMatthew G. Knepley 383a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 38430602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 38530602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 386a4e36a32SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 387a4e36a32SMatthew G. Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 388a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 389a4e36a32SMatthew G. Knepley } 390a4e36a32SMatthew G. Knepley 391a4e36a32SMatthew G. Knepley static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 392a4e36a32SMatthew G. Knepley { 393a4e36a32SMatthew G. Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 394a4e36a32SMatthew G. Knepley PetscErrorCode (*exactFuncs_t[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 395a4e36a32SMatthew G. Knepley PetscDS prob; 39645480ffeSMatthew G. Knepley DMLabel label; 397a4e36a32SMatthew G. Knepley Parameter *ctx; 398a4e36a32SMatthew G. Knepley PetscInt id; 399a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 400a4e36a32SMatthew G. Knepley 401a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 40245480ffeSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 403a4e36a32SMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 404a4e36a32SMatthew G. Knepley switch(user->solType) { 405a4e36a32SMatthew G. Knepley case SOL_TRIG_TRIG: 406a4e36a32SMatthew G. Knepley ierr = PetscDSSetResidual(prob, 0, f0_trig_trig_v, f1_v);CHKERRQ(ierr); 407a4e36a32SMatthew G. Knepley ierr = PetscDSSetResidual(prob, 2, f0_trig_trig_w, f1_w);CHKERRQ(ierr); 408a4e36a32SMatthew G. Knepley 409a4e36a32SMatthew G. Knepley exactFuncs[0] = trig_trig_u; 410a4e36a32SMatthew G. Knepley exactFuncs[1] = trig_trig_p; 411a4e36a32SMatthew G. Knepley exactFuncs[2] = trig_trig_T; 412a4e36a32SMatthew G. Knepley exactFuncs_t[0] = trig_trig_u_t; 413a4e36a32SMatthew G. Knepley exactFuncs_t[1] = NULL; 414a4e36a32SMatthew G. Knepley exactFuncs_t[2] = trig_trig_T_t; 415a4e36a32SMatthew G. Knepley break; 416a4e36a32SMatthew G. Knepley default: SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 417a4e36a32SMatthew G. Knepley } 418a4e36a32SMatthew G. Knepley 419a4e36a32SMatthew G. Knepley ierr = PetscDSSetResidual(prob, 1, f0_q, NULL);CHKERRQ(ierr); 420a4e36a32SMatthew G. Knepley 421a4e36a32SMatthew G. Knepley ierr = PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu);CHKERRQ(ierr); 422a4e36a32SMatthew G. Knepley ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL);CHKERRQ(ierr); 423a4e36a32SMatthew G. Knepley ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL);CHKERRQ(ierr); 424a4e36a32SMatthew G. Knepley ierr = PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL);CHKERRQ(ierr); 425a4e36a32SMatthew G. Knepley ierr = PetscDSSetJacobian(prob, 2, 2, g0_wT, g1_wT, NULL, g3_wT);CHKERRQ(ierr); 426a4e36a32SMatthew G. Knepley /* Setup constants */ 427a4e36a32SMatthew G. Knepley { 428a4e36a32SMatthew G. Knepley Parameter *param; 429a4e36a32SMatthew G. Knepley PetscScalar constants[4]; 430a4e36a32SMatthew G. Knepley 431a4e36a32SMatthew G. Knepley ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 432a4e36a32SMatthew G. Knepley 433a4e36a32SMatthew G. Knepley constants[0] = param->nu; 434a4e36a32SMatthew G. Knepley constants[1] = param->alpha; 435a4e36a32SMatthew G. Knepley constants[2] = param->T_in; 436a4e36a32SMatthew G. Knepley constants[3] = param->omega; 437a4e36a32SMatthew G. Knepley ierr = PetscDSSetConstants(prob, 4, constants);CHKERRQ(ierr); 438a4e36a32SMatthew G. Knepley } 439a4e36a32SMatthew G. Knepley /* Setup Boundary Conditions */ 440a4e36a32SMatthew G. Knepley ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr); 441a4e36a32SMatthew G. Knepley id = 3; 44245480ffeSMatthew G. Knepley ierr = 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);CHKERRQ(ierr); 443a4e36a32SMatthew G. Knepley id = 1; 44445480ffeSMatthew G. Knepley ierr = 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);CHKERRQ(ierr); 445a4e36a32SMatthew G. Knepley id = 2; 44645480ffeSMatthew G. Knepley ierr = 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);CHKERRQ(ierr); 447a4e36a32SMatthew G. Knepley id = 4; 44845480ffeSMatthew G. Knepley ierr = 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);CHKERRQ(ierr); 449a4e36a32SMatthew G. Knepley id = 3; 45045480ffeSMatthew G. Knepley ierr = 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);CHKERRQ(ierr); 451a4e36a32SMatthew G. Knepley id = 1; 45245480ffeSMatthew G. Knepley ierr = 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);CHKERRQ(ierr); 453a4e36a32SMatthew G. Knepley id = 2; 45445480ffeSMatthew G. Knepley ierr = 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);CHKERRQ(ierr); 455a4e36a32SMatthew G. Knepley id = 4; 45645480ffeSMatthew G. Knepley ierr = 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);CHKERRQ(ierr); 457a4e36a32SMatthew G. Knepley 458a4e36a32SMatthew G. Knepley /*setup exact solution.*/ 459a4e36a32SMatthew G. Knepley ierr = PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx);CHKERRQ(ierr); 460a4e36a32SMatthew G. Knepley ierr = PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx);CHKERRQ(ierr); 461a4e36a32SMatthew G. Knepley ierr = PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx);CHKERRQ(ierr); 462a4e36a32SMatthew G. Knepley ierr = PetscDSSetExactSolutionTimeDerivative(prob, 0, exactFuncs_t[0], ctx);CHKERRQ(ierr); 463a4e36a32SMatthew G. Knepley ierr = PetscDSSetExactSolutionTimeDerivative(prob, 1, exactFuncs_t[1], ctx);CHKERRQ(ierr); 464a4e36a32SMatthew G. Knepley ierr = PetscDSSetExactSolutionTimeDerivative(prob, 2, exactFuncs_t[2], ctx);CHKERRQ(ierr); 465a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 466a4e36a32SMatthew G. Knepley } 467a4e36a32SMatthew G. Knepley 468a4e36a32SMatthew G. Knepley /* x_t = v 469a4e36a32SMatthew G. Knepley 470a4e36a32SMatthew G. Knepley Note that here we use the velocity field at t_{n+1} to advect the particles from 471a4e36a32SMatthew G. Knepley t_n to t_{n+1}. If we use both of these fields, we could use Crank-Nicholson or 472a4e36a32SMatthew G. Knepley the method of characteristics. 473a4e36a32SMatthew G. Knepley */ 474a4e36a32SMatthew G. Knepley static PetscErrorCode FreeStreaming(TS ts, PetscReal t, Vec X, Vec F, void *ctx) 475a4e36a32SMatthew G. Knepley { 476a4e36a32SMatthew G. Knepley AdvCtx *adv = (AdvCtx *) ctx; 477a4e36a32SMatthew G. Knepley Vec u = adv->ui; 478a4e36a32SMatthew G. Knepley DM sdm, dm, vdm; 479a4e36a32SMatthew G. Knepley Vec vel, locvel, pvel; 480a4e36a32SMatthew G. Knepley IS vis; 481a4e36a32SMatthew G. Knepley DMInterpolationInfo ictx; 482a4e36a32SMatthew G. Knepley const PetscScalar *coords, *v; 483a4e36a32SMatthew G. Knepley PetscScalar *f; 484a4e36a32SMatthew G. Knepley PetscInt vf[1] = {0}; 485a4e36a32SMatthew G. Knepley PetscInt dim, Np; 486a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 487a4e36a32SMatthew G. Knepley 488a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 489a4e36a32SMatthew G. Knepley ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr); 490a4e36a32SMatthew G. Knepley ierr = DMSwarmGetCellDM(sdm, &dm);CHKERRQ(ierr); 491a4e36a32SMatthew G. Knepley ierr = DMGetGlobalVector(sdm, &pvel);CHKERRQ(ierr); 492a4e36a32SMatthew G. Knepley ierr = DMSwarmGetLocalSize(sdm, &Np);CHKERRQ(ierr); 493a4e36a32SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 494a4e36a32SMatthew G. Knepley /* Get local velocity */ 495a4e36a32SMatthew G. Knepley ierr = DMCreateSubDM(dm, 1, vf, &vis, &vdm);CHKERRQ(ierr); 496a4e36a32SMatthew G. Knepley ierr = VecGetSubVector(u, vis, &vel);CHKERRQ(ierr); 497a4e36a32SMatthew G. Knepley ierr = DMGetLocalVector(vdm, &locvel);CHKERRQ(ierr); 498a4e36a32SMatthew G. Knepley ierr = DMPlexInsertBoundaryValues(vdm, PETSC_TRUE, locvel, adv->ti, NULL, NULL, NULL);CHKERRQ(ierr); 499a4e36a32SMatthew G. Knepley ierr = DMGlobalToLocalBegin(vdm, vel, INSERT_VALUES, locvel);CHKERRQ(ierr); 500a4e36a32SMatthew G. Knepley ierr = DMGlobalToLocalEnd(vdm, vel, INSERT_VALUES, locvel);CHKERRQ(ierr); 501a4e36a32SMatthew G. Knepley ierr = VecRestoreSubVector(u, vis, &vel);CHKERRQ(ierr); 502a4e36a32SMatthew G. Knepley ierr = ISDestroy(&vis);CHKERRQ(ierr); 503a4e36a32SMatthew G. Knepley /* Interpolate velocity */ 504a4e36a32SMatthew G. Knepley ierr = DMInterpolationCreate(PETSC_COMM_SELF, &ictx);CHKERRQ(ierr); 505a4e36a32SMatthew G. Knepley ierr = DMInterpolationSetDim(ictx, dim);CHKERRQ(ierr); 506a4e36a32SMatthew G. Knepley ierr = DMInterpolationSetDof(ictx, dim);CHKERRQ(ierr); 507a4e36a32SMatthew G. Knepley ierr = VecGetArrayRead(X, &coords);CHKERRQ(ierr); 508a4e36a32SMatthew G. Knepley ierr = DMInterpolationAddPoints(ictx, Np, (PetscReal *) coords);CHKERRQ(ierr); 509a4e36a32SMatthew G. Knepley ierr = VecRestoreArrayRead(X, &coords);CHKERRQ(ierr); 51052aa1562SMatthew G. Knepley /* Particles that lie outside the domain should be dropped, 51152aa1562SMatthew G. Knepley whereas particles that move to another partition should trigger a migration */ 51252aa1562SMatthew G. Knepley ierr = DMInterpolationSetUp(ictx, vdm, PETSC_FALSE, PETSC_TRUE);CHKERRQ(ierr); 51352aa1562SMatthew G. Knepley ierr = VecSet(pvel, 0.);CHKERRQ(ierr); 514a4e36a32SMatthew G. Knepley ierr = DMInterpolationEvaluate(ictx, vdm, locvel, pvel);CHKERRQ(ierr); 515a4e36a32SMatthew G. Knepley ierr = DMInterpolationDestroy(&ictx);CHKERRQ(ierr); 516a4e36a32SMatthew G. Knepley ierr = DMRestoreLocalVector(vdm, &locvel);CHKERRQ(ierr); 517a4e36a32SMatthew G. Knepley ierr = DMDestroy(&vdm);CHKERRQ(ierr); 518a4e36a32SMatthew G. Knepley 519a4e36a32SMatthew G. Knepley ierr = VecGetArray(F, &f);CHKERRQ(ierr); 520a4e36a32SMatthew G. Knepley ierr = VecGetArrayRead(pvel, &v);CHKERRQ(ierr); 521a4e36a32SMatthew G. Knepley ierr = PetscArraycpy(f, v, Np*dim);CHKERRQ(ierr); 522a4e36a32SMatthew G. Knepley ierr = VecRestoreArrayRead(pvel, &v);CHKERRQ(ierr); 523a4e36a32SMatthew G. Knepley ierr = VecRestoreArray(F, &f);CHKERRQ(ierr); 524a4e36a32SMatthew G. Knepley ierr = DMRestoreGlobalVector(sdm, &pvel);CHKERRQ(ierr); 525a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 526a4e36a32SMatthew G. Knepley } 527a4e36a32SMatthew G. Knepley 528a4e36a32SMatthew G. Knepley static PetscErrorCode SetInitialParticleConditions(TS ts, Vec u) 529a4e36a32SMatthew G. Knepley { 530a4e36a32SMatthew G. Knepley AppCtx *user; 531a4e36a32SMatthew G. Knepley void *ctx; 532a4e36a32SMatthew G. Knepley DM dm; 533a4e36a32SMatthew G. Knepley PetscScalar *coords; 534a4e36a32SMatthew G. Knepley PetscReal x[3], dx[3]; 535a4e36a32SMatthew G. Knepley PetscInt n[3]; 536a4e36a32SMatthew G. Knepley PetscInt Np, dim, d, i, j, k; 537a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 538a4e36a32SMatthew G. Knepley 539a4e36a32SMatthew G. Knepley PetscFunctionBegin; 540a4e36a32SMatthew G. Knepley ierr = TSGetApplicationContext(ts, &ctx);CHKERRQ(ierr); 541a4e36a32SMatthew G. Knepley user = ((AdvCtx *) ctx)->ctx; 542a4e36a32SMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 543a4e36a32SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 544a4e36a32SMatthew G. Knepley switch (user->partLayout) { 545a4e36a32SMatthew G. Knepley case PART_LAYOUT_CELL: 546a4e36a32SMatthew G. Knepley ierr = DMSwarmSetPointCoordinatesRandom(dm, user->Npc);CHKERRQ(ierr); 547a4e36a32SMatthew G. Knepley break; 548a4e36a32SMatthew G. Knepley case PART_LAYOUT_BOX: 549a4e36a32SMatthew G. Knepley Np = 1; 550a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) { 551a4e36a32SMatthew G. Knepley n[d] = user->Npb; 552a4e36a32SMatthew G. Knepley dx[d] = (user->partUpper[d] - user->partLower[d])/PetscMax(1, n[d] - 1); 553a4e36a32SMatthew G. Knepley Np *= n[d]; 554a4e36a32SMatthew G. Knepley } 555a4e36a32SMatthew G. Knepley ierr = VecGetArray(u, &coords);CHKERRQ(ierr); 556a4e36a32SMatthew G. Knepley switch (dim) { 557a4e36a32SMatthew G. Knepley case 2: 558a4e36a32SMatthew G. Knepley x[0] = user->partLower[0]; 559a4e36a32SMatthew G. Knepley for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 560a4e36a32SMatthew G. Knepley x[1] = user->partLower[1]; 561a4e36a32SMatthew G. Knepley for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 562a4e36a32SMatthew G. Knepley const PetscInt p = j*n[0] + i; 563a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 564a4e36a32SMatthew G. Knepley } 565a4e36a32SMatthew G. Knepley } 566a4e36a32SMatthew G. Knepley break; 567a4e36a32SMatthew G. Knepley case 3: 568a4e36a32SMatthew G. Knepley x[0] = user->partLower[0]; 569a4e36a32SMatthew G. Knepley for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 570a4e36a32SMatthew G. Knepley x[1] = user->partLower[1]; 571a4e36a32SMatthew G. Knepley for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 572a4e36a32SMatthew G. Knepley x[2] = user->partLower[2]; 573a4e36a32SMatthew G. Knepley for (k = 0; k < n[2]; ++k, x[2] += dx[2]) { 574a4e36a32SMatthew G. Knepley const PetscInt p = (k*n[1] + j)*n[0] + i; 575a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 576a4e36a32SMatthew G. Knepley } 577a4e36a32SMatthew G. Knepley } 578a4e36a32SMatthew G. Knepley } 579a4e36a32SMatthew G. Knepley break; 580a4e36a32SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_SUP, "Do not support particle layout in dimension %D", dim); 581a4e36a32SMatthew G. Knepley } 582a4e36a32SMatthew G. Knepley ierr = VecRestoreArray(u, &coords);CHKERRQ(ierr); 583a4e36a32SMatthew G. Knepley break; 584a4e36a32SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]); 585a4e36a32SMatthew G. Knepley } 586a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 587a4e36a32SMatthew G. Knepley } 588a4e36a32SMatthew G. Knepley 589a4e36a32SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, DM sdm, AppCtx *user) 590a4e36a32SMatthew G. Knepley { 591a4e36a32SMatthew G. Knepley DM cdm = dm; 592a4e36a32SMatthew G. Knepley PetscFE fe[3]; 593a4e36a32SMatthew G. Knepley Parameter *param; 594a4e36a32SMatthew G. Knepley PetscInt *cellid, n[3]; 595a4e36a32SMatthew G. Knepley PetscReal x[3], dx[3]; 596a4e36a32SMatthew G. Knepley PetscScalar *coords; 597a4e36a32SMatthew G. Knepley DMPolytopeType ct; 598a4e36a32SMatthew G. Knepley PetscInt dim, d, cStart, cEnd, c, Np, p, i, j, k; 599a4e36a32SMatthew G. Knepley PetscBool simplex; 600a4e36a32SMatthew G. Knepley MPI_Comm comm; 601a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 602a4e36a32SMatthew G. Knepley 603a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 604a4e36a32SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 605a4e36a32SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 606a4e36a32SMatthew G. Knepley ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 607a4e36a32SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 608a4e36a32SMatthew G. Knepley /* Create finite element */ 609a4e36a32SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 610a4e36a32SMatthew G. Knepley ierr = PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 611a4e36a32SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr); 612a4e36a32SMatthew G. Knepley 613a4e36a32SMatthew G. Knepley ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 614a4e36a32SMatthew G. Knepley ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 615a4e36a32SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr); 616a4e36a32SMatthew G. Knepley 617a4e36a32SMatthew G. Knepley ierr = PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]);CHKERRQ(ierr); 618a4e36a32SMatthew G. Knepley ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr); 619a4e36a32SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) fe[2], "temperature");CHKERRQ(ierr); 620a4e36a32SMatthew G. Knepley 621a4e36a32SMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */ 622a4e36a32SMatthew G. Knepley ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 623a4e36a32SMatthew G. Knepley ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr); 624a4e36a32SMatthew G. Knepley ierr = DMSetField(dm, 2, NULL, (PetscObject) fe[2]);CHKERRQ(ierr); 625a4e36a32SMatthew G. Knepley ierr = DMCreateDS(dm);CHKERRQ(ierr); 626a4e36a32SMatthew G. Knepley ierr = SetupProblem(dm, user);CHKERRQ(ierr); 627a4e36a32SMatthew G. Knepley ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 628a4e36a32SMatthew G. Knepley while (cdm) { 629a4e36a32SMatthew G. Knepley ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 630a4e36a32SMatthew G. Knepley ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 631a4e36a32SMatthew G. Knepley } 632a4e36a32SMatthew G. Knepley ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); 633a4e36a32SMatthew G. Knepley ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); 634a4e36a32SMatthew G. Knepley ierr = PetscFEDestroy(&fe[2]);CHKERRQ(ierr); 635a4e36a32SMatthew G. Knepley 636a4e36a32SMatthew G. Knepley { 637a4e36a32SMatthew G. Knepley PetscObject pressure; 638a4e36a32SMatthew G. Knepley MatNullSpace nullspacePres; 639a4e36a32SMatthew G. Knepley 640a4e36a32SMatthew G. Knepley ierr = DMGetField(dm, 1, NULL, &pressure);CHKERRQ(ierr); 641a4e36a32SMatthew G. Knepley ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres);CHKERRQ(ierr); 642a4e36a32SMatthew G. Knepley ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres);CHKERRQ(ierr); 643a4e36a32SMatthew G. Knepley ierr = MatNullSpaceDestroy(&nullspacePres);CHKERRQ(ierr); 644a4e36a32SMatthew G. Knepley } 645a4e36a32SMatthew G. Knepley 646a4e36a32SMatthew G. Knepley /* Setup particle information */ 647a4e36a32SMatthew G. Knepley ierr = DMSwarmSetType(sdm, DMSWARM_PIC);CHKERRQ(ierr); 648a4e36a32SMatthew G. Knepley ierr = DMSwarmRegisterPetscDatatypeField(sdm, "mass", 1, PETSC_REAL);CHKERRQ(ierr); 649a4e36a32SMatthew G. Knepley ierr = DMSwarmFinalizeFieldRegister(sdm);CHKERRQ(ierr); 650a4e36a32SMatthew G. Knepley switch (user->partLayout) { 651a4e36a32SMatthew G. Knepley case PART_LAYOUT_CELL: 652a4e36a32SMatthew G. Knepley ierr = DMSwarmSetLocalSizes(sdm, (cEnd - cStart) * user->Npc, 0);CHKERRQ(ierr); 653a4e36a32SMatthew G. Knepley ierr = DMSetFromOptions(sdm);CHKERRQ(ierr); 654a4e36a32SMatthew G. Knepley ierr = DMSwarmGetField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 655a4e36a32SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 656a4e36a32SMatthew G. Knepley for (p = 0; p < user->Npc; ++p) { 657a4e36a32SMatthew G. Knepley const PetscInt n = c*user->Npc + p; 658a4e36a32SMatthew G. Knepley 659a4e36a32SMatthew G. Knepley cellid[n] = c; 660a4e36a32SMatthew G. Knepley } 661a4e36a32SMatthew G. Knepley } 662a4e36a32SMatthew G. Knepley ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 663a4e36a32SMatthew G. Knepley ierr = DMSwarmSetPointCoordinatesRandom(sdm, user->Npc);CHKERRQ(ierr); 664a4e36a32SMatthew G. Knepley break; 665a4e36a32SMatthew G. Knepley case PART_LAYOUT_BOX: 666a4e36a32SMatthew G. Knepley Np = 1; 667a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) { 668a4e36a32SMatthew G. Knepley n[d] = user->Npb; 669a4e36a32SMatthew G. Knepley dx[d] = (user->partUpper[d] - user->partLower[d])/PetscMax(1, n[d] - 1); 670a4e36a32SMatthew G. Knepley Np *= n[d]; 671a4e36a32SMatthew G. Knepley } 672a4e36a32SMatthew G. Knepley ierr = DMSwarmSetLocalSizes(sdm, Np, 0);CHKERRQ(ierr); 673a4e36a32SMatthew G. Knepley ierr = DMSetFromOptions(sdm);CHKERRQ(ierr); 674a4e36a32SMatthew G. Knepley ierr = DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 675a4e36a32SMatthew G. Knepley switch (dim) { 676a4e36a32SMatthew G. Knepley case 2: 677a4e36a32SMatthew G. Knepley x[0] = user->partLower[0]; 678a4e36a32SMatthew G. Knepley for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 679a4e36a32SMatthew G. Knepley x[1] = user->partLower[1]; 680a4e36a32SMatthew G. Knepley for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 681a4e36a32SMatthew G. Knepley const PetscInt p = j*n[0] + i; 682a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 683a4e36a32SMatthew G. Knepley } 684a4e36a32SMatthew G. Knepley } 685a4e36a32SMatthew G. Knepley break; 686a4e36a32SMatthew G. Knepley case 3: 687a4e36a32SMatthew G. Knepley x[0] = user->partLower[0]; 688a4e36a32SMatthew G. Knepley for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 689a4e36a32SMatthew G. Knepley x[1] = user->partLower[1]; 690a4e36a32SMatthew G. Knepley for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 691a4e36a32SMatthew G. Knepley x[2] = user->partLower[2]; 692a4e36a32SMatthew G. Knepley for (k = 0; k < n[2]; ++k, x[2] += dx[2]) { 693a4e36a32SMatthew G. Knepley const PetscInt p = (k*n[1] + j)*n[0] + i; 694a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 695a4e36a32SMatthew G. Knepley } 696a4e36a32SMatthew G. Knepley } 697a4e36a32SMatthew G. Knepley } 698a4e36a32SMatthew G. Knepley break; 699a4e36a32SMatthew G. Knepley default: SETERRQ1(comm, PETSC_ERR_SUP, "Do not support particle layout in dimension %D", dim); 700a4e36a32SMatthew G. Knepley } 701a4e36a32SMatthew G. Knepley ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 702a4e36a32SMatthew G. Knepley ierr = DMSwarmGetField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 703a4e36a32SMatthew G. Knepley for (p = 0; p < Np; ++p) cellid[p] = 0; 704a4e36a32SMatthew G. Knepley ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 705a4e36a32SMatthew G. Knepley ierr = DMSwarmMigrate(sdm, PETSC_TRUE);CHKERRQ(ierr); 706a4e36a32SMatthew G. Knepley break; 707a4e36a32SMatthew G. Knepley default: SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]); 708a4e36a32SMatthew G. Knepley } 709a4e36a32SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) sdm, "Particles");CHKERRQ(ierr); 710a4e36a32SMatthew G. Knepley ierr = DMViewFromOptions(sdm, NULL, "-dm_view");CHKERRQ(ierr); 711a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 712a4e36a32SMatthew G. Knepley } 713a4e36a32SMatthew G. Knepley 714a4e36a32SMatthew G. Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) 715a4e36a32SMatthew G. Knepley { 716a4e36a32SMatthew G. Knepley Vec vec; 717a4e36a32SMatthew G. Knepley PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 718a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 719a4e36a32SMatthew G. Knepley 720a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 721a4e36a32SMatthew G. Knepley if (ofield != 1) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield); 722a4e36a32SMatthew G. Knepley funcs[nfield] = constant; 723a4e36a32SMatthew G. Knepley ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr); 724a4e36a32SMatthew G. Knepley ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr); 725a4e36a32SMatthew G. Knepley ierr = VecNormalize(vec, NULL);CHKERRQ(ierr); 726a4e36a32SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr); 727a4e36a32SMatthew G. Knepley ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr); 728a4e36a32SMatthew G. Knepley ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace);CHKERRQ(ierr); 729a4e36a32SMatthew G. Knepley ierr = VecDestroy(&vec);CHKERRQ(ierr); 730a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 731a4e36a32SMatthew G. Knepley } 732a4e36a32SMatthew G. Knepley 733a4e36a32SMatthew G. Knepley static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u) 734a4e36a32SMatthew G. Knepley { 735a4e36a32SMatthew G. Knepley DM dm; 736a4e36a32SMatthew G. Knepley MatNullSpace nullsp; 737a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 738a4e36a32SMatthew G. Knepley 739a4e36a32SMatthew G. Knepley PetscFunctionBegin; 740a4e36a32SMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 741a4e36a32SMatthew G. Knepley ierr = CreatePressureNullSpace(dm, 1, 1, &nullsp);CHKERRQ(ierr); 742a4e36a32SMatthew G. Knepley ierr = MatNullSpaceRemove(nullsp, u);CHKERRQ(ierr); 743a4e36a32SMatthew G. Knepley ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 744a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 745a4e36a32SMatthew G. Knepley } 746a4e36a32SMatthew G. Knepley 747a4e36a32SMatthew G. Knepley /* Make the discrete pressure discretely divergence free */ 748a4e36a32SMatthew G. Knepley static PetscErrorCode RemoveDiscretePressureNullspace(TS ts) 749a4e36a32SMatthew G. Knepley { 750a4e36a32SMatthew G. Knepley Vec u; 751a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 752a4e36a32SMatthew G. Knepley 753a4e36a32SMatthew G. Knepley PetscFunctionBegin; 754a4e36a32SMatthew G. Knepley ierr = TSGetSolution(ts, &u);CHKERRQ(ierr); 755a4e36a32SMatthew G. Knepley ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr); 756a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 757a4e36a32SMatthew G. Knepley } 758a4e36a32SMatthew G. Knepley 759a4e36a32SMatthew G. Knepley static PetscErrorCode SetInitialConditions(TS ts, Vec u) 760a4e36a32SMatthew G. Knepley { 761a4e36a32SMatthew G. Knepley DM dm; 762a4e36a32SMatthew G. Knepley PetscReal t; 763a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 764a4e36a32SMatthew G. Knepley 765a4e36a32SMatthew G. Knepley PetscFunctionBegin; 766a4e36a32SMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 767a4e36a32SMatthew G. Knepley ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 768a4e36a32SMatthew G. Knepley ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr); 769a4e36a32SMatthew G. Knepley ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr); 770a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 771a4e36a32SMatthew G. Knepley } 772a4e36a32SMatthew G. Knepley 773a4e36a32SMatthew G. Knepley static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) 774a4e36a32SMatthew G. Knepley { 775a4e36a32SMatthew G. Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 776a4e36a32SMatthew G. Knepley void *ctxs[3]; 777a4e36a32SMatthew G. Knepley DM dm; 778a4e36a32SMatthew G. Knepley PetscDS ds; 779a4e36a32SMatthew G. Knepley Vec v; 780a4e36a32SMatthew G. Knepley PetscReal ferrors[3]; 781b22b7e2eSMatthew G. Knepley PetscInt tl, l, f; 782a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 783a4e36a32SMatthew G. Knepley 784a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 785a4e36a32SMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 786a4e36a32SMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 787a4e36a32SMatthew G. Knepley 788a4e36a32SMatthew G. Knepley for (f = 0; f < 3; ++f) {ierr = PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]);CHKERRQ(ierr);} 789a4e36a32SMatthew G. Knepley ierr = DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors);CHKERRQ(ierr); 790b22b7e2eSMatthew G. Knepley ierr = PetscObjectGetTabLevel((PetscObject) ts, &tl);CHKERRQ(ierr); 791b22b7e2eSMatthew G. Knepley for (l = 0; l < tl; ++l) {ierr = PetscPrintf(PETSC_COMM_WORLD, "\t");CHKERRQ(ierr);} 792a4e36a32SMatthew G. Knepley ierr = 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]);CHKERRQ(ierr); 793a4e36a32SMatthew G. Knepley 794a4e36a32SMatthew G. Knepley ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr); 795a4e36a32SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr); 796a4e36a32SMatthew G. Knepley ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 797a4e36a32SMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr); 798a4e36a32SMatthew G. Knepley 799a4e36a32SMatthew G. Knepley ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr); 800a4e36a32SMatthew G. Knepley ierr = DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, v);CHKERRQ(ierr); 801a4e36a32SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) v, "Exact Solution");CHKERRQ(ierr); 802a4e36a32SMatthew G. Knepley ierr = VecViewFromOptions(v, NULL, "-exact_vec_view");CHKERRQ(ierr); 803a4e36a32SMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr); 804a4e36a32SMatthew G. Knepley 805a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 806a4e36a32SMatthew G. Knepley } 807a4e36a32SMatthew G. Knepley 808d2ae0b49SMatthew G. Knepley /* Note that adv->x0 will not be correct after migration */ 809a4e36a32SMatthew G. Knepley static PetscErrorCode ComputeParticleError(TS ts, Vec u, Vec e) 810a4e36a32SMatthew G. Knepley { 811a4e36a32SMatthew G. Knepley AdvCtx *adv; 812a4e36a32SMatthew G. Knepley DM sdm; 813a4e36a32SMatthew G. Knepley Parameter *param; 814a4e36a32SMatthew G. Knepley const PetscScalar *xp0, *xp; 815a4e36a32SMatthew G. Knepley PetscScalar *ep; 816a4e36a32SMatthew G. Knepley PetscReal time; 817a4e36a32SMatthew G. Knepley PetscInt dim, Np, p; 818a4e36a32SMatthew G. Knepley MPI_Comm comm; 819a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 820a4e36a32SMatthew G. Knepley 821a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 822a4e36a32SMatthew G. Knepley ierr = TSGetTime(ts, &time);CHKERRQ(ierr); 823a4e36a32SMatthew G. Knepley ierr = TSGetApplicationContext(ts, (void **) &adv);CHKERRQ(ierr); 824a4e36a32SMatthew G. Knepley ierr = PetscBagGetData(adv->ctx->bag, (void **) ¶m);CHKERRQ(ierr); 825a4e36a32SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 826a4e36a32SMatthew G. Knepley ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr); 827a4e36a32SMatthew G. Knepley ierr = DMGetDimension(sdm, &dim);CHKERRQ(ierr); 828a4e36a32SMatthew G. Knepley ierr = DMSwarmGetLocalSize(sdm, &Np);CHKERRQ(ierr); 829a4e36a32SMatthew G. Knepley ierr = VecGetArrayRead(adv->x0, &xp0);CHKERRQ(ierr); 830a4e36a32SMatthew G. Knepley ierr = VecGetArrayRead(u, &xp);CHKERRQ(ierr); 831a4e36a32SMatthew G. Knepley ierr = VecGetArrayWrite(e, &ep);CHKERRQ(ierr); 832a4e36a32SMatthew G. Knepley for (p = 0; p < Np; ++p) { 833a4e36a32SMatthew G. Knepley PetscScalar x[3]; 834a4e36a32SMatthew G. Knepley PetscReal x0[3]; 835a4e36a32SMatthew G. Knepley PetscInt d; 836a4e36a32SMatthew G. Knepley 837a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p*dim+d]); 838a4e36a32SMatthew G. Knepley ierr = adv->exact(dim, time, x0, 1, x, param);CHKERRQ(ierr); 839a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) ep[p*dim+d] += x[d] - xp[p*dim+d]; 840a4e36a32SMatthew G. Knepley } 841a4e36a32SMatthew G. Knepley ierr = VecRestoreArrayRead(adv->x0, &xp0);CHKERRQ(ierr); 842a4e36a32SMatthew G. Knepley ierr = VecRestoreArrayRead(u, &xp);CHKERRQ(ierr); 843a4e36a32SMatthew G. Knepley ierr = VecRestoreArrayWrite(e, &ep);CHKERRQ(ierr); 844a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 845a4e36a32SMatthew G. Knepley } 846a4e36a32SMatthew G. Knepley 847a4e36a32SMatthew G. Knepley static PetscErrorCode MonitorParticleError(TS ts, PetscInt step, PetscReal time, Vec u, void *ctx) 848a4e36a32SMatthew G. Knepley { 849a4e36a32SMatthew G. Knepley AdvCtx *adv = (AdvCtx *) ctx; 850a4e36a32SMatthew G. Knepley DM sdm; 851a4e36a32SMatthew G. Knepley Parameter *param; 852a4e36a32SMatthew G. Knepley const PetscScalar *xp0, *xp; 853a4e36a32SMatthew G. Knepley PetscReal error = 0.0; 854b22b7e2eSMatthew G. Knepley PetscInt dim, tl, l, Np, p; 855a4e36a32SMatthew G. Knepley MPI_Comm comm; 856a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 857a4e36a32SMatthew G. Knepley 858a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 859a4e36a32SMatthew G. Knepley ierr = PetscBagGetData(adv->ctx->bag, (void **) ¶m);CHKERRQ(ierr); 860a4e36a32SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 861a4e36a32SMatthew G. Knepley ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr); 862a4e36a32SMatthew G. Knepley ierr = DMGetDimension(sdm, &dim);CHKERRQ(ierr); 863a4e36a32SMatthew G. Knepley ierr = DMSwarmGetLocalSize(sdm, &Np);CHKERRQ(ierr); 864a4e36a32SMatthew G. Knepley ierr = VecGetArrayRead(adv->x0, &xp0);CHKERRQ(ierr); 865a4e36a32SMatthew G. Knepley ierr = VecGetArrayRead(u, &xp);CHKERRQ(ierr); 866a4e36a32SMatthew G. Knepley for (p = 0; p < Np; ++p) { 867a4e36a32SMatthew G. Knepley PetscScalar x[3]; 868a4e36a32SMatthew G. Knepley PetscReal x0[3]; 869a4e36a32SMatthew G. Knepley PetscReal perror = 0.0; 870a4e36a32SMatthew G. Knepley PetscInt d; 871a4e36a32SMatthew G. Knepley 872a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p*dim+d]); 873a4e36a32SMatthew G. Knepley ierr = adv->exact(dim, time, x0, 1, x, param);CHKERRQ(ierr); 874a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) perror += PetscSqr(PetscRealPart(x[d] - xp[p*dim+d])); 875a4e36a32SMatthew G. Knepley error += perror; 876a4e36a32SMatthew G. Knepley } 877a4e36a32SMatthew G. Knepley ierr = VecRestoreArrayRead(adv->x0, &xp0);CHKERRQ(ierr); 878a4e36a32SMatthew G. Knepley ierr = VecRestoreArrayRead(u, &xp);CHKERRQ(ierr); 879b22b7e2eSMatthew G. Knepley ierr = PetscObjectGetTabLevel((PetscObject) ts, &tl);CHKERRQ(ierr); 880b22b7e2eSMatthew G. Knepley for (l = 0; l < tl; ++l) {ierr = PetscPrintf(PETSC_COMM_WORLD, "\t");CHKERRQ(ierr);} 881a4e36a32SMatthew G. Knepley ierr = PetscPrintf(comm, "Timestep: %04d time = %-8.4g \t L_2 Particle Error: [%2.3g]\n", (int) step, (double) time, (double) error);CHKERRQ(ierr); 882a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 883a4e36a32SMatthew G. Knepley } 884a4e36a32SMatthew G. Knepley 885a4e36a32SMatthew G. Knepley static PetscErrorCode AdvectParticles(TS ts) 886a4e36a32SMatthew G. Knepley { 887a4e36a32SMatthew G. Knepley TS sts; 888a4e36a32SMatthew G. Knepley DM sdm; 889ad95fc09SMatthew G. Knepley Vec coordinates; 890a4e36a32SMatthew G. Knepley AdvCtx *adv; 891a4e36a32SMatthew G. Knepley PetscReal time; 892d2ae0b49SMatthew G. Knepley PetscBool lreset, reset; 893d2ae0b49SMatthew G. Knepley PetscInt dim, n, N, newn, newN; 894a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 895a4e36a32SMatthew G. Knepley 896a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 897a4e36a32SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) ts, "_SwarmTS", (PetscObject *) &sts);CHKERRQ(ierr); 898a4e36a32SMatthew G. Knepley ierr = TSGetDM(sts, &sdm);CHKERRQ(ierr); 899ad95fc09SMatthew G. Knepley ierr = TSGetRHSFunction(sts, NULL, NULL, (void **) &adv);CHKERRQ(ierr); 900d2ae0b49SMatthew G. Knepley ierr = DMGetDimension(sdm, &dim);CHKERRQ(ierr); 901ad95fc09SMatthew G. Knepley ierr = DMSwarmGetSize(sdm, &N);CHKERRQ(ierr); 902ad95fc09SMatthew G. Knepley ierr = DMSwarmGetLocalSize(sdm, &n);CHKERRQ(ierr); 903ad95fc09SMatthew G. Knepley ierr = DMSwarmCreateGlobalVectorFromField(sdm, DMSwarmPICField_coor, &coordinates);CHKERRQ(ierr); 904a4e36a32SMatthew G. Knepley ierr = TSGetTime(ts, &time);CHKERRQ(ierr); 905a4e36a32SMatthew G. Knepley ierr = TSSetMaxTime(sts, time);CHKERRQ(ierr); 906a4e36a32SMatthew G. Knepley adv->tf = time; 907ad95fc09SMatthew G. Knepley ierr = TSSolve(sts, coordinates);CHKERRQ(ierr); 908ad95fc09SMatthew G. Knepley ierr = DMSwarmDestroyGlobalVectorFromField(sdm, DMSwarmPICField_coor, &coordinates);CHKERRQ(ierr); 909a4e36a32SMatthew G. Knepley ierr = VecCopy(adv->uf, adv->ui);CHKERRQ(ierr); 910a4e36a32SMatthew G. Knepley adv->ti = adv->tf; 911a4e36a32SMatthew G. Knepley 912a4e36a32SMatthew G. Knepley ierr = DMSwarmMigrate(sdm, PETSC_TRUE);CHKERRQ(ierr); 913ad95fc09SMatthew G. Knepley ierr = DMSwarmGetSize(sdm, &newN);CHKERRQ(ierr); 914ad95fc09SMatthew G. Knepley ierr = DMSwarmGetLocalSize(sdm, &newn);CHKERRQ(ierr); 915d2ae0b49SMatthew G. Knepley lreset = (n != newn || N != newN) ? PETSC_TRUE : PETSC_FALSE; 916ad95fc09SMatthew G. Knepley ierr = MPI_Allreduce(&lreset, &reset, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject) sts));CHKERRMPI(ierr); 917ad95fc09SMatthew G. Knepley if (reset) { 918ad95fc09SMatthew G. Knepley ierr = TSReset(sts);CHKERRQ(ierr); 919ad95fc09SMatthew G. Knepley ierr = DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor);CHKERRQ(ierr); 920ad95fc09SMatthew G. Knepley } 921a4e36a32SMatthew G. Knepley ierr = DMViewFromOptions(sdm, NULL, "-dm_view");CHKERRQ(ierr); 922a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 923a4e36a32SMatthew G. Knepley } 924a4e36a32SMatthew G. Knepley 925a4e36a32SMatthew G. Knepley int main(int argc, char **argv) 926a4e36a32SMatthew G. Knepley { 927a4e36a32SMatthew G. Knepley DM dm, sdm; 928a4e36a32SMatthew G. Knepley TS ts, sts; 929d2ae0b49SMatthew G. Knepley Vec u, xtmp; 930a4e36a32SMatthew G. Knepley AppCtx user; 931a4e36a32SMatthew G. Knepley AdvCtx adv; 932a4e36a32SMatthew G. Knepley PetscReal t; 933a4e36a32SMatthew G. Knepley PetscInt dim; 934a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 935a4e36a32SMatthew G. Knepley 936a4e36a32SMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 937a4e36a32SMatthew G. Knepley ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 938a4e36a32SMatthew G. Knepley ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr); 939a4e36a32SMatthew G. Knepley ierr = SetupParameters(&user);CHKERRQ(ierr); 940a4e36a32SMatthew G. Knepley ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 941a4e36a32SMatthew G. Knepley ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 942a4e36a32SMatthew G. Knepley ierr = TSSetDM(ts, dm);CHKERRQ(ierr); 943a4e36a32SMatthew G. Knepley ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr); 944a4e36a32SMatthew G. Knepley /* Discretize chemical species */ 945a4e36a32SMatthew G. Knepley ierr = DMCreate(PETSC_COMM_WORLD, &sdm);CHKERRQ(ierr); 946a4e36a32SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) sdm, "part_");CHKERRQ(ierr); 947a4e36a32SMatthew G. Knepley ierr = DMSetType(sdm, DMSWARM);CHKERRQ(ierr); 948a4e36a32SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 949a4e36a32SMatthew G. Knepley ierr = DMSetDimension(sdm, dim);CHKERRQ(ierr); 950a4e36a32SMatthew G. Knepley ierr = DMSwarmSetCellDM(sdm, dm);CHKERRQ(ierr); 951a4e36a32SMatthew G. Knepley /* Setup problem */ 952a4e36a32SMatthew G. Knepley ierr = SetupDiscretization(dm, sdm, &user);CHKERRQ(ierr); 953a4e36a32SMatthew G. Knepley ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr); 954a4e36a32SMatthew G. Knepley 955a4e36a32SMatthew G. Knepley ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 956a4e36a32SMatthew G. Knepley ierr = DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace);CHKERRQ(ierr); 957a4e36a32SMatthew G. Knepley 958a4e36a32SMatthew G. Knepley ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr); 959a4e36a32SMatthew G. Knepley ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr); 960a4e36a32SMatthew G. Knepley ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr); 961a4e36a32SMatthew G. Knepley ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 962a4e36a32SMatthew G. Knepley ierr = TSSetPreStep(ts, RemoveDiscretePressureNullspace);CHKERRQ(ierr); 963a4e36a32SMatthew G. Knepley ierr = TSMonitorSet(ts, MonitorError, &user, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 964a4e36a32SMatthew G. Knepley ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 965a4e36a32SMatthew G. Knepley 966a4e36a32SMatthew G. Knepley ierr = TSSetComputeInitialCondition(ts, SetInitialConditions);CHKERRQ(ierr); /* Must come after SetFromOptions() */ 967a4e36a32SMatthew G. Knepley ierr = SetInitialConditions(ts, u);CHKERRQ(ierr); 968a4e36a32SMatthew G. Knepley ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 969a4e36a32SMatthew G. Knepley ierr = DMSetOutputSequenceNumber(dm, 0, t);CHKERRQ(ierr); 970a4e36a32SMatthew G. Knepley ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 971a4e36a32SMatthew G. Knepley 972a4e36a32SMatthew G. Knepley /* Setup particle position integrator */ 973a4e36a32SMatthew G. Knepley ierr = TSCreate(PETSC_COMM_WORLD, &sts);CHKERRQ(ierr); 974a4e36a32SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) sts, "part_");CHKERRQ(ierr); 975b22b7e2eSMatthew G. Knepley ierr = PetscObjectIncrementTabLevel((PetscObject) sts, (PetscObject) ts, 1);CHKERRQ(ierr); 976a4e36a32SMatthew G. Knepley ierr = TSSetDM(sts, sdm);CHKERRQ(ierr); 977a4e36a32SMatthew G. Knepley ierr = TSSetProblemType(sts, TS_NONLINEAR);CHKERRQ(ierr); 978a4e36a32SMatthew G. Knepley ierr = TSSetExactFinalTime(sts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 979a4e36a32SMatthew G. Knepley ierr = TSMonitorSet(sts, MonitorParticleError, &adv, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 980a4e36a32SMatthew G. Knepley ierr = TSSetFromOptions(sts);CHKERRQ(ierr); 981a4e36a32SMatthew G. Knepley ierr = TSSetApplicationContext(sts, &adv);CHKERRQ(ierr); 982a4e36a32SMatthew G. Knepley ierr = TSSetComputeExactError(sts, ComputeParticleError);CHKERRQ(ierr); 983a4e36a32SMatthew G. Knepley ierr = TSSetComputeInitialCondition(sts, SetInitialParticleConditions);CHKERRQ(ierr); 984a4e36a32SMatthew G. Knepley adv.ti = t; 985a4e36a32SMatthew G. Knepley adv.uf = u; 986*1e1ea65dSPierre Jolivet ierr = VecDuplicate(adv.uf, &adv.ui);CHKERRQ(ierr); 987a4e36a32SMatthew G. Knepley ierr = VecCopy(u, adv.ui);CHKERRQ(ierr); 988a4e36a32SMatthew G. Knepley ierr = TSSetRHSFunction(sts, NULL, FreeStreaming, &adv);CHKERRQ(ierr); 989a4e36a32SMatthew G. Knepley ierr = TSSetPostStep(ts, AdvectParticles);CHKERRQ(ierr); 990a4e36a32SMatthew G. Knepley ierr = PetscObjectCompose((PetscObject) ts, "_SwarmTS", (PetscObject) sts);CHKERRQ(ierr); 991a4e36a32SMatthew G. Knepley ierr = DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor);CHKERRQ(ierr); 992a4e36a32SMatthew G. Knepley ierr = DMCreateGlobalVector(sdm, &adv.x0);CHKERRQ(ierr); 993a4e36a32SMatthew G. Knepley ierr = DMSwarmCreateGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp);CHKERRQ(ierr); 994a4e36a32SMatthew G. Knepley ierr = VecCopy(xtmp, adv.x0);CHKERRQ(ierr); 995a4e36a32SMatthew G. Knepley ierr = DMSwarmDestroyGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp);CHKERRQ(ierr); 996a4e36a32SMatthew G. Knepley switch(user.solType) { 997a4e36a32SMatthew G. Knepley case SOL_TRIG_TRIG: adv.exact = trig_trig_x;break; 998a4e36a32SMatthew G. Knepley default: SETERRQ2(PetscObjectComm((PetscObject) sdm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user.solType, NUM_SOL_TYPES)], user.solType); 999a4e36a32SMatthew G. Knepley } 1000a4e36a32SMatthew G. Knepley adv.ctx = &user; 1001a4e36a32SMatthew G. Knepley 1002a4e36a32SMatthew G. Knepley ierr = TSSolve(ts, u);CHKERRQ(ierr); 1003a4e36a32SMatthew G. Knepley ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 1004a4e36a32SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr); 1005a4e36a32SMatthew G. Knepley 1006a4e36a32SMatthew G. Knepley ierr = VecDestroy(&u);CHKERRQ(ierr); 1007a4e36a32SMatthew G. Knepley ierr = VecDestroy(&adv.x0);CHKERRQ(ierr); 1008a4e36a32SMatthew G. Knepley ierr = VecDestroy(&adv.ui);CHKERRQ(ierr); 1009a4e36a32SMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 1010a4e36a32SMatthew G. Knepley ierr = DMDestroy(&sdm);CHKERRQ(ierr); 1011a4e36a32SMatthew G. Knepley ierr = TSDestroy(&ts);CHKERRQ(ierr); 1012a4e36a32SMatthew G. Knepley ierr = TSDestroy(&sts);CHKERRQ(ierr); 1013a4e36a32SMatthew G. Knepley ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr); 1014a4e36a32SMatthew G. Knepley ierr = PetscFinalize(); 1015a4e36a32SMatthew G. Knepley return ierr; 1016a4e36a32SMatthew G. Knepley } 1017a4e36a32SMatthew G. Knepley 1018a4e36a32SMatthew G. Knepley /*TEST 1019a4e36a32SMatthew G. Knepley 1020a4e36a32SMatthew G. Knepley # Swarm does not work with complex 1021a4e36a32SMatthew G. Knepley test: 1022a4e36a32SMatthew G. Knepley suffix: 2d_tri_p2_p1_p1_tconvp 1023a4e36a32SMatthew G. Knepley requires: triangle !single !complex 1024a4e36a32SMatthew G. Knepley args: -dm_plex_separate_marker -sol_type trig_trig -dm_refine 2 \ 1025a4e36a32SMatthew G. Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1026a4e36a32SMatthew G. Knepley -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 -ts_monitor_cancel \ 1027a4e36a32SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1028a4e36a32SMatthew 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 \ 1029a4e36a32SMatthew G. Knepley -fieldsplit_0_pc_type lu \ 1030a4e36a32SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \ 1031a4e36a32SMatthew G. Knepley -omega 0.5 -part_layout_type box -part_lower 0.25,0.25 -part_upper 0.75,0.75 -Npb 5 \ 1032a4e36a32SMatthew G. Knepley -part_ts_max_steps 2 -part_ts_dt 0.05 -part_ts_convergence_estimate -convest_num_refine 1 -part_ts_monitor_cancel 1033b22b7e2eSMatthew G. Knepley test: 1034b22b7e2eSMatthew G. Knepley suffix: 2d_tri_p2_p1_p1_exit 1035b22b7e2eSMatthew G. Knepley requires: triangle !single !complex 1036b22b7e2eSMatthew G. Knepley args: -dm_plex_separate_marker -sol_type trig_trig -dm_refine 1 \ 1037b22b7e2eSMatthew G. Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1038b22b7e2eSMatthew G. Knepley -dmts_check .001 -ts_max_steps 10 -ts_dt 0.1 \ 1039b22b7e2eSMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1040b22b7e2eSMatthew 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 \ 1041b22b7e2eSMatthew G. Knepley -fieldsplit_0_pc_type lu \ 1042b22b7e2eSMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \ 1043b22b7e2eSMatthew G. Knepley -omega 0.5 -part_layout_type box -part_lower 0.25,0.25 -part_upper 0.75,0.75 -Npb 5 \ 1044b22b7e2eSMatthew G. Knepley -part_ts_max_steps 20 -part_ts_dt 0.05 1045a4e36a32SMatthew G. Knepley 1046a4e36a32SMatthew G. Knepley TEST*/ 1047