1*649ef022SMatthew Knepley static char help[] = "Time-dependent Low Mach Flow in 2d and 3d channels with finite elements.\n\ 2*649ef022SMatthew Knepley We solve the Low Mach flow problem in a rectangular\n\ 3*649ef022SMatthew Knepley domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n"; 4*649ef022SMatthew Knepley 5*649ef022SMatthew Knepley /*F 6*649ef022SMatthew Knepley This Low Mach flow is time-dependent isoviscous Navier-Stokes flow. We discretize using the 7*649ef022SMatthew Knepley finite element method on an unstructured mesh. The weak form equations are 8*649ef022SMatthew Knepley 9*649ef022SMatthew Knepley \begin{align*} 10*649ef022SMatthew Knepley < q, \nabla\cdot u > = 0 11*649ef022SMatthew Knepley <v, du/dt> + <v, u \cdot \nabla u> + < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p > - < v, f > = 0 12*649ef022SMatthew Knepley < w, u \cdot \nabla T > + < \nabla w, \alpha \nabla T > - < w, Q > = 0 13*649ef022SMatthew Knepley \end{align*} 14*649ef022SMatthew Knepley 15*649ef022SMatthew Knepley where $\nu$ is the kinematic viscosity and $\alpha$ is thermal diffusivity. 16*649ef022SMatthew Knepley 17*649ef022SMatthew Knepley For visualization, use 18*649ef022SMatthew Knepley 19*649ef022SMatthew Knepley -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append 20*649ef022SMatthew Knepley F*/ 21*649ef022SMatthew Knepley 22*649ef022SMatthew Knepley #include <petscdmplex.h> 23*649ef022SMatthew Knepley #include <petscsnes.h> 24*649ef022SMatthew Knepley #include <petscts.h> 25*649ef022SMatthew Knepley #include <petscds.h> 26*649ef022SMatthew Knepley #include <petscbag.h> 27*649ef022SMatthew Knepley 28*649ef022SMatthew Knepley typedef enum {SOL_QUADRATIC, SOL_CUBIC, SOL_CUBIC_TRIG, NUM_SOL_TYPES} SolType; 29*649ef022SMatthew Knepley const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic", "cubic_trig", "unknown"}; 30*649ef022SMatthew Knepley 31*649ef022SMatthew Knepley typedef struct { 32*649ef022SMatthew Knepley PetscReal nu; /* Kinematic viscosity */ 33*649ef022SMatthew Knepley PetscReal alpha; /* Thermal diffusivity */ 34*649ef022SMatthew Knepley PetscReal T_in; /* Inlet temperature*/ 35*649ef022SMatthew Knepley } Parameter; 36*649ef022SMatthew Knepley 37*649ef022SMatthew Knepley typedef struct { 38*649ef022SMatthew Knepley /* Problem definition */ 39*649ef022SMatthew Knepley PetscBag bag; /* Holds problem parameters */ 40*649ef022SMatthew Knepley SolType solType; /* MMS solution type */ 41*649ef022SMatthew Knepley } AppCtx; 42*649ef022SMatthew Knepley 43*649ef022SMatthew Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 44*649ef022SMatthew Knepley { 45*649ef022SMatthew Knepley PetscInt d; 46*649ef022SMatthew Knepley for (d = 0; d < Nc; ++d) u[d] = 0.0; 47*649ef022SMatthew Knepley return 0; 48*649ef022SMatthew Knepley } 49*649ef022SMatthew Knepley 50*649ef022SMatthew Knepley static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 51*649ef022SMatthew Knepley { 52*649ef022SMatthew Knepley PetscInt d; 53*649ef022SMatthew Knepley for (d = 0; d < Nc; ++d) u[d] = 1.0; 54*649ef022SMatthew Knepley return 0; 55*649ef022SMatthew Knepley } 56*649ef022SMatthew Knepley 57*649ef022SMatthew Knepley /* 58*649ef022SMatthew Knepley CASE: quadratic 59*649ef022SMatthew Knepley In 2D we use exact solution: 60*649ef022SMatthew Knepley 61*649ef022SMatthew Knepley u = t + x^2 + y^2 62*649ef022SMatthew Knepley v = t + 2x^2 - 2xy 63*649ef022SMatthew Knepley p = x + y - 1 64*649ef022SMatthew Knepley T = t + x + y 65*649ef022SMatthew Knepley f = <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2 -4\nu + 2, t (2x - 2y) + 4xy^2 + 2x^2y - 2y^3 -4\nu + 2> 66*649ef022SMatthew Knepley Q = 1 + 2t + 3x^2 - 2xy + y^2 67*649ef022SMatthew Knepley 68*649ef022SMatthew Knepley so that 69*649ef022SMatthew Knepley 70*649ef022SMatthew Knepley \nabla \cdot u = 2x - 2x = 0 71*649ef022SMatthew Knepley 72*649ef022SMatthew Knepley f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p 73*649ef022SMatthew Knepley = <1, 1> + <t + x^2 + y^2, t + 2x^2 - 2xy> . <<2x, 4x - 2y>, <2y, -2x>> - \nu <4, 4> + <1, 1> 74*649ef022SMatthew Knepley = <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2, t (2x - 2y) + 2x^2y + 4xy^2 - 2y^3> + <-4 \nu + 2, -4\nu + 2> 75*649ef022SMatthew Knepley = <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2 - 4\nu + 2, t (2x - 2y) + 4xy^2 + 2x^2y - 2y^3 - 4\nu + 2> 76*649ef022SMatthew Knepley 77*649ef022SMatthew Knepley Q = dT/dt + u \cdot \nabla T - \alpha \Delta T 78*649ef022SMatthew Knepley = 1 + <t + x^2 + y^2, t + 2x^2 - 2xy> . <1, 1> - \alpha 0 79*649ef022SMatthew Knepley = 1 + 2t + 3x^2 - 2xy + y^2 80*649ef022SMatthew Knepley */ 81*649ef022SMatthew Knepley 82*649ef022SMatthew Knepley static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 83*649ef022SMatthew Knepley { 84*649ef022SMatthew Knepley u[0] = time + X[0]*X[0] + X[1]*X[1]; 85*649ef022SMatthew Knepley u[1] = time + 2.0*X[0]*X[0] - 2.0*X[0]*X[1]; 86*649ef022SMatthew Knepley return 0; 87*649ef022SMatthew Knepley } 88*649ef022SMatthew Knepley static PetscErrorCode quadratic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 89*649ef022SMatthew Knepley { 90*649ef022SMatthew Knepley u[0] = 1.0; 91*649ef022SMatthew Knepley u[1] = 1.0; 92*649ef022SMatthew Knepley return 0; 93*649ef022SMatthew Knepley } 94*649ef022SMatthew Knepley 95*649ef022SMatthew Knepley static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 96*649ef022SMatthew Knepley { 97*649ef022SMatthew Knepley p[0] = X[0] + X[1] - 1.0; 98*649ef022SMatthew Knepley return 0; 99*649ef022SMatthew Knepley } 100*649ef022SMatthew Knepley 101*649ef022SMatthew Knepley static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 102*649ef022SMatthew Knepley { 103*649ef022SMatthew Knepley T[0] = time + X[0] + X[1]; 104*649ef022SMatthew Knepley return 0; 105*649ef022SMatthew Knepley } 106*649ef022SMatthew Knepley static PetscErrorCode quadratic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 107*649ef022SMatthew Knepley { 108*649ef022SMatthew Knepley T[0] = 1.0; 109*649ef022SMatthew Knepley return 0; 110*649ef022SMatthew Knepley } 111*649ef022SMatthew Knepley 112*649ef022SMatthew Knepley /* f0_v = du/dt - f */ 113*649ef022SMatthew Knepley static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 114*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 115*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 116*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 117*649ef022SMatthew Knepley { 118*649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 119*649ef022SMatthew Knepley PetscInt Nc = dim; 120*649ef022SMatthew Knepley PetscInt c, d; 121*649ef022SMatthew Knepley 122*649ef022SMatthew Knepley for (d = 0; d<dim; ++d) f0[d] = u_t[uOff[0]+d]; 123*649ef022SMatthew Knepley 124*649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 125*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 126*649ef022SMatthew Knepley } 127*649ef022SMatthew Knepley f0[0] -= (t*(2*X[0] + 2*X[1]) + 2*X[0]*X[0]*X[0] + 4*X[0]*X[0]*X[1] - 2*X[0]*X[1]*X[1] - 4.0*nu + 2); 128*649ef022SMatthew Knepley f0[1] -= (t*(2*X[0] - 2*X[1]) + 4*X[0]*X[1]*X[1] + 2*X[0]*X[0]*X[1] - 2*X[1]*X[1]*X[1] - 4.0*nu + 2); 129*649ef022SMatthew Knepley } 130*649ef022SMatthew Knepley 131*649ef022SMatthew Knepley /* f0_w = dT/dt + u.grad(T) - Q */ 132*649ef022SMatthew Knepley static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 133*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 134*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 135*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 136*649ef022SMatthew Knepley { 137*649ef022SMatthew Knepley PetscInt d; 138*649ef022SMatthew Knepley f0[0] = 0; 139*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 140*649ef022SMatthew Knepley f0[0] += u_t[uOff[2]] - (2*t + 1 + 3*X[0]*X[0] - 2*X[0]*X[1] + X[1]*X[1]); 141*649ef022SMatthew Knepley } 142*649ef022SMatthew Knepley 143*649ef022SMatthew Knepley /* 144*649ef022SMatthew Knepley CASE: cubic 145*649ef022SMatthew Knepley In 2D we use exact solution: 146*649ef022SMatthew Knepley 147*649ef022SMatthew Knepley u = t + x^3 + y^3 148*649ef022SMatthew Knepley v = t + 2x^3 - 3x^2y 149*649ef022SMatthew Knepley p = 3/2 x^2 + 3/2 y^2 - 1 150*649ef022SMatthew Knepley T = t + 1/2 x^2 + 1/2 y^2 151*649ef022SMatthew Knepley f = < t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x + 1, 152*649ef022SMatthew Knepley t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y + 1> 153*649ef022SMatthew Knepley Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 + xt + yt - 2\alpha + 1 154*649ef022SMatthew Knepley 155*649ef022SMatthew Knepley so that 156*649ef022SMatthew Knepley 157*649ef022SMatthew Knepley \nabla \cdot u = 3x^2 - 3x^2 = 0 158*649ef022SMatthew Knepley 159*649ef022SMatthew Knepley du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p - f 160*649ef022SMatthew Knepley = <1,1> + <t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3, t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4> - \nu<6x + 6y, 12x - 6y> + <3x, 3y> - <t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x + 1, t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y + 1> = 0 161*649ef022SMatthew Knepley 162*649ef022SMatthew Knepley dT/dt + u \cdot \nabla T - \alpha \Delta T - Q = 1 + (x^3 + y^3) x + (2x^3 - 3x^2y) y - 2*\alpha - (x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2*\alpha +1) = 0 163*649ef022SMatthew Knepley */ 164*649ef022SMatthew Knepley static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 165*649ef022SMatthew Knepley { 166*649ef022SMatthew Knepley u[0] = time + X[0]*X[0]*X[0] + X[1]*X[1]*X[1]; 167*649ef022SMatthew Knepley u[1] = time + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1]; 168*649ef022SMatthew Knepley return 0; 169*649ef022SMatthew Knepley } 170*649ef022SMatthew Knepley static PetscErrorCode cubic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 171*649ef022SMatthew Knepley { 172*649ef022SMatthew Knepley u[0] = 1.0; 173*649ef022SMatthew Knepley u[1] = 1.0; 174*649ef022SMatthew Knepley return 0; 175*649ef022SMatthew Knepley } 176*649ef022SMatthew Knepley 177*649ef022SMatthew Knepley static PetscErrorCode cubic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 178*649ef022SMatthew Knepley { 179*649ef022SMatthew Knepley p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0; 180*649ef022SMatthew Knepley return 0; 181*649ef022SMatthew Knepley } 182*649ef022SMatthew Knepley 183*649ef022SMatthew Knepley static PetscErrorCode cubic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 184*649ef022SMatthew Knepley { 185*649ef022SMatthew Knepley T[0] = time + X[0]*X[0]/2.0 + X[1]*X[1]/2.0; 186*649ef022SMatthew Knepley return 0; 187*649ef022SMatthew Knepley } 188*649ef022SMatthew Knepley static PetscErrorCode cubic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 189*649ef022SMatthew Knepley { 190*649ef022SMatthew Knepley T[0] = 1.0; 191*649ef022SMatthew Knepley return 0; 192*649ef022SMatthew Knepley } 193*649ef022SMatthew Knepley 194*649ef022SMatthew Knepley 195*649ef022SMatthew Knepley static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 196*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 197*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 198*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 199*649ef022SMatthew Knepley { 200*649ef022SMatthew Knepley PetscInt c, d; 201*649ef022SMatthew Knepley PetscInt Nc = dim; 202*649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 203*649ef022SMatthew Knepley 204*649ef022SMatthew Knepley for (d=0; d<dim; ++d) f0[d] = u_t[uOff[0]+d]; 205*649ef022SMatthew Knepley 206*649ef022SMatthew Knepley for (c=0; c<Nc; ++c) { 207*649ef022SMatthew Knepley for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 208*649ef022SMatthew Knepley } 209*649ef022SMatthew Knepley f0[0] -= (t*(3*X[0]*X[0] + 3*X[1]*X[1]) + 3*X[0]*X[0]*X[0]*X[0]*X[0] + 6*X[0]*X[0]*X[0]*X[1]*X[1] - 6*X[0]*X[0]*X[1]*X[1]*X[1] - ( 6*X[0] + 6*X[1])*nu + 3*X[0] + 1); 210*649ef022SMatthew Knepley f0[1] -= (t*(3*X[0]*X[0] - 6*X[0]*X[1]) + 3*X[0]*X[0]*X[0]*X[0]*X[1] + 6*X[0]*X[0]*X[1]*X[1]*X[1] - 6*X[0]*X[1]*X[1]*X[1]*X[1] - (12*X[0] - 6*X[1])*nu + 3*X[1] + 1); 211*649ef022SMatthew Knepley } 212*649ef022SMatthew Knepley 213*649ef022SMatthew Knepley static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 214*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 215*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 216*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 217*649ef022SMatthew Knepley { 218*649ef022SMatthew Knepley PetscInt d; 219*649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 220*649ef022SMatthew Knepley 221*649ef022SMatthew Knepley for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 222*649ef022SMatthew Knepley f0[0] += u_t[uOff[2]] - (X[0]*X[0]*X[0]*X[0] + 2.0*X[0]*X[0]*X[0]*X[1] - 3.0*X[0]*X[0]*X[1]*X[1] + X[0]*X[1]*X[1]*X[1] + X[0]*t + X[1]*t - 2.0*alpha + 1); 223*649ef022SMatthew Knepley } 224*649ef022SMatthew Knepley 225*649ef022SMatthew Knepley /* 226*649ef022SMatthew Knepley CASE: cubic-trigonometric 227*649ef022SMatthew Knepley In 2D we use exact solution: 228*649ef022SMatthew Knepley 229*649ef022SMatthew Knepley u = beta cos t + x^3 + y^3 230*649ef022SMatthew Knepley v = beta sin t + 2x^3 - 3x^2y 231*649ef022SMatthew Knepley p = 3/2 x^2 + 3/2 y^2 - 1 232*649ef022SMatthew Knepley T = 20 cos t + 1/2 x^2 + 1/2 y^2 233*649ef022SMatthew Knepley f = < beta cos t 3x^2 + beta sin t (3y^2 - 1) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x, 234*649ef022SMatthew Knepley beta cos t (6x^2 - 6xy) - beta sin t (3x^2) + 3x^4y + 6x^2y^3 - 6xy^4 - \nu(12x - 6y) + 3y> 235*649ef022SMatthew Knepley Q = beta cos t x + beta sin t (y - 1) + x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2\alpha 236*649ef022SMatthew Knepley 237*649ef022SMatthew Knepley so that 238*649ef022SMatthew Knepley 239*649ef022SMatthew Knepley \nabla \cdot u = 3x^2 - 3x^2 = 0 240*649ef022SMatthew Knepley 241*649ef022SMatthew Knepley f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p 242*649ef022SMatthew Knepley = <-sin t, cos t> + <cos t + x^3 + y^3, sin t + 2x^3 - 3x^2y> <<3x^2, 6x^2 - 6xy>, <3y^2, -3x^2>> - \nu <6x + 6y, 12x - 6y> + <3x, 3y> 243*649ef022SMatthew Knepley = <-sin t, cos t> + <cos t 3x^2 + 3x^5 + 3x^2y^3 + sin t 3y^2 + 6x^3y^2 - 9x^2y^3, cos t (6x^2 - 6xy) + 6x^5 - 6x^4y + 6x^2y^3 - 6xy^4 + sin t (-3x^2) - 6x^5 + 9x^4y> - \nu <6x + 6y, 12x - 6y> + <3x, 3y> 244*649ef022SMatthew Knepley = <cos t (3x^2) + sin t (3y^2 - 1) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu (6x + 6y) + 3x, 245*649ef022SMatthew Knepley cos t (6x^2 - 6xy) - sin t (3x^2) + 3x^4y + 6x^2y^3 - 6xy^4 - \nu (12x - 6y) + 3y> 246*649ef022SMatthew Knepley 247*649ef022SMatthew Knepley Q = dT/dt + u \cdot \nabla T - \alpha \Delta T 248*649ef022SMatthew Knepley = -sin t + <cos t + x^3 + y^3, sin t + 2x^3 - 3x^2y> . <x, y> - 2 \alpha 249*649ef022SMatthew Knepley = -sin t + cos t (x) + x^4 + xy^3 + sin t (y) + 2x^3y - 3x^2y^2 - 2 \alpha 250*649ef022SMatthew Knepley = cos t x + sin t (y - 1) + (x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2 \alpha) 251*649ef022SMatthew Knepley */ 252*649ef022SMatthew Knepley static PetscErrorCode cubic_trig_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 253*649ef022SMatthew Knepley { 254*649ef022SMatthew Knepley u[0] = 100.*PetscCosReal(time) + X[0]*X[0]*X[0] + X[1]*X[1]*X[1]; 255*649ef022SMatthew Knepley u[1] = 100.*PetscSinReal(time) + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1]; 256*649ef022SMatthew Knepley return 0; 257*649ef022SMatthew Knepley } 258*649ef022SMatthew Knepley static PetscErrorCode cubic_trig_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 259*649ef022SMatthew Knepley { 260*649ef022SMatthew Knepley u[0] = -100.*PetscSinReal(time); 261*649ef022SMatthew Knepley u[1] = 100.*PetscCosReal(time); 262*649ef022SMatthew Knepley return 0; 263*649ef022SMatthew Knepley } 264*649ef022SMatthew Knepley 265*649ef022SMatthew Knepley static PetscErrorCode cubic_trig_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 266*649ef022SMatthew Knepley { 267*649ef022SMatthew Knepley p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0; 268*649ef022SMatthew Knepley return 0; 269*649ef022SMatthew Knepley } 270*649ef022SMatthew Knepley 271*649ef022SMatthew Knepley static PetscErrorCode cubic_trig_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 272*649ef022SMatthew Knepley { 273*649ef022SMatthew Knepley T[0] = 100.*PetscCosReal(time) + X[0]*X[0]/2.0 + X[1]*X[1]/2.0; 274*649ef022SMatthew Knepley return 0; 275*649ef022SMatthew Knepley } 276*649ef022SMatthew Knepley static PetscErrorCode cubic_trig_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 277*649ef022SMatthew Knepley { 278*649ef022SMatthew Knepley T[0] = -100.*PetscSinReal(time); 279*649ef022SMatthew Knepley return 0; 280*649ef022SMatthew Knepley } 281*649ef022SMatthew Knepley 282*649ef022SMatthew Knepley static void f0_cubic_trig_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 283*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 284*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 285*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 286*649ef022SMatthew Knepley { 287*649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 288*649ef022SMatthew Knepley PetscInt Nc = dim; 289*649ef022SMatthew Knepley PetscInt c, d; 290*649ef022SMatthew Knepley 291*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f0[d] = u_t[uOff[0]+d]; 292*649ef022SMatthew Knepley 293*649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 294*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 295*649ef022SMatthew Knepley } 296*649ef022SMatthew Knepley f0[0] -= 100.*PetscCosReal(t)*(3*X[0]*X[0]) + 100.*PetscSinReal(t)*(3*X[1]*X[1] - 1.) + 3*X[0]*X[0]*X[0]*X[0]*X[0] + 6*X[0]*X[0]*X[0]*X[1]*X[1] - 6*X[0]*X[0]*X[1]*X[1]*X[1] - ( 6*X[0] + 6*X[1])*nu + 3*X[0]; 297*649ef022SMatthew Knepley f0[1] -= 100.*PetscCosReal(t)*(6*X[0]*X[0] - 6*X[0]*X[1]) - 100.*PetscSinReal(t)*(3*X[0]*X[0]) + 3*X[0]*X[0]*X[0]*X[0]*X[1] + 6*X[0]*X[0]*X[1]*X[1]*X[1] - 6*X[0]*X[1]*X[1]*X[1]*X[1] - (12*X[0] - 6*X[1])*nu + 3*X[1]; 298*649ef022SMatthew Knepley } 299*649ef022SMatthew Knepley 300*649ef022SMatthew Knepley static void f0_cubic_trig_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 301*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 302*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 303*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 304*649ef022SMatthew Knepley { 305*649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 306*649ef022SMatthew Knepley PetscInt d; 307*649ef022SMatthew Knepley 308*649ef022SMatthew Knepley for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 309*649ef022SMatthew Knepley f0[0] += u_t[uOff[2]] - (100.*PetscCosReal(t)*X[0] + 100.*PetscSinReal(t)*(X[1] - 1.) + X[0]*X[0]*X[0]*X[0] + 2.0*X[0]*X[0]*X[0]*X[1] - 3.0*X[0]*X[0]*X[1]*X[1] + X[0]*X[1]*X[1]*X[1] - 2.0*alpha); 310*649ef022SMatthew Knepley } 311*649ef022SMatthew Knepley 312*649ef022SMatthew Knepley static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 313*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 314*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 315*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 316*649ef022SMatthew Knepley { 317*649ef022SMatthew Knepley PetscInt d; 318*649ef022SMatthew Knepley for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; 319*649ef022SMatthew Knepley } 320*649ef022SMatthew Knepley 321*649ef022SMatthew Knepley /*f1_v = \nu[grad(u) + grad(u)^T] - pI */ 322*649ef022SMatthew Knepley static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 323*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 324*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 325*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 326*649ef022SMatthew Knepley { 327*649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 328*649ef022SMatthew Knepley const PetscInt Nc = dim; 329*649ef022SMatthew Knepley PetscInt c, d; 330*649ef022SMatthew Knepley 331*649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 332*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 333*649ef022SMatthew Knepley f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]); 334*649ef022SMatthew Knepley //f1[c*dim+d] = nu*u_x[c*dim+d]; 335*649ef022SMatthew Knepley } 336*649ef022SMatthew Knepley f1[c*dim+c] -= u[uOff[1]]; 337*649ef022SMatthew Knepley } 338*649ef022SMatthew Knepley } 339*649ef022SMatthew Knepley 340*649ef022SMatthew Knepley static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 341*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 342*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 343*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 344*649ef022SMatthew Knepley { 345*649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 346*649ef022SMatthew Knepley PetscInt d; 347*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d]; 348*649ef022SMatthew Knepley } 349*649ef022SMatthew Knepley 350*649ef022SMatthew Knepley /*Jacobians*/ 351*649ef022SMatthew Knepley static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 352*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 353*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 354*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 355*649ef022SMatthew Knepley { 356*649ef022SMatthew Knepley PetscInt d; 357*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; 358*649ef022SMatthew Knepley } 359*649ef022SMatthew Knepley 360*649ef022SMatthew Knepley static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 361*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 362*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 363*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 364*649ef022SMatthew Knepley { 365*649ef022SMatthew Knepley PetscInt c, d; 366*649ef022SMatthew Knepley const PetscInt Nc = dim; 367*649ef022SMatthew Knepley 368*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g0[d*dim+d] = u_tShift; 369*649ef022SMatthew Knepley 370*649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 371*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 372*649ef022SMatthew Knepley g0[c*Nc+d] += u_x[ c*Nc+d]; 373*649ef022SMatthew Knepley } 374*649ef022SMatthew Knepley } 375*649ef022SMatthew Knepley } 376*649ef022SMatthew Knepley 377*649ef022SMatthew Knepley static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 378*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 379*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 380*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 381*649ef022SMatthew Knepley { 382*649ef022SMatthew Knepley PetscInt NcI = dim; 383*649ef022SMatthew Knepley PetscInt NcJ = dim; 384*649ef022SMatthew Knepley PetscInt c, d, e; 385*649ef022SMatthew Knepley 386*649ef022SMatthew Knepley for (c = 0; c < NcI; ++c) { 387*649ef022SMatthew Knepley for (d = 0; d < NcJ; ++d) { 388*649ef022SMatthew Knepley for (e = 0; e < dim; ++e) { 389*649ef022SMatthew Knepley if (c == d) { 390*649ef022SMatthew Knepley g1[(c*NcJ+d)*dim+e] += u[e]; 391*649ef022SMatthew Knepley } 392*649ef022SMatthew Knepley } 393*649ef022SMatthew Knepley } 394*649ef022SMatthew Knepley } 395*649ef022SMatthew Knepley } 396*649ef022SMatthew Knepley 397*649ef022SMatthew Knepley 398*649ef022SMatthew Knepley static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 399*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 400*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 401*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 402*649ef022SMatthew Knepley { 403*649ef022SMatthew Knepley PetscInt d; 404*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; 405*649ef022SMatthew Knepley } 406*649ef022SMatthew Knepley 407*649ef022SMatthew Knepley static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 408*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 409*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 410*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 411*649ef022SMatthew Knepley { 412*649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 413*649ef022SMatthew Knepley const PetscInt Nc = dim; 414*649ef022SMatthew Knepley PetscInt c, d; 415*649ef022SMatthew Knepley 416*649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 417*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 418*649ef022SMatthew Knepley g3[((c*Nc+c)*dim+d)*dim+d] += nu; // gradU 419*649ef022SMatthew Knepley g3[((c*Nc+d)*dim+d)*dim+c] += nu; // gradU transpose 420*649ef022SMatthew Knepley } 421*649ef022SMatthew Knepley } 422*649ef022SMatthew Knepley } 423*649ef022SMatthew Knepley 424*649ef022SMatthew Knepley static void g0_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 425*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 426*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 427*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 428*649ef022SMatthew Knepley { 429*649ef022SMatthew Knepley PetscInt d; 430*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g0[d] = u_tShift; 431*649ef022SMatthew Knepley } 432*649ef022SMatthew Knepley 433*649ef022SMatthew Knepley static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 434*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 435*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 436*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 437*649ef022SMatthew Knepley { 438*649ef022SMatthew Knepley PetscInt d; 439*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d]; 440*649ef022SMatthew Knepley } 441*649ef022SMatthew Knepley 442*649ef022SMatthew Knepley static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 443*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 444*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 445*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 446*649ef022SMatthew Knepley { 447*649ef022SMatthew Knepley PetscInt d; 448*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d]; 449*649ef022SMatthew Knepley } 450*649ef022SMatthew Knepley 451*649ef022SMatthew Knepley static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 452*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 453*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 454*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 455*649ef022SMatthew Knepley { 456*649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 457*649ef022SMatthew Knepley PetscInt d; 458*649ef022SMatthew Knepley 459*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha; 460*649ef022SMatthew Knepley } 461*649ef022SMatthew Knepley 462*649ef022SMatthew Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 463*649ef022SMatthew Knepley { 464*649ef022SMatthew Knepley PetscInt sol; 465*649ef022SMatthew Knepley PetscErrorCode ierr; 466*649ef022SMatthew Knepley 467*649ef022SMatthew Knepley 468*649ef022SMatthew Knepley PetscFunctionBeginUser; 469*649ef022SMatthew Knepley options->solType = SOL_QUADRATIC; 470*649ef022SMatthew Knepley 471*649ef022SMatthew Knepley ierr = PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX");CHKERRQ(ierr); 472*649ef022SMatthew Knepley sol = options->solType; 473*649ef022SMatthew Knepley ierr = PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr); 474*649ef022SMatthew Knepley options->solType = (SolType) sol; 475*649ef022SMatthew Knepley ierr = PetscOptionsEnd(); 476*649ef022SMatthew Knepley PetscFunctionReturn(0); 477*649ef022SMatthew Knepley } 478*649ef022SMatthew Knepley 479*649ef022SMatthew Knepley static PetscErrorCode SetupParameters(AppCtx *user) 480*649ef022SMatthew Knepley { 481*649ef022SMatthew Knepley PetscBag bag; 482*649ef022SMatthew Knepley Parameter *p; 483*649ef022SMatthew Knepley PetscErrorCode ierr; 484*649ef022SMatthew Knepley 485*649ef022SMatthew Knepley PetscFunctionBeginUser; 486*649ef022SMatthew Knepley /* setup PETSc parameter bag */ 487*649ef022SMatthew Knepley ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr); 488*649ef022SMatthew Knepley ierr = PetscBagSetName(user->bag, "par", "Low Mach flow parameters");CHKERRQ(ierr); 489*649ef022SMatthew Knepley bag = user->bag; 490*649ef022SMatthew Knepley ierr = PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity");CHKERRQ(ierr); 491*649ef022SMatthew Knepley ierr = PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity");CHKERRQ(ierr); 492*649ef022SMatthew Knepley ierr = PetscBagRegisterReal(bag, &p->T_in, 1.0, "T_in", "Inlet temperature");CHKERRQ(ierr); 493*649ef022SMatthew Knepley PetscFunctionReturn(0); 494*649ef022SMatthew Knepley } 495*649ef022SMatthew Knepley 496*649ef022SMatthew Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 497*649ef022SMatthew Knepley { 498*649ef022SMatthew Knepley PetscErrorCode ierr; 499*649ef022SMatthew Knepley 500*649ef022SMatthew Knepley PetscFunctionBeginUser; 501*649ef022SMatthew Knepley ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 502*649ef022SMatthew Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 503*649ef022SMatthew Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 504*649ef022SMatthew Knepley PetscFunctionReturn(0); 505*649ef022SMatthew Knepley } 506*649ef022SMatthew Knepley 507*649ef022SMatthew Knepley static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 508*649ef022SMatthew Knepley { 509*649ef022SMatthew Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 510*649ef022SMatthew Knepley PetscErrorCode (*exactFuncs_t[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 511*649ef022SMatthew Knepley PetscDS prob; 512*649ef022SMatthew Knepley Parameter *ctx; 513*649ef022SMatthew Knepley PetscInt id; 514*649ef022SMatthew Knepley PetscErrorCode ierr; 515*649ef022SMatthew Knepley 516*649ef022SMatthew Knepley PetscFunctionBeginUser; 517*649ef022SMatthew Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 518*649ef022SMatthew Knepley switch(user->solType){ 519*649ef022SMatthew Knepley case SOL_QUADRATIC: 520*649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v);CHKERRQ(ierr); 521*649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w);CHKERRQ(ierr); 522*649ef022SMatthew Knepley 523*649ef022SMatthew Knepley exactFuncs[0] = quadratic_u; 524*649ef022SMatthew Knepley exactFuncs[1] = quadratic_p; 525*649ef022SMatthew Knepley exactFuncs[2] = quadratic_T; 526*649ef022SMatthew Knepley exactFuncs_t[0] = quadratic_u_t; 527*649ef022SMatthew Knepley exactFuncs_t[1] = NULL; 528*649ef022SMatthew Knepley exactFuncs_t[2] = quadratic_T_t; 529*649ef022SMatthew Knepley break; 530*649ef022SMatthew Knepley case SOL_CUBIC: 531*649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v);CHKERRQ(ierr); 532*649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w);CHKERRQ(ierr); 533*649ef022SMatthew Knepley 534*649ef022SMatthew Knepley exactFuncs[0] = cubic_u; 535*649ef022SMatthew Knepley exactFuncs[1] = cubic_p; 536*649ef022SMatthew Knepley exactFuncs[2] = cubic_T; 537*649ef022SMatthew Knepley exactFuncs_t[0] = cubic_u_t; 538*649ef022SMatthew Knepley exactFuncs_t[1] = NULL; 539*649ef022SMatthew Knepley exactFuncs_t[2] = cubic_T_t; 540*649ef022SMatthew Knepley break; 541*649ef022SMatthew Knepley case SOL_CUBIC_TRIG: 542*649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 0, f0_cubic_trig_v, f1_v);CHKERRQ(ierr); 543*649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 2, f0_cubic_trig_w, f1_w);CHKERRQ(ierr); 544*649ef022SMatthew Knepley 545*649ef022SMatthew Knepley exactFuncs[0] = cubic_trig_u; 546*649ef022SMatthew Knepley exactFuncs[1] = cubic_trig_p; 547*649ef022SMatthew Knepley exactFuncs[2] = cubic_trig_T; 548*649ef022SMatthew Knepley exactFuncs_t[0] = cubic_trig_u_t; 549*649ef022SMatthew Knepley exactFuncs_t[1] = NULL; 550*649ef022SMatthew Knepley exactFuncs_t[2] = cubic_trig_T_t; 551*649ef022SMatthew Knepley break; 552*649ef022SMatthew Knepley default: SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 553*649ef022SMatthew Knepley } 554*649ef022SMatthew Knepley 555*649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 1, f0_q, NULL);CHKERRQ(ierr); 556*649ef022SMatthew Knepley 557*649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu);CHKERRQ(ierr); 558*649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL);CHKERRQ(ierr); 559*649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL);CHKERRQ(ierr); 560*649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL);CHKERRQ(ierr); 561*649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 2, 2, g0_wT, g1_wT, NULL, g3_wT);CHKERRQ(ierr); 562*649ef022SMatthew Knepley /* Setup constants */ 563*649ef022SMatthew Knepley { 564*649ef022SMatthew Knepley Parameter *param; 565*649ef022SMatthew Knepley PetscScalar constants[3]; 566*649ef022SMatthew Knepley 567*649ef022SMatthew Knepley ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 568*649ef022SMatthew Knepley 569*649ef022SMatthew Knepley constants[0] = param->nu; 570*649ef022SMatthew Knepley constants[1] = param->alpha; 571*649ef022SMatthew Knepley constants[2] = param->T_in; 572*649ef022SMatthew Knepley ierr = PetscDSSetConstants(prob, 3, constants);CHKERRQ(ierr); 573*649ef022SMatthew Knepley } 574*649ef022SMatthew Knepley /* Setup Boundary Conditions */ 575*649ef022SMatthew Knepley ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr); 576*649ef022SMatthew Knepley id = 3; 577*649ef022SMatthew 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); 578*649ef022SMatthew Knepley id = 1; 579*649ef022SMatthew 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); 580*649ef022SMatthew Knepley id = 2; 581*649ef022SMatthew 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); 582*649ef022SMatthew Knepley id = 4; 583*649ef022SMatthew 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); 584*649ef022SMatthew Knepley id = 3; 585*649ef022SMatthew 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); 586*649ef022SMatthew Knepley id = 1; 587*649ef022SMatthew 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); 588*649ef022SMatthew Knepley id = 2; 589*649ef022SMatthew 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); 590*649ef022SMatthew Knepley id = 4; 591*649ef022SMatthew 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); 592*649ef022SMatthew Knepley 593*649ef022SMatthew Knepley /*setup exact solution.*/ 594*649ef022SMatthew Knepley ierr = PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx);CHKERRQ(ierr); 595*649ef022SMatthew Knepley ierr = PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx);CHKERRQ(ierr); 596*649ef022SMatthew Knepley ierr = PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx);CHKERRQ(ierr); 597*649ef022SMatthew Knepley ierr = PetscDSSetExactSolutionTimeDerivative(prob, 0, exactFuncs_t[0], ctx);CHKERRQ(ierr); 598*649ef022SMatthew Knepley ierr = PetscDSSetExactSolutionTimeDerivative(prob, 1, exactFuncs_t[1], ctx);CHKERRQ(ierr); 599*649ef022SMatthew Knepley ierr = PetscDSSetExactSolutionTimeDerivative(prob, 2, exactFuncs_t[2], ctx);CHKERRQ(ierr); 600*649ef022SMatthew Knepley PetscFunctionReturn(0); 601*649ef022SMatthew Knepley } 602*649ef022SMatthew Knepley 603*649ef022SMatthew Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 604*649ef022SMatthew Knepley { 605*649ef022SMatthew Knepley DM cdm = dm; 606*649ef022SMatthew Knepley PetscFE fe[3]; 607*649ef022SMatthew Knepley Parameter *param; 608*649ef022SMatthew Knepley MPI_Comm comm; 609*649ef022SMatthew Knepley DMPolytopeType ct; 610*649ef022SMatthew Knepley PetscInt dim, cStart; 611*649ef022SMatthew Knepley PetscBool simplex; 612*649ef022SMatthew Knepley PetscErrorCode ierr; 613*649ef022SMatthew Knepley 614*649ef022SMatthew Knepley PetscFunctionBeginUser; 615*649ef022SMatthew Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 616*649ef022SMatthew Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 617*649ef022SMatthew Knepley ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 618*649ef022SMatthew Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 619*649ef022SMatthew Knepley /* Create finite element */ 620*649ef022SMatthew Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 621*649ef022SMatthew Knepley ierr = PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 622*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr); 623*649ef022SMatthew Knepley 624*649ef022SMatthew Knepley ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 625*649ef022SMatthew Knepley ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 626*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr); 627*649ef022SMatthew Knepley 628*649ef022SMatthew Knepley ierr = PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]);CHKERRQ(ierr); 629*649ef022SMatthew Knepley ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr); 630*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) fe[2], "temperature");CHKERRQ(ierr); 631*649ef022SMatthew Knepley 632*649ef022SMatthew Knepley /* Set discretization and boundary conditions for each mesh */ 633*649ef022SMatthew Knepley ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 634*649ef022SMatthew Knepley ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr); 635*649ef022SMatthew Knepley ierr = DMSetField(dm, 2, NULL, (PetscObject) fe[2]);CHKERRQ(ierr); 636*649ef022SMatthew Knepley ierr = DMCreateDS(dm);CHKERRQ(ierr); 637*649ef022SMatthew Knepley ierr = SetupProblem(dm, user);CHKERRQ(ierr); 638*649ef022SMatthew Knepley ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 639*649ef022SMatthew Knepley while (cdm) { 640*649ef022SMatthew Knepley ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 641*649ef022SMatthew Knepley ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 642*649ef022SMatthew Knepley } 643*649ef022SMatthew Knepley ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); 644*649ef022SMatthew Knepley ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); 645*649ef022SMatthew Knepley ierr = PetscFEDestroy(&fe[2]);CHKERRQ(ierr); 646*649ef022SMatthew Knepley 647*649ef022SMatthew Knepley { 648*649ef022SMatthew Knepley PetscObject pressure; 649*649ef022SMatthew Knepley MatNullSpace nullspacePres; 650*649ef022SMatthew Knepley 651*649ef022SMatthew Knepley ierr = DMGetField(dm, 1, NULL, &pressure);CHKERRQ(ierr); 652*649ef022SMatthew Knepley ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres);CHKERRQ(ierr); 653*649ef022SMatthew Knepley ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres);CHKERRQ(ierr); 654*649ef022SMatthew Knepley ierr = MatNullSpaceDestroy(&nullspacePres);CHKERRQ(ierr); 655*649ef022SMatthew Knepley } 656*649ef022SMatthew Knepley 657*649ef022SMatthew Knepley PetscFunctionReturn(0); 658*649ef022SMatthew Knepley } 659*649ef022SMatthew Knepley 660*649ef022SMatthew Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) 661*649ef022SMatthew Knepley { 662*649ef022SMatthew Knepley Vec vec; 663*649ef022SMatthew Knepley PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 664*649ef022SMatthew Knepley PetscErrorCode ierr; 665*649ef022SMatthew Knepley 666*649ef022SMatthew Knepley PetscFunctionBeginUser; 667*649ef022SMatthew Knepley if (ofield != 1) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield); 668*649ef022SMatthew Knepley funcs[nfield] = constant; 669*649ef022SMatthew Knepley ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr); 670*649ef022SMatthew Knepley ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr); 671*649ef022SMatthew Knepley ierr = VecNormalize(vec, NULL);CHKERRQ(ierr); 672*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr); 673*649ef022SMatthew Knepley ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr); 674*649ef022SMatthew Knepley ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace);CHKERRQ(ierr); 675*649ef022SMatthew Knepley ierr = VecDestroy(&vec);CHKERRQ(ierr); 676*649ef022SMatthew Knepley PetscFunctionReturn(0); 677*649ef022SMatthew Knepley } 678*649ef022SMatthew Knepley 679*649ef022SMatthew Knepley static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u) 680*649ef022SMatthew Knepley { 681*649ef022SMatthew Knepley DM dm; 682*649ef022SMatthew Knepley MatNullSpace nullsp; 683*649ef022SMatthew Knepley PetscErrorCode ierr; 684*649ef022SMatthew Knepley 685*649ef022SMatthew Knepley PetscFunctionBegin; 686*649ef022SMatthew Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 687*649ef022SMatthew Knepley ierr = CreatePressureNullSpace(dm, 1, 1, &nullsp);CHKERRQ(ierr); 688*649ef022SMatthew Knepley ierr = MatNullSpaceRemove(nullsp, u);CHKERRQ(ierr); 689*649ef022SMatthew Knepley ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 690*649ef022SMatthew Knepley PetscFunctionReturn(0); 691*649ef022SMatthew Knepley } 692*649ef022SMatthew Knepley 693*649ef022SMatthew Knepley /* Make the discrete pressure discretely divergence free */ 694*649ef022SMatthew Knepley static PetscErrorCode RemoveDiscretePressureNullspace(TS ts) 695*649ef022SMatthew Knepley { 696*649ef022SMatthew Knepley Vec u; 697*649ef022SMatthew Knepley PetscErrorCode ierr; 698*649ef022SMatthew Knepley 699*649ef022SMatthew Knepley PetscFunctionBegin; 700*649ef022SMatthew Knepley ierr = TSGetSolution(ts, &u);CHKERRQ(ierr); 701*649ef022SMatthew Knepley ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr); 702*649ef022SMatthew Knepley PetscFunctionReturn(0); 703*649ef022SMatthew Knepley } 704*649ef022SMatthew Knepley 705*649ef022SMatthew Knepley static PetscErrorCode SetInitialConditions(TS ts, Vec u) 706*649ef022SMatthew Knepley { 707*649ef022SMatthew Knepley DM dm; 708*649ef022SMatthew Knepley PetscReal t; 709*649ef022SMatthew Knepley PetscErrorCode ierr; 710*649ef022SMatthew Knepley 711*649ef022SMatthew Knepley PetscFunctionBegin; 712*649ef022SMatthew Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 713*649ef022SMatthew Knepley ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 714*649ef022SMatthew Knepley ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr); 715*649ef022SMatthew Knepley ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr); 716*649ef022SMatthew Knepley PetscFunctionReturn(0); 717*649ef022SMatthew Knepley } 718*649ef022SMatthew Knepley 719*649ef022SMatthew Knepley static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) 720*649ef022SMatthew Knepley { 721*649ef022SMatthew Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 722*649ef022SMatthew Knepley void *ctxs[3]; 723*649ef022SMatthew Knepley DM dm; 724*649ef022SMatthew Knepley PetscDS ds; 725*649ef022SMatthew Knepley Vec v; 726*649ef022SMatthew Knepley PetscReal ferrors[3]; 727*649ef022SMatthew Knepley PetscInt f; 728*649ef022SMatthew Knepley PetscErrorCode ierr; 729*649ef022SMatthew Knepley 730*649ef022SMatthew Knepley PetscFunctionBeginUser; 731*649ef022SMatthew Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 732*649ef022SMatthew Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 733*649ef022SMatthew Knepley 734*649ef022SMatthew Knepley for (f = 0; f < 3; ++f) {ierr = PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]);CHKERRQ(ierr);} 735*649ef022SMatthew Knepley ierr = DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors);CHKERRQ(ierr); 736*649ef022SMatthew 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); 737*649ef022SMatthew Knepley 738*649ef022SMatthew Knepley ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr); 739*649ef022SMatthew Knepley //ierr = TSGetSolution(ts, &u);CHKERRQ(ierr); 740*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr); 741*649ef022SMatthew Knepley ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 742*649ef022SMatthew Knepley ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr); 743*649ef022SMatthew Knepley 744*649ef022SMatthew Knepley ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr); 745*649ef022SMatthew Knepley // ierr = VecSet(v, 0.0);CHKERRQ(ierr); 746*649ef022SMatthew Knepley ierr = DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, v);CHKERRQ(ierr); 747*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) v, "Exact Solution");CHKERRQ(ierr); 748*649ef022SMatthew Knepley ierr = VecViewFromOptions(v, NULL, "-exact_vec_view");CHKERRQ(ierr); 749*649ef022SMatthew Knepley ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr); 750*649ef022SMatthew Knepley 751*649ef022SMatthew Knepley PetscFunctionReturn(0); 752*649ef022SMatthew Knepley } 753*649ef022SMatthew Knepley 754*649ef022SMatthew Knepley int main(int argc, char **argv) 755*649ef022SMatthew Knepley { 756*649ef022SMatthew Knepley DM dm; /* problem definition */ 757*649ef022SMatthew Knepley TS ts; /* timestepper */ 758*649ef022SMatthew Knepley Vec u; /* solution */ 759*649ef022SMatthew Knepley AppCtx user; /* user-defined work context */ 760*649ef022SMatthew Knepley PetscReal t; 761*649ef022SMatthew Knepley PetscErrorCode ierr; 762*649ef022SMatthew Knepley 763*649ef022SMatthew Knepley ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 764*649ef022SMatthew Knepley ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 765*649ef022SMatthew Knepley ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr); 766*649ef022SMatthew Knepley ierr = SetupParameters(&user);CHKERRQ(ierr); 767*649ef022SMatthew Knepley ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 768*649ef022SMatthew Knepley ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 769*649ef022SMatthew Knepley ierr = TSSetDM(ts, dm);CHKERRQ(ierr); 770*649ef022SMatthew Knepley ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr); 771*649ef022SMatthew Knepley /* Setup problem */ 772*649ef022SMatthew Knepley ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr); 773*649ef022SMatthew Knepley ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr); 774*649ef022SMatthew Knepley 775*649ef022SMatthew Knepley ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 776*649ef022SMatthew Knepley ierr = DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace);CHKERRQ(ierr); 777*649ef022SMatthew Knepley 778*649ef022SMatthew Knepley ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr); 779*649ef022SMatthew Knepley ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr); 780*649ef022SMatthew Knepley ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr); 781*649ef022SMatthew Knepley ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 782*649ef022SMatthew Knepley ierr = TSSetPreStep(ts, RemoveDiscretePressureNullspace);CHKERRQ(ierr); 783*649ef022SMatthew Knepley ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 784*649ef022SMatthew Knepley 785*649ef022SMatthew Knepley ierr = TSSetComputeInitialCondition(ts, SetInitialConditions);CHKERRQ(ierr); /* Must come after SetFromOptions() */ 786*649ef022SMatthew Knepley ierr = SetInitialConditions(ts, u);CHKERRQ(ierr); 787*649ef022SMatthew Knepley ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 788*649ef022SMatthew Knepley ierr = DMSetOutputSequenceNumber(dm, 0, t);CHKERRQ(ierr); 789*649ef022SMatthew Knepley ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 790*649ef022SMatthew Knepley ierr = TSMonitorSet(ts, MonitorError, &user, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 791*649ef022SMatthew Knepley 792*649ef022SMatthew Knepley ierr = TSSolve(ts, u);CHKERRQ(ierr); 793*649ef022SMatthew Knepley ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 794*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr); 795*649ef022SMatthew Knepley 796*649ef022SMatthew Knepley ierr = VecDestroy(&u);CHKERRQ(ierr); 797*649ef022SMatthew Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 798*649ef022SMatthew Knepley ierr = TSDestroy(&ts);CHKERRQ(ierr); 799*649ef022SMatthew Knepley ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr); 800*649ef022SMatthew Knepley ierr = PetscFinalize(); 801*649ef022SMatthew Knepley return ierr; 802*649ef022SMatthew Knepley } 803*649ef022SMatthew Knepley 804*649ef022SMatthew Knepley /*TEST 805*649ef022SMatthew Knepley 806*649ef022SMatthew Knepley test: 807*649ef022SMatthew Knepley suffix: 2d_tri_p2_p1_p1 808*649ef022SMatthew Knepley requires: triangle !single 809*649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \ 810*649ef022SMatthew Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 811*649ef022SMatthew Knepley -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 \ 812*649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 813*649ef022SMatthew 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 \ 814*649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 815*649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 816*649ef022SMatthew Knepley 817*649ef022SMatthew Knepley # TODO Need trig t for convergence in time, also need to refine in space 818*649ef022SMatthew Knepley test: 819*649ef022SMatthew Knepley # Using -dm_refine 5 -convest_num_refine 2 gives L_2 convergence rate: [0.89, 0.011, 1.0] 820*649ef022SMatthew Knepley suffix: 2d_tri_p2_p1_p1_tconv 821*649ef022SMatthew Knepley requires: triangle !single 822*649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type cubic_trig -dm_refine 0 \ 823*649ef022SMatthew Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 824*649ef022SMatthew Knepley -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1 \ 825*649ef022SMatthew Knepley -snes_error_if_not_converged -snes_convergence_test correct_pressure \ 826*649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 827*649ef022SMatthew 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 \ 828*649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 829*649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 830*649ef022SMatthew Knepley 831*649ef022SMatthew Knepley test: 832*649ef022SMatthew Knepley # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.5, 1.9] 833*649ef022SMatthew Knepley suffix: 2d_tri_p2_p1_p1_sconv 834*649ef022SMatthew Knepley requires: triangle !single 835*649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 836*649ef022SMatthew Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 837*649ef022SMatthew Knepley -ts_max_steps 1 -ts_dt 1e-4 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 838*649ef022SMatthew Knepley -snes_error_if_not_converged -snes_convergence_test correct_pressure \ 839*649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_atol 1e-16 -ksp_error_if_not_converged \ 840*649ef022SMatthew 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 \ 841*649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 842*649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 843*649ef022SMatthew Knepley 844*649ef022SMatthew Knepley test: 845*649ef022SMatthew Knepley suffix: 2d_tri_p3_p2_p2 846*649ef022SMatthew Knepley requires: triangle !single 847*649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 848*649ef022SMatthew Knepley -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \ 849*649ef022SMatthew Knepley -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 \ 850*649ef022SMatthew Knepley -snes_convergence_test correct_pressure \ 851*649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 852*649ef022SMatthew 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 \ 853*649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 854*649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 855*649ef022SMatthew Knepley 856*649ef022SMatthew Knepley TEST*/ 857