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; 346*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-sol_type", "The solution type", "ex77.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL)); 347a4e36a32SMatthew G. Knepley options->solType = (SolType) sol; 348a4e36a32SMatthew G. Knepley pl = options->partLayout; 349*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-part_layout_type", "The particle layout type", "ex77.c", partLayoutTypes, NUM_PART_LAYOUT_TYPES, partLayoutTypes[options->partLayout], &pl, NULL)); 350a4e36a32SMatthew G. Knepley options->partLayout = (PartLayoutType) pl; 351*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-Npc", "The initial number of particles per cell", "ex77.c", options->Npc, &options->Npc, NULL)); 352a4e36a32SMatthew G. Knepley n = 3; 353*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsRealArray("-part_lower", "The lower left corner of the particle box", "ex77.c", options->partLower, &n, NULL)); 354a4e36a32SMatthew G. Knepley n = 3; 355*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsRealArray("-part_upper", "The upper right corner of the particle box", "ex77.c", options->partUpper, &n, NULL)); 356*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-Npb", "The initial number of particles per box dimension", "ex77.c", options->Npb, &options->Npb, NULL)); 3571e1ea65dSPierre 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 366a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 367a4e36a32SMatthew G. Knepley /* setup PETSc parameter bag */ 368*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(user->bag, (void **) &p)); 369*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagSetName(user->bag, "par", "Low Mach flow parameters")); 370a4e36a32SMatthew G. Knepley bag = user->bag; 371*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity")); 372*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity")); 373*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagRegisterReal(bag, &p->T_in, 1.0, "T_in", "Inlet temperature")); 374*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagRegisterReal(bag, &p->omega, 1.0, "omega", "Rotation speed in MMS benchmark")); 375a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 376a4e36a32SMatthew G. Knepley } 377a4e36a32SMatthew G. Knepley 378a4e36a32SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 379a4e36a32SMatthew G. Knepley { 380a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 381*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 382*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 383*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 384*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 385a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 386a4e36a32SMatthew G. Knepley } 387a4e36a32SMatthew G. Knepley 388a4e36a32SMatthew G. Knepley static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 389a4e36a32SMatthew G. Knepley { 390a4e36a32SMatthew G. Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 391a4e36a32SMatthew G. Knepley PetscErrorCode (*exactFuncs_t[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 392a4e36a32SMatthew G. Knepley PetscDS prob; 39345480ffeSMatthew G. Knepley DMLabel label; 394a4e36a32SMatthew G. Knepley Parameter *ctx; 395a4e36a32SMatthew G. Knepley PetscInt id; 396a4e36a32SMatthew G. Knepley 397a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 398*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 399*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &prob)); 400a4e36a32SMatthew G. Knepley switch(user->solType) { 401a4e36a32SMatthew G. Knepley case SOL_TRIG_TRIG: 402*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 0, f0_trig_trig_v, f1_v)); 403*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 2, f0_trig_trig_w, f1_w)); 404a4e36a32SMatthew G. Knepley 405a4e36a32SMatthew G. Knepley exactFuncs[0] = trig_trig_u; 406a4e36a32SMatthew G. Knepley exactFuncs[1] = trig_trig_p; 407a4e36a32SMatthew G. Knepley exactFuncs[2] = trig_trig_T; 408a4e36a32SMatthew G. Knepley exactFuncs_t[0] = trig_trig_u_t; 409a4e36a32SMatthew G. Knepley exactFuncs_t[1] = NULL; 410a4e36a32SMatthew G. Knepley exactFuncs_t[2] = trig_trig_T_t; 411a4e36a32SMatthew G. Knepley break; 41298921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 413a4e36a32SMatthew G. Knepley } 414a4e36a32SMatthew G. Knepley 415*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 1, f0_q, NULL)); 416a4e36a32SMatthew G. Knepley 417*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu)); 418*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL)); 419*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL)); 420*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL)); 421*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 2, 2, g0_wT, g1_wT, NULL, g3_wT)); 422a4e36a32SMatthew G. Knepley /* Setup constants */ 423a4e36a32SMatthew G. Knepley { 424a4e36a32SMatthew G. Knepley Parameter *param; 425a4e36a32SMatthew G. Knepley PetscScalar constants[4]; 426a4e36a32SMatthew G. Knepley 427*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(user->bag, (void **) ¶m)); 428a4e36a32SMatthew G. Knepley 429a4e36a32SMatthew G. Knepley constants[0] = param->nu; 430a4e36a32SMatthew G. Knepley constants[1] = param->alpha; 431a4e36a32SMatthew G. Knepley constants[2] = param->T_in; 432a4e36a32SMatthew G. Knepley constants[3] = param->omega; 433*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetConstants(prob, 4, constants)); 434a4e36a32SMatthew G. Knepley } 435a4e36a32SMatthew G. Knepley /* Setup Boundary Conditions */ 436*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(user->bag, (void **) &ctx)); 437a4e36a32SMatthew G. Knepley id = 3; 438*5f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 439a4e36a32SMatthew G. Knepley id = 1; 440*5f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 441a4e36a32SMatthew G. Knepley id = 2; 442*5f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 443a4e36a32SMatthew G. Knepley id = 4; 444*5f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 445a4e36a32SMatthew G. Knepley id = 3; 446*5f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 447a4e36a32SMatthew G. Knepley id = 1; 448*5f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 449a4e36a32SMatthew G. Knepley id = 2; 450*5f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 451a4e36a32SMatthew G. Knepley id = 4; 452*5f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 453a4e36a32SMatthew G. Knepley 454a4e36a32SMatthew G. Knepley /*setup exact solution.*/ 455*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx)); 456*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx)); 457*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx)); 458*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolutionTimeDerivative(prob, 0, exactFuncs_t[0], ctx)); 459*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolutionTimeDerivative(prob, 1, exactFuncs_t[1], ctx)); 460*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolutionTimeDerivative(prob, 2, exactFuncs_t[2], ctx)); 461a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 462a4e36a32SMatthew G. Knepley } 463a4e36a32SMatthew G. Knepley 464a4e36a32SMatthew G. Knepley /* x_t = v 465a4e36a32SMatthew G. Knepley 466a4e36a32SMatthew G. Knepley Note that here we use the velocity field at t_{n+1} to advect the particles from 467a4e36a32SMatthew G. Knepley t_n to t_{n+1}. If we use both of these fields, we could use Crank-Nicholson or 468a4e36a32SMatthew G. Knepley the method of characteristics. 469a4e36a32SMatthew G. Knepley */ 470a4e36a32SMatthew G. Knepley static PetscErrorCode FreeStreaming(TS ts, PetscReal t, Vec X, Vec F, void *ctx) 471a4e36a32SMatthew G. Knepley { 472a4e36a32SMatthew G. Knepley AdvCtx *adv = (AdvCtx *) ctx; 473a4e36a32SMatthew G. Knepley Vec u = adv->ui; 474a4e36a32SMatthew G. Knepley DM sdm, dm, vdm; 475a4e36a32SMatthew G. Knepley Vec vel, locvel, pvel; 476a4e36a32SMatthew G. Knepley IS vis; 477a4e36a32SMatthew G. Knepley DMInterpolationInfo ictx; 478a4e36a32SMatthew G. Knepley const PetscScalar *coords, *v; 479a4e36a32SMatthew G. Knepley PetscScalar *f; 480a4e36a32SMatthew G. Knepley PetscInt vf[1] = {0}; 481a4e36a32SMatthew G. Knepley PetscInt dim, Np; 482a4e36a32SMatthew G. Knepley 483a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 484*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &sdm)); 485*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(sdm, &dm)); 486*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(sdm, &pvel)); 487*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetLocalSize(sdm, &Np)); 488*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 489a4e36a32SMatthew G. Knepley /* Get local velocity */ 490*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateSubDM(dm, 1, vf, &vis, &vdm)); 491*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSubVector(u, vis, &vel)); 492*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(vdm, &locvel)); 493*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInsertBoundaryValues(vdm, PETSC_TRUE, locvel, adv->ti, NULL, NULL, NULL)); 494*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(vdm, vel, INSERT_VALUES, locvel)); 495*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(vdm, vel, INSERT_VALUES, locvel)); 496*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreSubVector(u, vis, &vel)); 497*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&vis)); 498a4e36a32SMatthew G. Knepley /* Interpolate velocity */ 499*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationCreate(PETSC_COMM_SELF, &ictx)); 500*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationSetDim(ictx, dim)); 501*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationSetDof(ictx, dim)); 502*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X, &coords)); 503*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationAddPoints(ictx, Np, (PetscReal *) coords)); 504*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X, &coords)); 50552aa1562SMatthew G. Knepley /* Particles that lie outside the domain should be dropped, 50652aa1562SMatthew G. Knepley whereas particles that move to another partition should trigger a migration */ 507*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationSetUp(ictx, vdm, PETSC_FALSE, PETSC_TRUE)); 508*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(pvel, 0.)); 509*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationEvaluate(ictx, vdm, locvel, pvel)); 510*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationDestroy(&ictx)); 511*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(vdm, &locvel)); 512*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&vdm)); 513a4e36a32SMatthew G. Knepley 514*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(F, &f)); 515*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(pvel, &v)); 516*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(f, v, Np*dim)); 517*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(pvel, &v)); 518*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(F, &f)); 519*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(sdm, &pvel)); 520a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 521a4e36a32SMatthew G. Knepley } 522a4e36a32SMatthew G. Knepley 523a4e36a32SMatthew G. Knepley static PetscErrorCode SetInitialParticleConditions(TS ts, Vec u) 524a4e36a32SMatthew G. Knepley { 525a4e36a32SMatthew G. Knepley AppCtx *user; 526a4e36a32SMatthew G. Knepley void *ctx; 527a4e36a32SMatthew G. Knepley DM dm; 528a4e36a32SMatthew G. Knepley PetscScalar *coords; 529a4e36a32SMatthew G. Knepley PetscReal x[3], dx[3]; 530a4e36a32SMatthew G. Knepley PetscInt n[3]; 531*5f80ce2aSJacob Faibussowitsch PetscInt dim, d, i, j, k; 532a4e36a32SMatthew G. Knepley 533a4e36a32SMatthew G. Knepley PetscFunctionBegin; 534*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetApplicationContext(ts, &ctx)); 535a4e36a32SMatthew G. Knepley user = ((AdvCtx *) ctx)->ctx; 536*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dm)); 537*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 538a4e36a32SMatthew G. Knepley switch (user->partLayout) { 539a4e36a32SMatthew G. Knepley case PART_LAYOUT_CELL: 540*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetPointCoordinatesRandom(dm, user->Npc)); 541a4e36a32SMatthew G. Knepley break; 542a4e36a32SMatthew G. Knepley case PART_LAYOUT_BOX: 543a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) { 544a4e36a32SMatthew G. Knepley n[d] = user->Npb; 545a4e36a32SMatthew G. Knepley dx[d] = (user->partUpper[d] - user->partLower[d])/PetscMax(1, n[d] - 1); 546a4e36a32SMatthew G. Knepley } 547*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(u, &coords)); 548a4e36a32SMatthew G. Knepley switch (dim) { 549a4e36a32SMatthew G. Knepley case 2: 550a4e36a32SMatthew G. Knepley x[0] = user->partLower[0]; 551a4e36a32SMatthew G. Knepley for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 552a4e36a32SMatthew G. Knepley x[1] = user->partLower[1]; 553a4e36a32SMatthew G. Knepley for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 554a4e36a32SMatthew G. Knepley const PetscInt p = j*n[0] + i; 555a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 556a4e36a32SMatthew G. Knepley } 557a4e36a32SMatthew G. Knepley } 558a4e36a32SMatthew G. Knepley break; 559a4e36a32SMatthew G. Knepley case 3: 560a4e36a32SMatthew G. Knepley x[0] = user->partLower[0]; 561a4e36a32SMatthew G. Knepley for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 562a4e36a32SMatthew G. Knepley x[1] = user->partLower[1]; 563a4e36a32SMatthew G. Knepley for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 564a4e36a32SMatthew G. Knepley x[2] = user->partLower[2]; 565a4e36a32SMatthew G. Knepley for (k = 0; k < n[2]; ++k, x[2] += dx[2]) { 566a4e36a32SMatthew G. Knepley const PetscInt p = (k*n[1] + j)*n[0] + i; 567a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 568a4e36a32SMatthew G. Knepley } 569a4e36a32SMatthew G. Knepley } 570a4e36a32SMatthew G. Knepley } 571a4e36a32SMatthew G. Knepley break; 57298921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ts), PETSC_ERR_SUP, "Do not support particle layout in dimension %D", dim); 573a4e36a32SMatthew G. Knepley } 574*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(u, &coords)); 575a4e36a32SMatthew G. Knepley break; 57698921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]); 577a4e36a32SMatthew G. Knepley } 578a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 579a4e36a32SMatthew G. Knepley } 580a4e36a32SMatthew G. Knepley 581a4e36a32SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, DM sdm, AppCtx *user) 582a4e36a32SMatthew G. Knepley { 583a4e36a32SMatthew G. Knepley DM cdm = dm; 584a4e36a32SMatthew G. Knepley PetscFE fe[3]; 585a4e36a32SMatthew G. Knepley Parameter *param; 586a4e36a32SMatthew G. Knepley PetscInt *cellid, n[3]; 587a4e36a32SMatthew G. Knepley PetscReal x[3], dx[3]; 588a4e36a32SMatthew G. Knepley PetscScalar *coords; 589a4e36a32SMatthew G. Knepley DMPolytopeType ct; 590a4e36a32SMatthew G. Knepley PetscInt dim, d, cStart, cEnd, c, Np, p, i, j, k; 591a4e36a32SMatthew G. Knepley PetscBool simplex; 592a4e36a32SMatthew G. Knepley MPI_Comm comm; 593a4e36a32SMatthew G. Knepley 594a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 595*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 596*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 597*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dm, cStart, &ct)); 598a4e36a32SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 599a4e36a32SMatthew G. Knepley /* Create finite element */ 600*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 601*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0])); 602*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe[0], "velocity")); 603a4e36a32SMatthew G. Knepley 604*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1])); 605*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECopyQuadrature(fe[0], fe[1])); 606*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe[1], "pressure")); 607a4e36a32SMatthew G. Knepley 608*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2])); 609*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECopyQuadrature(fe[0], fe[2])); 610*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe[2], "temperature")); 611a4e36a32SMatthew G. Knepley 612a4e36a32SMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */ 613*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe[0])); 614*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 1, NULL, (PetscObject) fe[1])); 615*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 2, NULL, (PetscObject) fe[2])); 616*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 617*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupProblem(dm, user)); 618*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(user->bag, (void **) ¶m)); 619a4e36a32SMatthew G. Knepley while (cdm) { 620*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCopyDisc(dm, cdm)); 621*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 622a4e36a32SMatthew G. Knepley } 623*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe[0])); 624*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe[1])); 625*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe[2])); 626a4e36a32SMatthew G. Knepley 627a4e36a32SMatthew G. Knepley { 628a4e36a32SMatthew G. Knepley PetscObject pressure; 629a4e36a32SMatthew G. Knepley MatNullSpace nullspacePres; 630a4e36a32SMatthew G. Knepley 631*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetField(dm, 1, NULL, &pressure)); 632*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres)); 633*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres)); 634*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceDestroy(&nullspacePres)); 635a4e36a32SMatthew G. Knepley } 636a4e36a32SMatthew G. Knepley 637a4e36a32SMatthew G. Knepley /* Setup particle information */ 638*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetType(sdm, DMSWARM_PIC)); 639*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRegisterPetscDatatypeField(sdm, "mass", 1, PETSC_REAL)); 640*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmFinalizeFieldRegister(sdm)); 641a4e36a32SMatthew G. Knepley switch (user->partLayout) { 642a4e36a32SMatthew G. Knepley case PART_LAYOUT_CELL: 643*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetLocalSizes(sdm, (cEnd - cStart) * user->Npc, 0)); 644*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(sdm)); 645*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 646a4e36a32SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 647a4e36a32SMatthew G. Knepley for (p = 0; p < user->Npc; ++p) { 648a4e36a32SMatthew G. Knepley const PetscInt n = c*user->Npc + p; 649a4e36a32SMatthew G. Knepley 650a4e36a32SMatthew G. Knepley cellid[n] = c; 651a4e36a32SMatthew G. Knepley } 652a4e36a32SMatthew G. Knepley } 653*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 654*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetPointCoordinatesRandom(sdm, user->Npc)); 655a4e36a32SMatthew G. Knepley break; 656a4e36a32SMatthew G. Knepley case PART_LAYOUT_BOX: 657a4e36a32SMatthew G. Knepley Np = 1; 658a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) { 659a4e36a32SMatthew G. Knepley n[d] = user->Npb; 660a4e36a32SMatthew G. Knepley dx[d] = (user->partUpper[d] - user->partLower[d])/PetscMax(1, n[d] - 1); 661a4e36a32SMatthew G. Knepley Np *= n[d]; 662a4e36a32SMatthew G. Knepley } 663*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetLocalSizes(sdm, Np, 0)); 664*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(sdm)); 665*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 666a4e36a32SMatthew G. Knepley switch (dim) { 667a4e36a32SMatthew G. Knepley case 2: 668a4e36a32SMatthew G. Knepley x[0] = user->partLower[0]; 669a4e36a32SMatthew G. Knepley for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 670a4e36a32SMatthew G. Knepley x[1] = user->partLower[1]; 671a4e36a32SMatthew G. Knepley for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 672a4e36a32SMatthew G. Knepley const PetscInt p = j*n[0] + i; 673a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 674a4e36a32SMatthew G. Knepley } 675a4e36a32SMatthew G. Knepley } 676a4e36a32SMatthew G. Knepley break; 677a4e36a32SMatthew G. Knepley case 3: 678a4e36a32SMatthew G. Knepley x[0] = user->partLower[0]; 679a4e36a32SMatthew G. Knepley for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 680a4e36a32SMatthew G. Knepley x[1] = user->partLower[1]; 681a4e36a32SMatthew G. Knepley for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 682a4e36a32SMatthew G. Knepley x[2] = user->partLower[2]; 683a4e36a32SMatthew G. Knepley for (k = 0; k < n[2]; ++k, x[2] += dx[2]) { 684a4e36a32SMatthew G. Knepley const PetscInt p = (k*n[1] + j)*n[0] + i; 685a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 686a4e36a32SMatthew G. Knepley } 687a4e36a32SMatthew G. Knepley } 688a4e36a32SMatthew G. Knepley } 689a4e36a32SMatthew G. Knepley break; 69098921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_SUP, "Do not support particle layout in dimension %D", dim); 691a4e36a32SMatthew G. Knepley } 692*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 693*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 694a4e36a32SMatthew G. Knepley for (p = 0; p < Np; ++p) cellid[p] = 0; 695*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 696*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmMigrate(sdm, PETSC_TRUE)); 697a4e36a32SMatthew G. Knepley break; 69898921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]); 699a4e36a32SMatthew G. Knepley } 700*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) sdm, "Particles")); 701*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(sdm, NULL, "-dm_view")); 702a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 703a4e36a32SMatthew G. Knepley } 704a4e36a32SMatthew G. Knepley 705a4e36a32SMatthew G. Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) 706a4e36a32SMatthew G. Knepley { 707a4e36a32SMatthew G. Knepley Vec vec; 708a4e36a32SMatthew G. Knepley PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 709a4e36a32SMatthew G. Knepley 710a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 7113c633725SBarry Smith PetscCheck(ofield == 1,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield); 712a4e36a32SMatthew G. Knepley funcs[nfield] = constant; 713*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &vec)); 714*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec)); 715*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNormalize(vec, NULL)); 716*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) vec, "Pressure Null Space")); 717*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(vec, NULL, "-pressure_nullspace_view")); 718*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace)); 719*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&vec)); 720a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 721a4e36a32SMatthew G. Knepley } 722a4e36a32SMatthew G. Knepley 723a4e36a32SMatthew G. Knepley static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u) 724a4e36a32SMatthew G. Knepley { 725a4e36a32SMatthew G. Knepley DM dm; 726a4e36a32SMatthew G. Knepley MatNullSpace nullsp; 727a4e36a32SMatthew G. Knepley 728a4e36a32SMatthew G. Knepley PetscFunctionBegin; 729*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dm)); 730*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreatePressureNullSpace(dm, 1, 1, &nullsp)); 731*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceRemove(nullsp, u)); 732*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceDestroy(&nullsp)); 733a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 734a4e36a32SMatthew G. Knepley } 735a4e36a32SMatthew G. Knepley 736a4e36a32SMatthew G. Knepley /* Make the discrete pressure discretely divergence free */ 737a4e36a32SMatthew G. Knepley static PetscErrorCode RemoveDiscretePressureNullspace(TS ts) 738a4e36a32SMatthew G. Knepley { 739a4e36a32SMatthew G. Knepley Vec u; 740a4e36a32SMatthew G. Knepley 741a4e36a32SMatthew G. Knepley PetscFunctionBegin; 742*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolution(ts, &u)); 743*5f80ce2aSJacob Faibussowitsch CHKERRQ(RemoveDiscretePressureNullspace_Private(ts, u)); 744a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 745a4e36a32SMatthew G. Knepley } 746a4e36a32SMatthew G. Knepley 747a4e36a32SMatthew G. Knepley static PetscErrorCode SetInitialConditions(TS ts, Vec u) 748a4e36a32SMatthew G. Knepley { 749a4e36a32SMatthew G. Knepley DM dm; 750a4e36a32SMatthew G. Knepley PetscReal t; 751a4e36a32SMatthew G. Knepley 752a4e36a32SMatthew G. Knepley PetscFunctionBegin; 753*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dm)); 754*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTime(ts, &t)); 755*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMComputeExactSolution(dm, t, u, NULL)); 756*5f80ce2aSJacob Faibussowitsch CHKERRQ(RemoveDiscretePressureNullspace_Private(ts, u)); 757a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 758a4e36a32SMatthew G. Knepley } 759a4e36a32SMatthew G. Knepley 760a4e36a32SMatthew G. Knepley static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) 761a4e36a32SMatthew G. Knepley { 762a4e36a32SMatthew G. Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 763a4e36a32SMatthew G. Knepley void *ctxs[3]; 764a4e36a32SMatthew G. Knepley DM dm; 765a4e36a32SMatthew G. Knepley PetscDS ds; 766a4e36a32SMatthew G. Knepley Vec v; 767a4e36a32SMatthew G. Knepley PetscReal ferrors[3]; 768b22b7e2eSMatthew G. Knepley PetscInt tl, l, f; 769a4e36a32SMatthew G. Knepley 770a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 771*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dm)); 772*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 773a4e36a32SMatthew G. Knepley 774*5f80ce2aSJacob Faibussowitsch for (f = 0; f < 3; ++f) CHKERRQ(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f])); 775*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors)); 776*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetTabLevel((PetscObject) ts, &tl)); 777*5f80ce2aSJacob Faibussowitsch for (l = 0; l < tl; ++l) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\t")); 778*5f80ce2aSJacob Faibussowitsch CHKERRQ(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])); 779a4e36a32SMatthew G. Knepley 780*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dm, &u)); 781*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) u, "Numerical Solution")); 782*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-sol_vec_view")); 783*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dm, &u)); 784a4e36a32SMatthew G. Knepley 785*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dm, &v)); 786*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, v)); 787*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) v, "Exact Solution")); 788*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(v, NULL, "-exact_vec_view")); 789*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dm, &v)); 790a4e36a32SMatthew G. Knepley 791a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 792a4e36a32SMatthew G. Knepley } 793a4e36a32SMatthew G. Knepley 794d2ae0b49SMatthew G. Knepley /* Note that adv->x0 will not be correct after migration */ 795a4e36a32SMatthew G. Knepley static PetscErrorCode ComputeParticleError(TS ts, Vec u, Vec e) 796a4e36a32SMatthew G. Knepley { 797a4e36a32SMatthew G. Knepley AdvCtx *adv; 798a4e36a32SMatthew G. Knepley DM sdm; 799a4e36a32SMatthew G. Knepley Parameter *param; 800a4e36a32SMatthew G. Knepley const PetscScalar *xp0, *xp; 801a4e36a32SMatthew G. Knepley PetscScalar *ep; 802a4e36a32SMatthew G. Knepley PetscReal time; 803a4e36a32SMatthew G. Knepley PetscInt dim, Np, p; 804a4e36a32SMatthew G. Knepley MPI_Comm comm; 805a4e36a32SMatthew G. Knepley 806a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 807*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTime(ts, &time)); 808*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetApplicationContext(ts, &adv)); 809*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(adv->ctx->bag, (void **) ¶m)); 810*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) ts, &comm)); 811*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &sdm)); 812*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(sdm, &dim)); 813*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetLocalSize(sdm, &Np)); 814*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(adv->x0, &xp0)); 815*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(u, &xp)); 816*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayWrite(e, &ep)); 817a4e36a32SMatthew G. Knepley for (p = 0; p < Np; ++p) { 818a4e36a32SMatthew G. Knepley PetscScalar x[3]; 819a4e36a32SMatthew G. Knepley PetscReal x0[3]; 820a4e36a32SMatthew G. Knepley PetscInt d; 821a4e36a32SMatthew G. Knepley 822a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p*dim+d]); 823*5f80ce2aSJacob Faibussowitsch CHKERRQ(adv->exact(dim, time, x0, 1, x, param)); 824a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) ep[p*dim+d] += x[d] - xp[p*dim+d]; 825a4e36a32SMatthew G. Knepley } 826*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(adv->x0, &xp0)); 827*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(u, &xp)); 828*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayWrite(e, &ep)); 829a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 830a4e36a32SMatthew G. Knepley } 831a4e36a32SMatthew G. Knepley 832a4e36a32SMatthew G. Knepley static PetscErrorCode MonitorParticleError(TS ts, PetscInt step, PetscReal time, Vec u, void *ctx) 833a4e36a32SMatthew G. Knepley { 834a4e36a32SMatthew G. Knepley AdvCtx *adv = (AdvCtx *) ctx; 835a4e36a32SMatthew G. Knepley DM sdm; 836a4e36a32SMatthew G. Knepley Parameter *param; 837a4e36a32SMatthew G. Knepley const PetscScalar *xp0, *xp; 838a4e36a32SMatthew G. Knepley PetscReal error = 0.0; 839b22b7e2eSMatthew G. Knepley PetscInt dim, tl, l, Np, p; 840a4e36a32SMatthew G. Knepley MPI_Comm comm; 841a4e36a32SMatthew G. Knepley 842a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 843*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(adv->ctx->bag, (void **) ¶m)); 844*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) ts, &comm)); 845*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &sdm)); 846*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(sdm, &dim)); 847*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetLocalSize(sdm, &Np)); 848*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(adv->x0, &xp0)); 849*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(u, &xp)); 850a4e36a32SMatthew G. Knepley for (p = 0; p < Np; ++p) { 851a4e36a32SMatthew G. Knepley PetscScalar x[3]; 852a4e36a32SMatthew G. Knepley PetscReal x0[3]; 853a4e36a32SMatthew G. Knepley PetscReal perror = 0.0; 854a4e36a32SMatthew G. Knepley PetscInt d; 855a4e36a32SMatthew G. Knepley 856a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p*dim+d]); 857*5f80ce2aSJacob Faibussowitsch CHKERRQ(adv->exact(dim, time, x0, 1, x, param)); 858a4e36a32SMatthew G. Knepley for (d = 0; d < dim; ++d) perror += PetscSqr(PetscRealPart(x[d] - xp[p*dim+d])); 859a4e36a32SMatthew G. Knepley error += perror; 860a4e36a32SMatthew G. Knepley } 861*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(adv->x0, &xp0)); 862*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(u, &xp)); 863*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetTabLevel((PetscObject) ts, &tl)); 864*5f80ce2aSJacob Faibussowitsch for (l = 0; l < tl; ++l) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\t")); 865*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(comm, "Timestep: %04d time = %-8.4g \t L_2 Particle Error: [%2.3g]\n", (int) step, (double) time, (double) error)); 866a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 867a4e36a32SMatthew G. Knepley } 868a4e36a32SMatthew G. Knepley 869a4e36a32SMatthew G. Knepley static PetscErrorCode AdvectParticles(TS ts) 870a4e36a32SMatthew G. Knepley { 871a4e36a32SMatthew G. Knepley TS sts; 872a4e36a32SMatthew G. Knepley DM sdm; 873ad95fc09SMatthew G. Knepley Vec coordinates; 874a4e36a32SMatthew G. Knepley AdvCtx *adv; 875a4e36a32SMatthew G. Knepley PetscReal time; 876d2ae0b49SMatthew G. Knepley PetscBool lreset, reset; 877d2ae0b49SMatthew G. Knepley PetscInt dim, n, N, newn, newN; 878a4e36a32SMatthew G. Knepley 879a4e36a32SMatthew G. Knepley PetscFunctionBeginUser; 880*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectQuery((PetscObject) ts, "_SwarmTS", (PetscObject *) &sts)); 881*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(sts, &sdm)); 882*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetRHSFunction(sts, NULL, NULL, (void **) &adv)); 883*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(sdm, &dim)); 884*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetSize(sdm, &N)); 885*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetLocalSize(sdm, &n)); 886*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmCreateGlobalVectorFromField(sdm, DMSwarmPICField_coor, &coordinates)); 887*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTime(ts, &time)); 888*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(sts, time)); 889a4e36a32SMatthew G. Knepley adv->tf = time; 890*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(sts, coordinates)); 891*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmDestroyGlobalVectorFromField(sdm, DMSwarmPICField_coor, &coordinates)); 892*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(adv->uf, adv->ui)); 893a4e36a32SMatthew G. Knepley adv->ti = adv->tf; 894a4e36a32SMatthew G. Knepley 895*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmMigrate(sdm, PETSC_TRUE)); 896*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetSize(sdm, &newN)); 897*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetLocalSize(sdm, &newn)); 898d2ae0b49SMatthew G. Knepley lreset = (n != newn || N != newN) ? PETSC_TRUE : PETSC_FALSE; 899*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(&lreset, &reset, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject) sts))); 900ad95fc09SMatthew G. Knepley if (reset) { 901*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSReset(sts)); 902*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor)); 903ad95fc09SMatthew G. Knepley } 904*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(sdm, NULL, "-dm_view")); 905a4e36a32SMatthew G. Knepley PetscFunctionReturn(0); 906a4e36a32SMatthew G. Knepley } 907a4e36a32SMatthew G. Knepley 908a4e36a32SMatthew G. Knepley int main(int argc, char **argv) 909a4e36a32SMatthew G. Knepley { 910a4e36a32SMatthew G. Knepley DM dm, sdm; 911a4e36a32SMatthew G. Knepley TS ts, sts; 912d2ae0b49SMatthew G. Knepley Vec u, xtmp; 913a4e36a32SMatthew G. Knepley AppCtx user; 914a4e36a32SMatthew G. Knepley AdvCtx adv; 915a4e36a32SMatthew G. Knepley PetscReal t; 916a4e36a32SMatthew G. Knepley PetscInt dim; 917a4e36a32SMatthew G. Knepley PetscErrorCode ierr; 918a4e36a32SMatthew G. Knepley 919a4e36a32SMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 920*5f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 921*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag)); 922*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupParameters(&user)); 923*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD, &ts)); 924*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 925*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts, dm)); 926*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(dm, &user)); 927a4e36a32SMatthew G. Knepley /* Discretize chemical species */ 928*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PETSC_COMM_WORLD, &sdm)); 929*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) sdm, "part_")); 930*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(sdm, DMSWARM)); 931*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 932*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(sdm, dim)); 933*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetCellDM(sdm, dm)); 934a4e36a32SMatthew G. Knepley /* Setup problem */ 935*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(dm, sdm, &user)); 936*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateClosureIndex(dm, NULL)); 937a4e36a32SMatthew G. Knepley 938*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &u)); 939*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace)); 940a4e36a32SMatthew G. Knepley 941*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 942*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 943*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 944*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 945*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetPreStep(ts, RemoveDiscretePressureNullspace)); 946*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(ts, MonitorError, &user, NULL)); 947*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 948a4e36a32SMatthew G. Knepley 949*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetComputeInitialCondition(ts, SetInitialConditions)); /* Must come after SetFromOptions() */ 950*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetInitialConditions(ts, u)); 951*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTime(ts, &t)); 952*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOutputSequenceNumber(dm, 0, t)); 953*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSCheckFromOptions(ts, u)); 954a4e36a32SMatthew G. Knepley 955a4e36a32SMatthew G. Knepley /* Setup particle position integrator */ 956*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD, &sts)); 957*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) sts, "part_")); 958*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectIncrementTabLevel((PetscObject) sts, (PetscObject) ts, 1)); 959*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(sts, sdm)); 960*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(sts, TS_NONLINEAR)); 961*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(sts, TS_EXACTFINALTIME_MATCHSTEP)); 962*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(sts, MonitorParticleError, &adv, NULL)); 963*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(sts)); 964*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetApplicationContext(sts, &adv)); 965*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetComputeExactError(sts, ComputeParticleError)); 966*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetComputeInitialCondition(sts, SetInitialParticleConditions)); 967a4e36a32SMatthew G. Knepley adv.ti = t; 968a4e36a32SMatthew G. Knepley adv.uf = u; 969*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(adv.uf, &adv.ui)); 970*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(u, adv.ui)); 971*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(sts, NULL, FreeStreaming, &adv)); 972*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetPostStep(ts, AdvectParticles)); 973*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectCompose((PetscObject) ts, "_SwarmTS", (PetscObject) sts)); 974*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor)); 975*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(sdm, &adv.x0)); 976*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmCreateGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp)); 977*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(xtmp, adv.x0)); 978*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmDestroyGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp)); 979a4e36a32SMatthew G. Knepley switch(user.solType) { 980a4e36a32SMatthew G. Knepley case SOL_TRIG_TRIG: adv.exact = trig_trig_x;break; 98198921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) sdm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user.solType, NUM_SOL_TYPES)], user.solType); 982a4e36a32SMatthew G. Knepley } 983a4e36a32SMatthew G. Knepley adv.ctx = &user; 984a4e36a32SMatthew G. Knepley 985*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts, u)); 986*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSCheckFromOptions(ts, u)); 987*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) u, "Numerical Solution")); 988a4e36a32SMatthew G. Knepley 989*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 990*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&adv.x0)); 991*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&adv.ui)); 992*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 993*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&sdm)); 994*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 995*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&sts)); 996*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagDestroy(&user.bag)); 997a4e36a32SMatthew G. Knepley ierr = PetscFinalize(); 998a4e36a32SMatthew G. Knepley return ierr; 999a4e36a32SMatthew G. Knepley } 1000a4e36a32SMatthew G. Knepley 1001a4e36a32SMatthew G. Knepley /*TEST 1002a4e36a32SMatthew G. Knepley 1003a4e36a32SMatthew G. Knepley # Swarm does not work with complex 1004a4e36a32SMatthew G. Knepley test: 1005a4e36a32SMatthew G. Knepley suffix: 2d_tri_p2_p1_p1_tconvp 1006a4e36a32SMatthew G. Knepley requires: triangle !single !complex 1007a4e36a32SMatthew G. Knepley args: -dm_plex_separate_marker -sol_type trig_trig -dm_refine 2 \ 1008a4e36a32SMatthew G. Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1009a4e36a32SMatthew G. Knepley -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 -ts_monitor_cancel \ 1010a4e36a32SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1011a4e36a32SMatthew 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 \ 1012a4e36a32SMatthew G. Knepley -fieldsplit_0_pc_type lu \ 1013a4e36a32SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \ 1014a4e36a32SMatthew G. Knepley -omega 0.5 -part_layout_type box -part_lower 0.25,0.25 -part_upper 0.75,0.75 -Npb 5 \ 1015a4e36a32SMatthew G. Knepley -part_ts_max_steps 2 -part_ts_dt 0.05 -part_ts_convergence_estimate -convest_num_refine 1 -part_ts_monitor_cancel 1016b22b7e2eSMatthew G. Knepley test: 1017b22b7e2eSMatthew G. Knepley suffix: 2d_tri_p2_p1_p1_exit 1018b22b7e2eSMatthew G. Knepley requires: triangle !single !complex 1019b22b7e2eSMatthew G. Knepley args: -dm_plex_separate_marker -sol_type trig_trig -dm_refine 1 \ 1020b22b7e2eSMatthew G. Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1021b22b7e2eSMatthew G. Knepley -dmts_check .001 -ts_max_steps 10 -ts_dt 0.1 \ 1022b22b7e2eSMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1023b22b7e2eSMatthew 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 \ 1024b22b7e2eSMatthew G. Knepley -fieldsplit_0_pc_type lu \ 1025b22b7e2eSMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \ 1026b22b7e2eSMatthew G. Knepley -omega 0.5 -part_layout_type box -part_lower 0.25,0.25 -part_upper 0.75,0.75 -Npb 5 \ 1027b22b7e2eSMatthew G. Knepley -part_ts_max_steps 20 -part_ts_dt 0.05 1028a4e36a32SMatthew G. Knepley 1029a4e36a32SMatthew G. Knepley TEST*/ 1030