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