1649ef022SMatthew Knepley static char help[] = "Time-dependent Low Mach Flow in 2d and 3d channels with finite elements.\n\ 2649ef022SMatthew Knepley We solve the Low Mach flow problem in a rectangular\n\ 3649ef022SMatthew Knepley domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n"; 4649ef022SMatthew Knepley 5649ef022SMatthew Knepley /*F 6649ef022SMatthew Knepley This Low Mach flow is time-dependent isoviscous Navier-Stokes flow. We discretize using the 7649ef022SMatthew Knepley finite element method on an unstructured mesh. The weak form equations are 8649ef022SMatthew Knepley 9649ef022SMatthew Knepley \begin{align*} 10649ef022SMatthew Knepley < q, \nabla\cdot u > = 0 11649ef022SMatthew Knepley <v, du/dt> + <v, u \cdot \nabla u> + < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p > - < v, f > = 0 12649ef022SMatthew Knepley < w, u \cdot \nabla T > + < \nabla w, \alpha \nabla T > - < w, Q > = 0 13649ef022SMatthew Knepley \end{align*} 14649ef022SMatthew Knepley 15649ef022SMatthew Knepley where $\nu$ is the kinematic viscosity and $\alpha$ is thermal diffusivity. 16649ef022SMatthew Knepley 17649ef022SMatthew Knepley For visualization, use 18649ef022SMatthew Knepley 19649ef022SMatthew Knepley -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append 20649ef022SMatthew Knepley F*/ 21649ef022SMatthew Knepley 22649ef022SMatthew Knepley #include <petscdmplex.h> 23649ef022SMatthew Knepley #include <petscsnes.h> 24649ef022SMatthew Knepley #include <petscts.h> 25649ef022SMatthew Knepley #include <petscds.h> 26649ef022SMatthew Knepley #include <petscbag.h> 27649ef022SMatthew Knepley 28606d57d4SMatthew G. Knepley typedef enum {SOL_QUADRATIC, SOL_CUBIC, SOL_CUBIC_TRIG, SOL_TAYLOR_GREEN, NUM_SOL_TYPES} SolType; 29606d57d4SMatthew G. Knepley const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic", "cubic_trig", "taylor_green", "unknown"}; 30649ef022SMatthew Knepley 31649ef022SMatthew Knepley typedef struct { 32649ef022SMatthew Knepley PetscReal nu; /* Kinematic viscosity */ 33649ef022SMatthew Knepley PetscReal alpha; /* Thermal diffusivity */ 34649ef022SMatthew Knepley PetscReal T_in; /* Inlet temperature*/ 35649ef022SMatthew Knepley } Parameter; 36649ef022SMatthew Knepley 37649ef022SMatthew Knepley typedef struct { 38649ef022SMatthew Knepley /* Problem definition */ 39649ef022SMatthew Knepley PetscBag bag; /* Holds problem parameters */ 40649ef022SMatthew Knepley SolType solType; /* MMS solution type */ 41*a712f3bbSMatthew G. Knepley /* Flow diagnostics */ 42*a712f3bbSMatthew G. Knepley DM dmCell; /* A DM with piecewise constant discretization */ 43649ef022SMatthew Knepley } AppCtx; 44649ef022SMatthew Knepley 45649ef022SMatthew Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 46649ef022SMatthew Knepley { 47649ef022SMatthew Knepley PetscInt d; 48649ef022SMatthew Knepley for (d = 0; d < Nc; ++d) u[d] = 0.0; 49649ef022SMatthew Knepley return 0; 50649ef022SMatthew Knepley } 51649ef022SMatthew Knepley 52649ef022SMatthew Knepley static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 53649ef022SMatthew Knepley { 54649ef022SMatthew Knepley PetscInt d; 55649ef022SMatthew Knepley for (d = 0; d < Nc; ++d) u[d] = 1.0; 56649ef022SMatthew Knepley return 0; 57649ef022SMatthew Knepley } 58649ef022SMatthew Knepley 59649ef022SMatthew Knepley /* 60649ef022SMatthew Knepley CASE: quadratic 61649ef022SMatthew Knepley In 2D we use exact solution: 62649ef022SMatthew Knepley 63649ef022SMatthew Knepley u = t + x^2 + y^2 64649ef022SMatthew Knepley v = t + 2x^2 - 2xy 65649ef022SMatthew Knepley p = x + y - 1 66649ef022SMatthew Knepley T = t + x + y 67649ef022SMatthew 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> 68649ef022SMatthew Knepley Q = 1 + 2t + 3x^2 - 2xy + y^2 69649ef022SMatthew Knepley 70649ef022SMatthew Knepley so that 71649ef022SMatthew Knepley 72649ef022SMatthew Knepley \nabla \cdot u = 2x - 2x = 0 73649ef022SMatthew Knepley 74649ef022SMatthew Knepley f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p 75649ef022SMatthew Knepley = <1, 1> + <t + x^2 + y^2, t + 2x^2 - 2xy> . <<2x, 4x - 2y>, <2y, -2x>> - \nu <4, 4> + <1, 1> 76649ef022SMatthew 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> 77649ef022SMatthew 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> 78649ef022SMatthew Knepley 79649ef022SMatthew Knepley Q = dT/dt + u \cdot \nabla T - \alpha \Delta T 80649ef022SMatthew Knepley = 1 + <t + x^2 + y^2, t + 2x^2 - 2xy> . <1, 1> - \alpha 0 81649ef022SMatthew Knepley = 1 + 2t + 3x^2 - 2xy + y^2 82649ef022SMatthew Knepley */ 83649ef022SMatthew Knepley 84649ef022SMatthew Knepley static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 85649ef022SMatthew Knepley { 86649ef022SMatthew Knepley u[0] = time + X[0]*X[0] + X[1]*X[1]; 87649ef022SMatthew Knepley u[1] = time + 2.0*X[0]*X[0] - 2.0*X[0]*X[1]; 88649ef022SMatthew Knepley return 0; 89649ef022SMatthew Knepley } 90649ef022SMatthew Knepley static PetscErrorCode quadratic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 91649ef022SMatthew Knepley { 92649ef022SMatthew Knepley u[0] = 1.0; 93649ef022SMatthew Knepley u[1] = 1.0; 94649ef022SMatthew Knepley return 0; 95649ef022SMatthew Knepley } 96649ef022SMatthew Knepley 97649ef022SMatthew Knepley static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 98649ef022SMatthew Knepley { 99649ef022SMatthew Knepley p[0] = X[0] + X[1] - 1.0; 100649ef022SMatthew Knepley return 0; 101649ef022SMatthew Knepley } 102649ef022SMatthew Knepley 103649ef022SMatthew Knepley static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 104649ef022SMatthew Knepley { 105649ef022SMatthew Knepley T[0] = time + X[0] + X[1]; 106649ef022SMatthew Knepley return 0; 107649ef022SMatthew Knepley } 108649ef022SMatthew Knepley static PetscErrorCode quadratic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 109649ef022SMatthew Knepley { 110649ef022SMatthew Knepley T[0] = 1.0; 111649ef022SMatthew Knepley return 0; 112649ef022SMatthew Knepley } 113649ef022SMatthew Knepley 114649ef022SMatthew Knepley /* f0_v = du/dt - f */ 115649ef022SMatthew Knepley static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 116649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 117649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 118649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 119649ef022SMatthew Knepley { 120649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 121649ef022SMatthew Knepley PetscInt Nc = dim; 122649ef022SMatthew Knepley PetscInt c, d; 123649ef022SMatthew Knepley 124649ef022SMatthew Knepley for (d = 0; d<dim; ++d) f0[d] = u_t[uOff[0]+d]; 125649ef022SMatthew Knepley 126649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 127649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 128649ef022SMatthew Knepley } 129649ef022SMatthew 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); 130649ef022SMatthew 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); 131649ef022SMatthew Knepley } 132649ef022SMatthew Knepley 133649ef022SMatthew Knepley /* f0_w = dT/dt + u.grad(T) - Q */ 134649ef022SMatthew Knepley static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 135649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 136649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 137649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 138649ef022SMatthew Knepley { 139649ef022SMatthew Knepley PetscInt d; 140649ef022SMatthew Knepley f0[0] = 0; 141649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 142649ef022SMatthew Knepley f0[0] += u_t[uOff[2]] - (2*t + 1 + 3*X[0]*X[0] - 2*X[0]*X[1] + X[1]*X[1]); 143649ef022SMatthew Knepley } 144649ef022SMatthew Knepley 145649ef022SMatthew Knepley /* 146649ef022SMatthew Knepley CASE: cubic 147649ef022SMatthew Knepley In 2D we use exact solution: 148649ef022SMatthew Knepley 149649ef022SMatthew Knepley u = t + x^3 + y^3 150649ef022SMatthew Knepley v = t + 2x^3 - 3x^2y 151649ef022SMatthew Knepley p = 3/2 x^2 + 3/2 y^2 - 1 152649ef022SMatthew Knepley T = t + 1/2 x^2 + 1/2 y^2 153649ef022SMatthew Knepley f = < t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x + 1, 154649ef022SMatthew Knepley t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y + 1> 155649ef022SMatthew Knepley Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 + xt + yt - 2\alpha + 1 156649ef022SMatthew Knepley 157649ef022SMatthew Knepley so that 158649ef022SMatthew Knepley 159649ef022SMatthew Knepley \nabla \cdot u = 3x^2 - 3x^2 = 0 160649ef022SMatthew Knepley 161649ef022SMatthew Knepley du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p - f 162649ef022SMatthew 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 163649ef022SMatthew Knepley 164649ef022SMatthew 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 165649ef022SMatthew Knepley */ 166649ef022SMatthew Knepley static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 167649ef022SMatthew Knepley { 168649ef022SMatthew Knepley u[0] = time + X[0]*X[0]*X[0] + X[1]*X[1]*X[1]; 169649ef022SMatthew Knepley u[1] = time + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1]; 170649ef022SMatthew Knepley return 0; 171649ef022SMatthew Knepley } 172649ef022SMatthew Knepley static PetscErrorCode cubic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 173649ef022SMatthew Knepley { 174649ef022SMatthew Knepley u[0] = 1.0; 175649ef022SMatthew Knepley u[1] = 1.0; 176649ef022SMatthew Knepley return 0; 177649ef022SMatthew Knepley } 178649ef022SMatthew Knepley 179649ef022SMatthew Knepley static PetscErrorCode cubic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 180649ef022SMatthew Knepley { 181649ef022SMatthew Knepley p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0; 182649ef022SMatthew Knepley return 0; 183649ef022SMatthew Knepley } 184649ef022SMatthew Knepley 185649ef022SMatthew Knepley static PetscErrorCode cubic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 186649ef022SMatthew Knepley { 187649ef022SMatthew Knepley T[0] = time + X[0]*X[0]/2.0 + X[1]*X[1]/2.0; 188649ef022SMatthew Knepley return 0; 189649ef022SMatthew Knepley } 190649ef022SMatthew Knepley static PetscErrorCode cubic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 191649ef022SMatthew Knepley { 192649ef022SMatthew Knepley T[0] = 1.0; 193649ef022SMatthew Knepley return 0; 194649ef022SMatthew Knepley } 195649ef022SMatthew Knepley 196649ef022SMatthew Knepley 197649ef022SMatthew Knepley static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 198649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 199649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 200649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 201649ef022SMatthew Knepley { 202649ef022SMatthew Knepley PetscInt c, d; 203649ef022SMatthew Knepley PetscInt Nc = dim; 204649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 205649ef022SMatthew Knepley 206649ef022SMatthew Knepley for (d=0; d<dim; ++d) f0[d] = u_t[uOff[0]+d]; 207649ef022SMatthew Knepley 208649ef022SMatthew Knepley for (c=0; c<Nc; ++c) { 209649ef022SMatthew Knepley for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 210649ef022SMatthew Knepley } 211649ef022SMatthew 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); 212649ef022SMatthew 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); 213649ef022SMatthew Knepley } 214649ef022SMatthew Knepley 215649ef022SMatthew Knepley static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 216649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 217649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 218649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 219649ef022SMatthew Knepley { 220649ef022SMatthew Knepley PetscInt d; 221649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 222649ef022SMatthew Knepley 223649ef022SMatthew Knepley for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 224649ef022SMatthew 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); 225649ef022SMatthew Knepley } 226649ef022SMatthew Knepley 227649ef022SMatthew Knepley /* 228649ef022SMatthew Knepley CASE: cubic-trigonometric 229649ef022SMatthew Knepley In 2D we use exact solution: 230649ef022SMatthew Knepley 231649ef022SMatthew Knepley u = beta cos t + x^3 + y^3 232649ef022SMatthew Knepley v = beta sin t + 2x^3 - 3x^2y 233649ef022SMatthew Knepley p = 3/2 x^2 + 3/2 y^2 - 1 234649ef022SMatthew Knepley T = 20 cos t + 1/2 x^2 + 1/2 y^2 235649ef022SMatthew 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, 236649ef022SMatthew Knepley beta cos t (6x^2 - 6xy) - beta sin t (3x^2) + 3x^4y + 6x^2y^3 - 6xy^4 - \nu(12x - 6y) + 3y> 237649ef022SMatthew Knepley Q = beta cos t x + beta sin t (y - 1) + x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2\alpha 238649ef022SMatthew Knepley 239649ef022SMatthew Knepley so that 240649ef022SMatthew Knepley 241649ef022SMatthew Knepley \nabla \cdot u = 3x^2 - 3x^2 = 0 242649ef022SMatthew Knepley 243649ef022SMatthew Knepley f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p 244649ef022SMatthew 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> 245649ef022SMatthew 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> 246649ef022SMatthew Knepley = <cos t (3x^2) + sin t (3y^2 - 1) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu (6x + 6y) + 3x, 247649ef022SMatthew Knepley cos t (6x^2 - 6xy) - sin t (3x^2) + 3x^4y + 6x^2y^3 - 6xy^4 - \nu (12x - 6y) + 3y> 248649ef022SMatthew Knepley 249649ef022SMatthew Knepley Q = dT/dt + u \cdot \nabla T - \alpha \Delta T 250649ef022SMatthew Knepley = -sin t + <cos t + x^3 + y^3, sin t + 2x^3 - 3x^2y> . <x, y> - 2 \alpha 251649ef022SMatthew Knepley = -sin t + cos t (x) + x^4 + xy^3 + sin t (y) + 2x^3y - 3x^2y^2 - 2 \alpha 252649ef022SMatthew Knepley = cos t x + sin t (y - 1) + (x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2 \alpha) 253649ef022SMatthew Knepley */ 254649ef022SMatthew Knepley static PetscErrorCode cubic_trig_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 255649ef022SMatthew Knepley { 256649ef022SMatthew Knepley u[0] = 100.*PetscCosReal(time) + X[0]*X[0]*X[0] + X[1]*X[1]*X[1]; 257649ef022SMatthew Knepley u[1] = 100.*PetscSinReal(time) + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1]; 258649ef022SMatthew Knepley return 0; 259649ef022SMatthew Knepley } 260649ef022SMatthew Knepley static PetscErrorCode cubic_trig_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 261649ef022SMatthew Knepley { 262649ef022SMatthew Knepley u[0] = -100.*PetscSinReal(time); 263649ef022SMatthew Knepley u[1] = 100.*PetscCosReal(time); 264649ef022SMatthew Knepley return 0; 265649ef022SMatthew Knepley } 266649ef022SMatthew Knepley 267649ef022SMatthew Knepley static PetscErrorCode cubic_trig_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 268649ef022SMatthew Knepley { 269649ef022SMatthew Knepley p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0; 270649ef022SMatthew Knepley return 0; 271649ef022SMatthew Knepley } 272649ef022SMatthew Knepley 273649ef022SMatthew Knepley static PetscErrorCode cubic_trig_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 274649ef022SMatthew Knepley { 275649ef022SMatthew Knepley T[0] = 100.*PetscCosReal(time) + X[0]*X[0]/2.0 + X[1]*X[1]/2.0; 276649ef022SMatthew Knepley return 0; 277649ef022SMatthew Knepley } 278649ef022SMatthew Knepley static PetscErrorCode cubic_trig_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 279649ef022SMatthew Knepley { 280649ef022SMatthew Knepley T[0] = -100.*PetscSinReal(time); 281649ef022SMatthew Knepley return 0; 282649ef022SMatthew Knepley } 283649ef022SMatthew Knepley 284649ef022SMatthew Knepley static void f0_cubic_trig_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 285649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 286649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 287649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 288649ef022SMatthew Knepley { 289649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 290649ef022SMatthew Knepley PetscInt Nc = dim; 291649ef022SMatthew Knepley PetscInt c, d; 292649ef022SMatthew Knepley 293649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f0[d] = u_t[uOff[0]+d]; 294649ef022SMatthew Knepley 295649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 296649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 297649ef022SMatthew Knepley } 298649ef022SMatthew 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]; 299649ef022SMatthew 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]; 300649ef022SMatthew Knepley } 301649ef022SMatthew Knepley 302649ef022SMatthew Knepley static void f0_cubic_trig_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 303649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 304649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 305649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 306649ef022SMatthew Knepley { 307649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 308649ef022SMatthew Knepley PetscInt d; 309649ef022SMatthew Knepley 310649ef022SMatthew Knepley for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 311649ef022SMatthew 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); 312649ef022SMatthew Knepley } 313649ef022SMatthew Knepley 314606d57d4SMatthew G. Knepley /* 315606d57d4SMatthew G. Knepley CASE: taylor-green vortex 316606d57d4SMatthew G. Knepley In 2D we use exact solution: 317606d57d4SMatthew G. Knepley 318606d57d4SMatthew G. Knepley u = 1 - cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t) 319606d57d4SMatthew G. Knepley v = 1 + sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t) 320606d57d4SMatthew G. Knepley p = -1/4 [cos(2 \pi(x - t)) + cos(2 \pi(y - t))] exp(-4 \pi^2 \nu t) 321606d57d4SMatthew G. Knepley T = t + x + y 322606d57d4SMatthew G. Knepley f = <\nu \pi^2 exp(-2\nu \pi^2 t) cos(\pi(x-t)) sin(\pi(y-t)), -\nu \pi^2 exp(-2\nu \pi^2 t) sin(\pi(x-t)) cos(\pi(y-t)) > 323606d57d4SMatthew G. Knepley Q = 3 + sin(\pi(x-y)) exp(-2\nu \pi^2 t) 324606d57d4SMatthew G. Knepley 325606d57d4SMatthew G. Knepley so that 326606d57d4SMatthew G. Knepley 327606d57d4SMatthew G. Knepley \nabla \cdot u = \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t) - \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t) = 0 328606d57d4SMatthew G. Knepley 329606d57d4SMatthew G. Knepley f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p 330606d57d4SMatthew G. Knepley = <-\pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t)) - 2\pi cos(\pi(x - t)) sin(\pi(y - t))) exp(-2 \pi^2 \nu t), 331606d57d4SMatthew G. Knepley \pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t)) - 2\pi sin(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t)> 332606d57d4SMatthew G. Knepley + < \pi (1 - cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)) sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t), 333606d57d4SMatthew G. Knepley \pi (1 - cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)) cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)> 334606d57d4SMatthew G. Knepley + <-\pi (1 + sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)) cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t), 335606d57d4SMatthew G. Knepley -\pi (1 + sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)) sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)> 336606d57d4SMatthew G. Knepley + <-2\pi^2 cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t), 337606d57d4SMatthew G. Knepley 2\pi^2 sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)> 338606d57d4SMatthew G. Knepley + < \pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t), 339606d57d4SMatthew G. Knepley \pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)> 340606d57d4SMatthew G. Knepley = <-\pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t)) - 2\pi cos(\pi(x - t)) sin(\pi(y - t))) exp(-2 \pi^2 \nu t), 341606d57d4SMatthew G. Knepley \pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t)) - 2\pi sin(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t)> 342606d57d4SMatthew G. Knepley + < \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t), 343606d57d4SMatthew G. Knepley \pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)> 344606d57d4SMatthew G. Knepley + <-\pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t), 345606d57d4SMatthew G. Knepley -\pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)> 346606d57d4SMatthew G. Knepley + <-\pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t), 347606d57d4SMatthew G. Knepley -\pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)> 348606d57d4SMatthew G. Knepley + <-2\pi^2 cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t), 349606d57d4SMatthew G. Knepley 2\pi^2 sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)> 350606d57d4SMatthew G. Knepley + < \pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t), 351606d57d4SMatthew G. Knepley \pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)> 352606d57d4SMatthew G. Knepley = <-\pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t), 353606d57d4SMatthew G. Knepley \pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t)> 354606d57d4SMatthew G. Knepley + < \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t), 355606d57d4SMatthew G. Knepley \pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)> 356606d57d4SMatthew G. Knepley + <-\pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t), 357606d57d4SMatthew G. Knepley -\pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)> 358606d57d4SMatthew G. Knepley = < \pi cos(\pi(x - t)) cos(\pi(y - t)), 359606d57d4SMatthew G. Knepley \pi sin(\pi(x - t)) sin(\pi(y - t))> 360606d57d4SMatthew G. Knepley + <-\pi cos(\pi(x - t)) cos(\pi(y - t)), 361606d57d4SMatthew G. Knepley -\pi sin(\pi(x - t)) sin(\pi(y - t))> = 0 362606d57d4SMatthew G. Knepley Q = dT/dt + u \cdot \nabla T - \alpha \Delta T 363606d57d4SMatthew G. Knepley = 1 + u \cdot <1, 1> - 0 364606d57d4SMatthew G. Knepley = 1 + u + v 365606d57d4SMatthew G. Knepley */ 366606d57d4SMatthew G. Knepley 367606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 368606d57d4SMatthew G. Knepley { 369606d57d4SMatthew G. Knepley u[0] = 1 - PetscCosReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time); 370606d57d4SMatthew G. Knepley u[1] = 1 + PetscSinReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time); 371606d57d4SMatthew G. Knepley return 0; 372606d57d4SMatthew G. Knepley } 373606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 374606d57d4SMatthew G. Knepley { 375606d57d4SMatthew G. Knepley u[0] = -PETSC_PI*(PetscSinReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time)) 376606d57d4SMatthew G. Knepley - PetscCosReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time)) 377606d57d4SMatthew G. Knepley - 2*PETSC_PI*PetscCosReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time)))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time); 378606d57d4SMatthew G. Knepley u[1] = PETSC_PI*(PetscSinReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time)) 379606d57d4SMatthew G. Knepley - PetscCosReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time)) 380606d57d4SMatthew G. Knepley - 2*PETSC_PI*PetscSinReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time)))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time); 381606d57d4SMatthew G. Knepley return 0; 382606d57d4SMatthew G. Knepley } 383606d57d4SMatthew G. Knepley 384606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 385606d57d4SMatthew G. Knepley { 386606d57d4SMatthew G. Knepley p[0] = -0.25*(PetscCosReal(2*PETSC_PI*(X[0]-time)) + PetscCosReal(2*PETSC_PI*(X[1]-time)))*PetscExpReal(-4*PETSC_PI*PETSC_PI*time); 387606d57d4SMatthew G. Knepley return 0; 388606d57d4SMatthew G. Knepley } 389606d57d4SMatthew G. Knepley 390606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_p_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 391606d57d4SMatthew G. Knepley { 392606d57d4SMatthew G. Knepley p[0] = PETSC_PI*(0.5*(PetscSinReal(2*PETSC_PI*(X[0]-time)) + PetscSinReal(2*PETSC_PI*(X[1]-time))) 393606d57d4SMatthew G. Knepley + PETSC_PI*(PetscCosReal(2*PETSC_PI*(X[0]-time)) + PetscCosReal(2*PETSC_PI*(X[1]-time))))*PetscExpReal(-4*PETSC_PI*PETSC_PI*time); 394606d57d4SMatthew G. Knepley return 0; 395606d57d4SMatthew G. Knepley } 396606d57d4SMatthew G. Knepley 397606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 398606d57d4SMatthew G. Knepley { 399606d57d4SMatthew G. Knepley T[0] = time + X[0] + X[1]; 400606d57d4SMatthew G. Knepley return 0; 401606d57d4SMatthew G. Knepley } 402606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 403606d57d4SMatthew G. Knepley { 404606d57d4SMatthew G. Knepley T[0] = 1.0; 405606d57d4SMatthew G. Knepley return 0; 406606d57d4SMatthew G. Knepley } 407606d57d4SMatthew G. Knepley 408606d57d4SMatthew G. Knepley static void f0_taylor_green_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 409606d57d4SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 410606d57d4SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 411606d57d4SMatthew G. Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 412606d57d4SMatthew G. Knepley { 413606d57d4SMatthew G. Knepley PetscInt Nc = dim; 414606d57d4SMatthew G. Knepley PetscInt c, d; 415606d57d4SMatthew G. Knepley 416606d57d4SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[d] = u_t[uOff[0]+d]; 417606d57d4SMatthew G. Knepley 418606d57d4SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 419606d57d4SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 420606d57d4SMatthew G. Knepley } 421606d57d4SMatthew G. Knepley } 422606d57d4SMatthew G. Knepley 423606d57d4SMatthew G. Knepley static void f0_taylor_green_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 424606d57d4SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 425606d57d4SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 426606d57d4SMatthew G. Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 427606d57d4SMatthew G. Knepley { 428606d57d4SMatthew G. Knepley PetscScalar vel[2]; 429606d57d4SMatthew G. Knepley PetscInt d; 430606d57d4SMatthew G. Knepley 431606d57d4SMatthew G. Knepley taylor_green_u(dim, t, X, Nf, vel, NULL); 432606d57d4SMatthew G. Knepley for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 433606d57d4SMatthew G. Knepley f0[0] += u_t[uOff[2]] - (1.0 + vel[0] + vel[1]); 434606d57d4SMatthew G. Knepley } 435606d57d4SMatthew G. Knepley 436649ef022SMatthew Knepley static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 437649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 438649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 439649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 440649ef022SMatthew Knepley { 441649ef022SMatthew Knepley PetscInt d; 442649ef022SMatthew Knepley for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; 443649ef022SMatthew Knepley } 444649ef022SMatthew Knepley 445649ef022SMatthew Knepley /*f1_v = \nu[grad(u) + grad(u)^T] - pI */ 446649ef022SMatthew Knepley static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 447649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 448649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 449649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 450649ef022SMatthew Knepley { 451649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 452649ef022SMatthew Knepley const PetscInt Nc = dim; 453649ef022SMatthew Knepley PetscInt c, d; 454649ef022SMatthew Knepley 455649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 456649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 457649ef022SMatthew Knepley f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]); 458649ef022SMatthew Knepley } 459649ef022SMatthew Knepley f1[c*dim+c] -= u[uOff[1]]; 460649ef022SMatthew Knepley } 461649ef022SMatthew Knepley } 462649ef022SMatthew Knepley 463649ef022SMatthew Knepley static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 464649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 465649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 466649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 467649ef022SMatthew Knepley { 468649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 469649ef022SMatthew Knepley PetscInt d; 470649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d]; 471649ef022SMatthew Knepley } 472649ef022SMatthew Knepley 473649ef022SMatthew Knepley /*Jacobians*/ 474649ef022SMatthew Knepley static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 475649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 476649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 477649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 478649ef022SMatthew Knepley { 479649ef022SMatthew Knepley PetscInt d; 480649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; 481649ef022SMatthew Knepley } 482649ef022SMatthew Knepley 483649ef022SMatthew Knepley static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 484649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 485649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 486649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 487649ef022SMatthew Knepley { 488649ef022SMatthew Knepley PetscInt c, d; 489649ef022SMatthew Knepley const PetscInt Nc = dim; 490649ef022SMatthew Knepley 491649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g0[d*dim+d] = u_tShift; 492649ef022SMatthew Knepley 493649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 494649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 495649ef022SMatthew Knepley g0[c*Nc+d] += u_x[ c*Nc+d]; 496649ef022SMatthew Knepley } 497649ef022SMatthew Knepley } 498649ef022SMatthew Knepley } 499649ef022SMatthew Knepley 500649ef022SMatthew Knepley static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 501649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 502649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 503649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 504649ef022SMatthew Knepley { 505649ef022SMatthew Knepley PetscInt NcI = dim; 506649ef022SMatthew Knepley PetscInt NcJ = dim; 507649ef022SMatthew Knepley PetscInt c, d, e; 508649ef022SMatthew Knepley 509649ef022SMatthew Knepley for (c = 0; c < NcI; ++c) { 510649ef022SMatthew Knepley for (d = 0; d < NcJ; ++d) { 511649ef022SMatthew Knepley for (e = 0; e < dim; ++e) { 512649ef022SMatthew Knepley if (c == d) { 513649ef022SMatthew Knepley g1[(c*NcJ+d)*dim+e] += u[e]; 514649ef022SMatthew Knepley } 515649ef022SMatthew Knepley } 516649ef022SMatthew Knepley } 517649ef022SMatthew Knepley } 518649ef022SMatthew Knepley } 519649ef022SMatthew Knepley 520649ef022SMatthew Knepley 521649ef022SMatthew Knepley static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 522649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 523649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 524649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 525649ef022SMatthew Knepley { 526649ef022SMatthew Knepley PetscInt d; 527649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; 528649ef022SMatthew Knepley } 529649ef022SMatthew Knepley 530649ef022SMatthew Knepley static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 531649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 532649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 533649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 534649ef022SMatthew Knepley { 535649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 536649ef022SMatthew Knepley const PetscInt Nc = dim; 537649ef022SMatthew Knepley PetscInt c, d; 538649ef022SMatthew Knepley 539649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 540649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 541606d57d4SMatthew G. Knepley g3[((c*Nc+c)*dim+d)*dim+d] += nu; 542606d57d4SMatthew G. Knepley g3[((c*Nc+d)*dim+d)*dim+c] += nu; 543649ef022SMatthew Knepley } 544649ef022SMatthew Knepley } 545649ef022SMatthew Knepley } 546649ef022SMatthew Knepley 547649ef022SMatthew Knepley static void g0_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 548649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 549649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 550649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 551649ef022SMatthew Knepley { 552*a712f3bbSMatthew G. Knepley g0[0] = u_tShift; 553649ef022SMatthew Knepley } 554649ef022SMatthew Knepley 555649ef022SMatthew Knepley static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 556649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 557649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 558649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 559649ef022SMatthew Knepley { 560649ef022SMatthew Knepley PetscInt d; 561649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d]; 562649ef022SMatthew Knepley } 563649ef022SMatthew Knepley 564649ef022SMatthew Knepley static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 565649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 566649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 567649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 568649ef022SMatthew Knepley { 569649ef022SMatthew Knepley PetscInt d; 570649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d]; 571649ef022SMatthew Knepley } 572649ef022SMatthew Knepley 573649ef022SMatthew Knepley static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 574649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 575649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 576649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 577649ef022SMatthew Knepley { 578649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 579649ef022SMatthew Knepley PetscInt d; 580649ef022SMatthew Knepley 581649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha; 582649ef022SMatthew Knepley } 583649ef022SMatthew Knepley 584649ef022SMatthew Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 585649ef022SMatthew Knepley { 586649ef022SMatthew Knepley PetscInt sol; 587649ef022SMatthew Knepley PetscErrorCode ierr; 588649ef022SMatthew Knepley 589649ef022SMatthew Knepley 590649ef022SMatthew Knepley PetscFunctionBeginUser; 591649ef022SMatthew Knepley options->solType = SOL_QUADRATIC; 592649ef022SMatthew Knepley 593649ef022SMatthew Knepley ierr = PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX");CHKERRQ(ierr); 594649ef022SMatthew Knepley sol = options->solType; 595649ef022SMatthew Knepley ierr = PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr); 596649ef022SMatthew Knepley options->solType = (SolType) sol; 597649ef022SMatthew Knepley ierr = PetscOptionsEnd(); 598649ef022SMatthew Knepley PetscFunctionReturn(0); 599649ef022SMatthew Knepley } 600649ef022SMatthew Knepley 601649ef022SMatthew Knepley static PetscErrorCode SetupParameters(AppCtx *user) 602649ef022SMatthew Knepley { 603649ef022SMatthew Knepley PetscBag bag; 604649ef022SMatthew Knepley Parameter *p; 605649ef022SMatthew Knepley PetscErrorCode ierr; 606649ef022SMatthew Knepley 607649ef022SMatthew Knepley PetscFunctionBeginUser; 608649ef022SMatthew Knepley /* setup PETSc parameter bag */ 609649ef022SMatthew Knepley ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr); 610649ef022SMatthew Knepley ierr = PetscBagSetName(user->bag, "par", "Low Mach flow parameters");CHKERRQ(ierr); 611649ef022SMatthew Knepley bag = user->bag; 612649ef022SMatthew Knepley ierr = PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity");CHKERRQ(ierr); 613649ef022SMatthew Knepley ierr = PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity");CHKERRQ(ierr); 614649ef022SMatthew Knepley ierr = PetscBagRegisterReal(bag, &p->T_in, 1.0, "T_in", "Inlet temperature");CHKERRQ(ierr); 615649ef022SMatthew Knepley PetscFunctionReturn(0); 616649ef022SMatthew Knepley } 617649ef022SMatthew Knepley 618649ef022SMatthew Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 619649ef022SMatthew Knepley { 620649ef022SMatthew Knepley PetscErrorCode ierr; 621649ef022SMatthew Knepley 622649ef022SMatthew Knepley PetscFunctionBeginUser; 623649ef022SMatthew Knepley ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 624649ef022SMatthew Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 625649ef022SMatthew Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 626649ef022SMatthew Knepley PetscFunctionReturn(0); 627649ef022SMatthew Knepley } 628649ef022SMatthew Knepley 629649ef022SMatthew Knepley static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 630649ef022SMatthew Knepley { 631649ef022SMatthew Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 632649ef022SMatthew Knepley PetscErrorCode (*exactFuncs_t[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 633649ef022SMatthew Knepley PetscDS prob; 634649ef022SMatthew Knepley Parameter *ctx; 635649ef022SMatthew Knepley PetscInt id; 636649ef022SMatthew Knepley PetscErrorCode ierr; 637649ef022SMatthew Knepley 638649ef022SMatthew Knepley PetscFunctionBeginUser; 639649ef022SMatthew Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 640649ef022SMatthew Knepley switch(user->solType){ 641649ef022SMatthew Knepley case SOL_QUADRATIC: 642649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v);CHKERRQ(ierr); 643649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w);CHKERRQ(ierr); 644649ef022SMatthew Knepley 645649ef022SMatthew Knepley exactFuncs[0] = quadratic_u; 646649ef022SMatthew Knepley exactFuncs[1] = quadratic_p; 647649ef022SMatthew Knepley exactFuncs[2] = quadratic_T; 648649ef022SMatthew Knepley exactFuncs_t[0] = quadratic_u_t; 649649ef022SMatthew Knepley exactFuncs_t[1] = NULL; 650649ef022SMatthew Knepley exactFuncs_t[2] = quadratic_T_t; 651649ef022SMatthew Knepley break; 652649ef022SMatthew Knepley case SOL_CUBIC: 653649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v);CHKERRQ(ierr); 654649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w);CHKERRQ(ierr); 655649ef022SMatthew Knepley 656649ef022SMatthew Knepley exactFuncs[0] = cubic_u; 657649ef022SMatthew Knepley exactFuncs[1] = cubic_p; 658649ef022SMatthew Knepley exactFuncs[2] = cubic_T; 659649ef022SMatthew Knepley exactFuncs_t[0] = cubic_u_t; 660649ef022SMatthew Knepley exactFuncs_t[1] = NULL; 661649ef022SMatthew Knepley exactFuncs_t[2] = cubic_T_t; 662649ef022SMatthew Knepley break; 663649ef022SMatthew Knepley case SOL_CUBIC_TRIG: 664649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 0, f0_cubic_trig_v, f1_v);CHKERRQ(ierr); 665649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 2, f0_cubic_trig_w, f1_w);CHKERRQ(ierr); 666649ef022SMatthew Knepley 667649ef022SMatthew Knepley exactFuncs[0] = cubic_trig_u; 668649ef022SMatthew Knepley exactFuncs[1] = cubic_trig_p; 669649ef022SMatthew Knepley exactFuncs[2] = cubic_trig_T; 670649ef022SMatthew Knepley exactFuncs_t[0] = cubic_trig_u_t; 671649ef022SMatthew Knepley exactFuncs_t[1] = NULL; 672649ef022SMatthew Knepley exactFuncs_t[2] = cubic_trig_T_t; 673649ef022SMatthew Knepley break; 674606d57d4SMatthew G. Knepley case SOL_TAYLOR_GREEN: 675606d57d4SMatthew G. Knepley ierr = PetscDSSetResidual(prob, 0, f0_taylor_green_v, f1_v);CHKERRQ(ierr); 676606d57d4SMatthew G. Knepley ierr = PetscDSSetResidual(prob, 2, f0_taylor_green_w, f1_w);CHKERRQ(ierr); 677606d57d4SMatthew G. Knepley 678606d57d4SMatthew G. Knepley exactFuncs[0] = taylor_green_u; 679606d57d4SMatthew G. Knepley exactFuncs[1] = taylor_green_p; 680606d57d4SMatthew G. Knepley exactFuncs[2] = taylor_green_T; 681606d57d4SMatthew G. Knepley exactFuncs_t[0] = taylor_green_u_t; 682606d57d4SMatthew G. Knepley exactFuncs_t[1] = taylor_green_p_t; 683606d57d4SMatthew G. Knepley exactFuncs_t[2] = taylor_green_T_t; 684606d57d4SMatthew G. Knepley break; 685649ef022SMatthew Knepley default: SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 686649ef022SMatthew Knepley } 687649ef022SMatthew Knepley 688649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 1, f0_q, NULL);CHKERRQ(ierr); 689649ef022SMatthew Knepley 690649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu);CHKERRQ(ierr); 691649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL);CHKERRQ(ierr); 692649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL);CHKERRQ(ierr); 693649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL);CHKERRQ(ierr); 694649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 2, 2, g0_wT, g1_wT, NULL, g3_wT);CHKERRQ(ierr); 695649ef022SMatthew Knepley /* Setup constants */ 696649ef022SMatthew Knepley { 697649ef022SMatthew Knepley Parameter *param; 698649ef022SMatthew Knepley PetscScalar constants[3]; 699649ef022SMatthew Knepley 700649ef022SMatthew Knepley ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 701649ef022SMatthew Knepley 702649ef022SMatthew Knepley constants[0] = param->nu; 703649ef022SMatthew Knepley constants[1] = param->alpha; 704649ef022SMatthew Knepley constants[2] = param->T_in; 705649ef022SMatthew Knepley ierr = PetscDSSetConstants(prob, 3, constants);CHKERRQ(ierr); 706649ef022SMatthew Knepley } 707649ef022SMatthew Knepley /* Setup Boundary Conditions */ 708649ef022SMatthew Knepley ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr); 709649ef022SMatthew Knepley id = 3; 710649ef022SMatthew 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); 711649ef022SMatthew Knepley id = 1; 712649ef022SMatthew 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); 713649ef022SMatthew Knepley id = 2; 714649ef022SMatthew 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); 715649ef022SMatthew Knepley id = 4; 716649ef022SMatthew 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); 717649ef022SMatthew Knepley id = 3; 718649ef022SMatthew 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); 719649ef022SMatthew Knepley id = 1; 720649ef022SMatthew 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); 721649ef022SMatthew Knepley id = 2; 722649ef022SMatthew 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); 723649ef022SMatthew Knepley id = 4; 724649ef022SMatthew 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); 725649ef022SMatthew Knepley 726649ef022SMatthew Knepley /*setup exact solution.*/ 727649ef022SMatthew Knepley ierr = PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx);CHKERRQ(ierr); 728649ef022SMatthew Knepley ierr = PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx);CHKERRQ(ierr); 729649ef022SMatthew Knepley ierr = PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx);CHKERRQ(ierr); 730649ef022SMatthew Knepley ierr = PetscDSSetExactSolutionTimeDerivative(prob, 0, exactFuncs_t[0], ctx);CHKERRQ(ierr); 731649ef022SMatthew Knepley ierr = PetscDSSetExactSolutionTimeDerivative(prob, 1, exactFuncs_t[1], ctx);CHKERRQ(ierr); 732649ef022SMatthew Knepley ierr = PetscDSSetExactSolutionTimeDerivative(prob, 2, exactFuncs_t[2], ctx);CHKERRQ(ierr); 733649ef022SMatthew Knepley PetscFunctionReturn(0); 734649ef022SMatthew Knepley } 735649ef022SMatthew Knepley 736649ef022SMatthew Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 737649ef022SMatthew Knepley { 738649ef022SMatthew Knepley DM cdm = dm; 739*a712f3bbSMatthew G. Knepley PetscFE fe[3], fediv; 740649ef022SMatthew Knepley Parameter *param; 741649ef022SMatthew Knepley DMPolytopeType ct; 742649ef022SMatthew Knepley PetscInt dim, cStart; 743649ef022SMatthew Knepley PetscBool simplex; 744649ef022SMatthew Knepley PetscErrorCode ierr; 745649ef022SMatthew Knepley 746649ef022SMatthew Knepley PetscFunctionBeginUser; 747649ef022SMatthew Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 748649ef022SMatthew Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 749649ef022SMatthew Knepley ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 750649ef022SMatthew Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 751649ef022SMatthew Knepley /* Create finite element */ 752*a712f3bbSMatthew G. Knepley ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 753649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr); 754649ef022SMatthew Knepley 755*a712f3bbSMatthew G. Knepley ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 756649ef022SMatthew Knepley ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 757649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr); 758649ef022SMatthew Knepley 759*a712f3bbSMatthew G. Knepley ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]);CHKERRQ(ierr); 760649ef022SMatthew Knepley ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr); 761649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) fe[2], "temperature");CHKERRQ(ierr); 762649ef022SMatthew Knepley 763*a712f3bbSMatthew G. Knepley ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "div_", PETSC_DEFAULT, &fediv);CHKERRQ(ierr); 764*a712f3bbSMatthew G. Knepley ierr = PetscFECopyQuadrature(fe[0], fediv);CHKERRQ(ierr); 765*a712f3bbSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) fediv, "divergence");CHKERRQ(ierr); 766*a712f3bbSMatthew G. Knepley 767649ef022SMatthew Knepley /* Set discretization and boundary conditions for each mesh */ 768649ef022SMatthew Knepley ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 769649ef022SMatthew Knepley ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr); 770649ef022SMatthew Knepley ierr = DMSetField(dm, 2, NULL, (PetscObject) fe[2]);CHKERRQ(ierr); 771649ef022SMatthew Knepley ierr = DMCreateDS(dm);CHKERRQ(ierr); 772649ef022SMatthew Knepley ierr = SetupProblem(dm, user);CHKERRQ(ierr); 773649ef022SMatthew Knepley ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 774649ef022SMatthew Knepley while (cdm) { 775649ef022SMatthew Knepley ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 776649ef022SMatthew Knepley ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 777649ef022SMatthew Knepley } 778649ef022SMatthew Knepley ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); 779649ef022SMatthew Knepley ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); 780649ef022SMatthew Knepley ierr = PetscFEDestroy(&fe[2]);CHKERRQ(ierr); 781649ef022SMatthew Knepley 782*a712f3bbSMatthew G. Knepley ierr = DMClone(dm, &user->dmCell);CHKERRQ(ierr); 783*a712f3bbSMatthew G. Knepley ierr = DMSetField(user->dmCell, 0, NULL, (PetscObject) fediv);CHKERRQ(ierr); 784*a712f3bbSMatthew G. Knepley ierr = DMCreateDS(user->dmCell);CHKERRQ(ierr); 785*a712f3bbSMatthew G. Knepley ierr = PetscFEDestroy(&fediv);CHKERRQ(ierr); 786*a712f3bbSMatthew G. Knepley 787649ef022SMatthew Knepley { 788649ef022SMatthew Knepley PetscObject pressure; 789649ef022SMatthew Knepley MatNullSpace nullspacePres; 790649ef022SMatthew Knepley 791649ef022SMatthew Knepley ierr = DMGetField(dm, 1, NULL, &pressure);CHKERRQ(ierr); 792649ef022SMatthew Knepley ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres);CHKERRQ(ierr); 793649ef022SMatthew Knepley ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres);CHKERRQ(ierr); 794649ef022SMatthew Knepley ierr = MatNullSpaceDestroy(&nullspacePres);CHKERRQ(ierr); 795649ef022SMatthew Knepley } 796649ef022SMatthew Knepley 797649ef022SMatthew Knepley PetscFunctionReturn(0); 798649ef022SMatthew Knepley } 799649ef022SMatthew Knepley 800649ef022SMatthew Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) 801649ef022SMatthew Knepley { 802649ef022SMatthew Knepley Vec vec; 803649ef022SMatthew Knepley PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 804649ef022SMatthew Knepley PetscErrorCode ierr; 805649ef022SMatthew Knepley 806649ef022SMatthew Knepley PetscFunctionBeginUser; 807649ef022SMatthew Knepley if (ofield != 1) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield); 808649ef022SMatthew Knepley funcs[nfield] = constant; 809649ef022SMatthew Knepley ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr); 810649ef022SMatthew Knepley ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr); 811649ef022SMatthew Knepley ierr = VecNormalize(vec, NULL);CHKERRQ(ierr); 812649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr); 813649ef022SMatthew Knepley ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr); 814649ef022SMatthew Knepley ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace);CHKERRQ(ierr); 815649ef022SMatthew Knepley ierr = VecDestroy(&vec);CHKERRQ(ierr); 816649ef022SMatthew Knepley PetscFunctionReturn(0); 817649ef022SMatthew Knepley } 818649ef022SMatthew Knepley 819649ef022SMatthew Knepley static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u) 820649ef022SMatthew Knepley { 821649ef022SMatthew Knepley DM dm; 822649ef022SMatthew Knepley MatNullSpace nullsp; 823649ef022SMatthew Knepley PetscErrorCode ierr; 824649ef022SMatthew Knepley 825649ef022SMatthew Knepley PetscFunctionBegin; 826649ef022SMatthew Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 827649ef022SMatthew Knepley ierr = CreatePressureNullSpace(dm, 1, 1, &nullsp);CHKERRQ(ierr); 828649ef022SMatthew Knepley ierr = MatNullSpaceRemove(nullsp, u);CHKERRQ(ierr); 829649ef022SMatthew Knepley ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 830649ef022SMatthew Knepley PetscFunctionReturn(0); 831649ef022SMatthew Knepley } 832649ef022SMatthew Knepley 833649ef022SMatthew Knepley /* Make the discrete pressure discretely divergence free */ 834649ef022SMatthew Knepley static PetscErrorCode RemoveDiscretePressureNullspace(TS ts) 835649ef022SMatthew Knepley { 836649ef022SMatthew Knepley Vec u; 837649ef022SMatthew Knepley PetscErrorCode ierr; 838649ef022SMatthew Knepley 839649ef022SMatthew Knepley PetscFunctionBegin; 840649ef022SMatthew Knepley ierr = TSGetSolution(ts, &u);CHKERRQ(ierr); 841649ef022SMatthew Knepley ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr); 842649ef022SMatthew Knepley PetscFunctionReturn(0); 843649ef022SMatthew Knepley } 844649ef022SMatthew Knepley 845649ef022SMatthew Knepley static PetscErrorCode SetInitialConditions(TS ts, Vec u) 846649ef022SMatthew Knepley { 847649ef022SMatthew Knepley DM dm; 848649ef022SMatthew Knepley PetscReal t; 849649ef022SMatthew Knepley PetscErrorCode ierr; 850649ef022SMatthew Knepley 851649ef022SMatthew Knepley PetscFunctionBegin; 852649ef022SMatthew Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 853649ef022SMatthew Knepley ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 854649ef022SMatthew Knepley ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr); 855649ef022SMatthew Knepley ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr); 856649ef022SMatthew Knepley PetscFunctionReturn(0); 857649ef022SMatthew Knepley } 858649ef022SMatthew Knepley 859*a712f3bbSMatthew G. Knepley static void divergence(PetscInt dim, PetscInt Nf, PetscInt NfAux, 860*a712f3bbSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 861*a712f3bbSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 862*a712f3bbSMatthew G. Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar divu[]) 863*a712f3bbSMatthew G. Knepley { 864*a712f3bbSMatthew G. Knepley PetscInt d; 865*a712f3bbSMatthew G. Knepley 866*a712f3bbSMatthew G. Knepley divu[0] = 0.; 867*a712f3bbSMatthew G. Knepley for (d = 0; d < dim; ++d) divu[0] += u_x[d*dim+d]; 868*a712f3bbSMatthew G. Knepley } 869*a712f3bbSMatthew G. Knepley 870649ef022SMatthew Knepley static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) 871649ef022SMatthew Knepley { 872649ef022SMatthew Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 873649ef022SMatthew Knepley void *ctxs[3]; 874*a712f3bbSMatthew G. Knepley PetscPointFunc diagnostics[1] = {divergence}; 875*a712f3bbSMatthew G. Knepley DM dm, dmCell = ((AppCtx *) ctx)->dmCell; 876649ef022SMatthew Knepley PetscDS ds; 877*a712f3bbSMatthew G. Knepley Vec v, divu; 878*a712f3bbSMatthew G. Knepley PetscReal ferrors[3], massFlux; 879649ef022SMatthew Knepley PetscInt f; 880649ef022SMatthew Knepley PetscErrorCode ierr; 881649ef022SMatthew Knepley 882649ef022SMatthew Knepley PetscFunctionBeginUser; 883649ef022SMatthew Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 884649ef022SMatthew Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 885649ef022SMatthew Knepley 886649ef022SMatthew Knepley for (f = 0; f < 3; ++f) {ierr = PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]);CHKERRQ(ierr);} 887649ef022SMatthew Knepley ierr = DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors);CHKERRQ(ierr); 888*a712f3bbSMatthew G. Knepley ierr = DMGetGlobalVector(dmCell, &divu);CHKERRQ(ierr); 889*a712f3bbSMatthew G. Knepley ierr = DMProjectField(dmCell, crtime, u, diagnostics, INSERT_VALUES, divu);CHKERRQ(ierr); 890*a712f3bbSMatthew G. Knepley ierr = VecViewFromOptions(divu, NULL, "-divu_vec_view");CHKERRQ(ierr); 891*a712f3bbSMatthew G. Knepley ierr = VecNorm(divu, NORM_2, &massFlux);CHKERRQ(ierr); 892*a712f3bbSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [%2.3g, %2.3g, %2.3g] ||div u||: %2.3g\n", (int) step, (double) crtime, (double) ferrors[0], (double) ferrors[1], (double) ferrors[2], (double) massFlux);CHKERRQ(ierr); 893649ef022SMatthew Knepley 894649ef022SMatthew Knepley ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 895649ef022SMatthew Knepley 896649ef022SMatthew Knepley ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr); 897*a712f3bbSMatthew G. Knepley ierr = DMProjectFunction(dm, crtime, exactFuncs, ctxs, INSERT_ALL_VALUES, v);CHKERRQ(ierr); 898649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) v, "Exact Solution");CHKERRQ(ierr); 899649ef022SMatthew Knepley ierr = VecViewFromOptions(v, NULL, "-exact_vec_view");CHKERRQ(ierr); 900649ef022SMatthew Knepley ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr); 901649ef022SMatthew Knepley 902*a712f3bbSMatthew G. Knepley ierr = VecViewFromOptions(divu, NULL, "-div_vec_view");CHKERRQ(ierr); 903*a712f3bbSMatthew G. Knepley ierr = DMRestoreGlobalVector(dmCell, &divu);CHKERRQ(ierr); 904*a712f3bbSMatthew G. Knepley 905649ef022SMatthew Knepley PetscFunctionReturn(0); 906649ef022SMatthew Knepley } 907649ef022SMatthew Knepley 908649ef022SMatthew Knepley int main(int argc, char **argv) 909649ef022SMatthew Knepley { 910649ef022SMatthew Knepley DM dm; /* problem definition */ 911649ef022SMatthew Knepley TS ts; /* timestepper */ 912649ef022SMatthew Knepley Vec u; /* solution */ 913649ef022SMatthew Knepley AppCtx user; /* user-defined work context */ 914649ef022SMatthew Knepley PetscReal t; 915649ef022SMatthew Knepley PetscErrorCode ierr; 916649ef022SMatthew Knepley 917649ef022SMatthew Knepley ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 918649ef022SMatthew Knepley ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 919649ef022SMatthew Knepley ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr); 920649ef022SMatthew Knepley ierr = SetupParameters(&user);CHKERRQ(ierr); 921649ef022SMatthew Knepley ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 922649ef022SMatthew Knepley ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 923649ef022SMatthew Knepley ierr = TSSetDM(ts, dm);CHKERRQ(ierr); 924649ef022SMatthew Knepley ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr); 925649ef022SMatthew Knepley /* Setup problem */ 926649ef022SMatthew Knepley ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr); 927649ef022SMatthew Knepley ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr); 928649ef022SMatthew Knepley 929649ef022SMatthew Knepley ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 930*a712f3bbSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr); 931649ef022SMatthew Knepley ierr = DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace);CHKERRQ(ierr); 932649ef022SMatthew Knepley 933649ef022SMatthew Knepley ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr); 934649ef022SMatthew Knepley ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr); 935649ef022SMatthew Knepley ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr); 936649ef022SMatthew Knepley ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 937649ef022SMatthew Knepley ierr = TSSetPreStep(ts, RemoveDiscretePressureNullspace);CHKERRQ(ierr); 938649ef022SMatthew Knepley ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 939649ef022SMatthew Knepley 940649ef022SMatthew Knepley ierr = TSSetComputeInitialCondition(ts, SetInitialConditions);CHKERRQ(ierr); /* Must come after SetFromOptions() */ 941649ef022SMatthew Knepley ierr = SetInitialConditions(ts, u);CHKERRQ(ierr); 942649ef022SMatthew Knepley ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 943649ef022SMatthew Knepley ierr = DMSetOutputSequenceNumber(dm, 0, t);CHKERRQ(ierr); 944649ef022SMatthew Knepley ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 945649ef022SMatthew Knepley ierr = TSMonitorSet(ts, MonitorError, &user, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 946649ef022SMatthew Knepley 947649ef022SMatthew Knepley ierr = TSSolve(ts, u);CHKERRQ(ierr); 948649ef022SMatthew Knepley ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 949649ef022SMatthew Knepley 950649ef022SMatthew Knepley ierr = VecDestroy(&u);CHKERRQ(ierr); 951*a712f3bbSMatthew G. Knepley ierr = DMDestroy(&user.dmCell);CHKERRQ(ierr); 952649ef022SMatthew Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 953649ef022SMatthew Knepley ierr = TSDestroy(&ts);CHKERRQ(ierr); 954649ef022SMatthew Knepley ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr); 955649ef022SMatthew Knepley ierr = PetscFinalize(); 956649ef022SMatthew Knepley return ierr; 957649ef022SMatthew Knepley } 958649ef022SMatthew Knepley 959649ef022SMatthew Knepley /*TEST 960649ef022SMatthew Knepley 961649ef022SMatthew Knepley test: 962649ef022SMatthew Knepley suffix: 2d_tri_p2_p1_p1 963649ef022SMatthew Knepley requires: triangle !single 964649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \ 965649ef022SMatthew Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 966*a712f3bbSMatthew G. Knepley -div_petscdualspace_lagrange_use_moments -div_petscdualspace_lagrange_moment_order 3 \ 967649ef022SMatthew Knepley -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 \ 968649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 969649ef022SMatthew 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 \ 970649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 971649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 972649ef022SMatthew Knepley 973649ef022SMatthew Knepley # TODO Need trig t for convergence in time, also need to refine in space 974649ef022SMatthew Knepley test: 975649ef022SMatthew Knepley # Using -dm_refine 5 -convest_num_refine 2 gives L_2 convergence rate: [0.89, 0.011, 1.0] 976649ef022SMatthew Knepley suffix: 2d_tri_p2_p1_p1_tconv 977649ef022SMatthew Knepley requires: triangle !single 978649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type cubic_trig -dm_refine 0 \ 979649ef022SMatthew Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 980*a712f3bbSMatthew G. Knepley -div_petscdualspace_lagrange_use_moments -div_petscdualspace_lagrange_moment_order 3 \ 981649ef022SMatthew Knepley -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1 \ 982649ef022SMatthew Knepley -snes_error_if_not_converged -snes_convergence_test correct_pressure \ 983649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 984649ef022SMatthew 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 \ 985649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 986649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 987649ef022SMatthew Knepley 988649ef022SMatthew Knepley test: 989649ef022SMatthew Knepley # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.5, 1.9] 990649ef022SMatthew Knepley suffix: 2d_tri_p2_p1_p1_sconv 991649ef022SMatthew Knepley requires: triangle !single 992649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 993649ef022SMatthew Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 994*a712f3bbSMatthew G. Knepley -div_petscdualspace_lagrange_use_moments -div_petscdualspace_lagrange_moment_order 3 \ 995649ef022SMatthew Knepley -ts_max_steps 1 -ts_dt 1e-4 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 996649ef022SMatthew Knepley -snes_error_if_not_converged -snes_convergence_test correct_pressure \ 997649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_atol 1e-16 -ksp_error_if_not_converged \ 998649ef022SMatthew 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 \ 999649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 1000649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 1001649ef022SMatthew Knepley 1002649ef022SMatthew Knepley test: 1003649ef022SMatthew Knepley suffix: 2d_tri_p3_p2_p2 1004649ef022SMatthew Knepley requires: triangle !single 1005649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 1006649ef022SMatthew Knepley -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \ 1007*a712f3bbSMatthew G. Knepley -div_petscdualspace_lagrange_use_moments -div_petscdualspace_lagrange_moment_order 3 \ 1008649ef022SMatthew Knepley -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 \ 1009649ef022SMatthew Knepley -snes_convergence_test correct_pressure \ 1010649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1011649ef022SMatthew 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 \ 1012649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 1013649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 1014649ef022SMatthew Knepley 1015606d57d4SMatthew G. Knepley test: 1016606d57d4SMatthew G. Knepley # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.1, 3.1] 1017606d57d4SMatthew G. Knepley suffix: 2d_tri_p2_p1_p1_tg_sconv 1018606d57d4SMatthew G. Knepley requires: triangle !single 1019606d57d4SMatthew G. Knepley args: -dm_plex_separate_marker -sol_type taylor_green -dm_refine 0 \ 1020606d57d4SMatthew G. Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1021*a712f3bbSMatthew G. Knepley -div_petscdualspace_lagrange_use_moments -div_petscdualspace_lagrange_moment_order 3 \ 1022606d57d4SMatthew G. Knepley -ts_max_steps 1 -ts_dt 1e-8 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 1023606d57d4SMatthew G. Knepley -snes_error_if_not_converged -snes_convergence_test correct_pressure \ 1024606d57d4SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_atol 1e-16 -ksp_error_if_not_converged \ 1025606d57d4SMatthew 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 \ 1026606d57d4SMatthew G. Knepley -fieldsplit_0_pc_type lu \ 1027606d57d4SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 1028606d57d4SMatthew G. Knepley 1029606d57d4SMatthew G. Knepley test: 1030606d57d4SMatthew G. Knepley # Using -dm_refine 3 -convest_num_refine 2 gives L_2 convergence rate: [1.2, 1.5, 1.2] 1031606d57d4SMatthew G. Knepley suffix: 2d_tri_p2_p1_p1_tg_tconv 1032606d57d4SMatthew G. Knepley requires: triangle !single 1033606d57d4SMatthew G. Knepley args: -dm_plex_separate_marker -sol_type taylor_green -dm_refine 0 \ 1034606d57d4SMatthew G. Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1035*a712f3bbSMatthew G. Knepley -div_petscdualspace_lagrange_use_moments -div_petscdualspace_lagrange_moment_order 3 \ 1036606d57d4SMatthew G. Knepley -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1 \ 1037606d57d4SMatthew G. Knepley -snes_error_if_not_converged -snes_convergence_test correct_pressure \ 1038606d57d4SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_atol 1e-16 -ksp_error_if_not_converged \ 1039606d57d4SMatthew 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 \ 1040606d57d4SMatthew G. Knepley -fieldsplit_0_pc_type lu \ 1041606d57d4SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 1042606d57d4SMatthew G. Knepley 1043649ef022SMatthew Knepley TEST*/ 1044