1*649ef022SMatthew Knepley static char help[] = "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 a steady-state 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, 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 <petscds.h> 25*649ef022SMatthew Knepley #include <petscbag.h> 26*649ef022SMatthew Knepley 27*649ef022SMatthew Knepley typedef enum {SOL_QUADRATIC, SOL_CUBIC, NUM_SOL_TYPES} SolType; 28*649ef022SMatthew Knepley const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic", "unknown"}; 29*649ef022SMatthew Knepley 30*649ef022SMatthew Knepley typedef struct { 31*649ef022SMatthew Knepley PetscReal nu; /* Kinematic viscosity */ 32*649ef022SMatthew Knepley PetscReal theta; /* Angle of pipe wall to x-axis */ 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 39*649ef022SMatthew Knepley PetscBool showError; /*showSolution */ 40*649ef022SMatthew Knepley /* Domain and mesh definition */ 41*649ef022SMatthew Knepley PetscInt dim; /* The topological mesh dimension */ 42*649ef022SMatthew Knepley PetscBool simplex; /* Use simplices or tensor product cells */ 43*649ef022SMatthew Knepley PetscInt cells[3]; /* The initial domain division */ 44*649ef022SMatthew Knepley /* Problem definition */ 45*649ef022SMatthew Knepley PetscBag bag; /* Holds problem parameters */ 46*649ef022SMatthew Knepley SolType solType; 47*649ef022SMatthew Knepley } AppCtx; 48*649ef022SMatthew Knepley 49*649ef022SMatthew Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 50*649ef022SMatthew Knepley { 51*649ef022SMatthew Knepley PetscInt d; 52*649ef022SMatthew Knepley for (d = 0; d < Nc; ++d) u[d] = 0.0; 53*649ef022SMatthew Knepley return 0; 54*649ef022SMatthew Knepley } 55*649ef022SMatthew Knepley 56*649ef022SMatthew Knepley static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 57*649ef022SMatthew Knepley { 58*649ef022SMatthew Knepley PetscInt d; 59*649ef022SMatthew Knepley for (d = 0; d < Nc; ++d) u[d] = 1.0; 60*649ef022SMatthew Knepley return 0; 61*649ef022SMatthew Knepley } 62*649ef022SMatthew Knepley 63*649ef022SMatthew Knepley /* 64*649ef022SMatthew Knepley CASE: quadratic 65*649ef022SMatthew Knepley In 2D we use exact solution: 66*649ef022SMatthew Knepley 67*649ef022SMatthew Knepley u = x^2 + y^2 68*649ef022SMatthew Knepley v = 2x^2 - 2xy 69*649ef022SMatthew Knepley p = x + y - 1 70*649ef022SMatthew Knepley T = x + y 71*649ef022SMatthew Knepley f = <2x^3 + 4x^2y - 2xy^2 -4\nu + 1, 4xy^2 + 2x^2y - 2y^3 -4\nu + 1> 72*649ef022SMatthew Knepley Q = 3x^2 + y^2 - 2xy 73*649ef022SMatthew Knepley 74*649ef022SMatthew Knepley so that 75*649ef022SMatthew Knepley 76*649ef022SMatthew Knepley (1) \nabla \cdot u = 2x - 2x = 0 77*649ef022SMatthew Knepley 78*649ef022SMatthew Knepley (2) u \cdot \nabla u - \nu \Delta u + \nabla p - f 79*649ef022SMatthew Knepley = <2x^3 + 4x^2y -2xy^2, 4xy^2 + 2x^2y - 2y^3> -\nu <4, 4> + <1, 1> - <2x^3 + 4x^2y - 2xy^2 -4\nu + 1, 4xy^2 + 2x^2y - 2y^3 - 4\nu + 1> = 0 80*649ef022SMatthew Knepley 81*649ef022SMatthew Knepley (3) u \cdot \nabla T - \alpha \Delta T - Q = 3x^2 + y^2 - 2xy - \alpha*0 - 3x^2 - y^2 + 2xy = 0 82*649ef022SMatthew Knepley */ 83*649ef022SMatthew Knepley 84*649ef022SMatthew Knepley static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 85*649ef022SMatthew Knepley { 86*649ef022SMatthew Knepley u[0] = X[0]*X[0] + X[1]*X[1]; 87*649ef022SMatthew Knepley u[1] = 2.0*X[0]*X[0] - 2.0*X[0]*X[1]; 88*649ef022SMatthew Knepley return 0; 89*649ef022SMatthew Knepley } 90*649ef022SMatthew Knepley 91*649ef022SMatthew Knepley static PetscErrorCode linear_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 92*649ef022SMatthew Knepley { 93*649ef022SMatthew Knepley p[0] = X[0] + X[1] - 1.0; 94*649ef022SMatthew Knepley return 0; 95*649ef022SMatthew Knepley } 96*649ef022SMatthew Knepley 97*649ef022SMatthew Knepley static PetscErrorCode linear_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 98*649ef022SMatthew Knepley { 99*649ef022SMatthew Knepley T[0] = X[0] + X[1]; 100*649ef022SMatthew Knepley return 0; 101*649ef022SMatthew Knepley } 102*649ef022SMatthew Knepley 103*649ef022SMatthew Knepley static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 104*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 105*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 106*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 107*649ef022SMatthew Knepley { 108*649ef022SMatthew Knepley PetscInt c, d; 109*649ef022SMatthew Knepley PetscInt Nc = dim; 110*649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 111*649ef022SMatthew Knepley 112*649ef022SMatthew Knepley for (c=0; c<Nc; ++c) { 113*649ef022SMatthew Knepley for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 114*649ef022SMatthew Knepley } 115*649ef022SMatthew Knepley f0[0] -= (2*X[0]*X[0]*X[0] + 4*X[0]*X[0]*X[1] - 2*X[0]*X[1]*X[1] - 4.0*nu + 1); 116*649ef022SMatthew Knepley f0[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 + 1); 117*649ef022SMatthew Knepley } 118*649ef022SMatthew Knepley 119*649ef022SMatthew Knepley static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 120*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 121*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 122*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 123*649ef022SMatthew Knepley { 124*649ef022SMatthew Knepley PetscInt d; 125*649ef022SMatthew Knepley for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 126*649ef022SMatthew Knepley f0[0] -= (3*X[0]*X[0] + X[1]*X[1] - 2*X[0]*X[1]); 127*649ef022SMatthew Knepley } 128*649ef022SMatthew Knepley 129*649ef022SMatthew Knepley 130*649ef022SMatthew Knepley /* 131*649ef022SMatthew Knepley CASE: cubic 132*649ef022SMatthew Knepley In 2D we use exact solution: 133*649ef022SMatthew Knepley 134*649ef022SMatthew Knepley u = x^3 + y^3 135*649ef022SMatthew Knepley v = 2x^3 - 3x^2y 136*649ef022SMatthew Knepley p = 3/2 x^2 + 3/2 y^2 - 1 137*649ef022SMatthew Knepley T = 1/2 x^2 + 1/2 y^2 138*649ef022SMatthew Knepley f = <3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y), 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y> 139*649ef022SMatthew Knepley Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2 140*649ef022SMatthew Knepley 141*649ef022SMatthew Knepley so that 142*649ef022SMatthew Knepley 143*649ef022SMatthew Knepley \nabla \cdot u = 3x^2 - 3x^2 = 0 144*649ef022SMatthew Knepley 145*649ef022SMatthew Knepley u \cdot \nabla u - \nu \Delta u + \nabla p - f 146*649ef022SMatthew Knepley = <3x^5 + 6x^3y^2 - 6x^2y^3, 6x^2y^3 + 3x^4y - 6xy^4> - \nu<6x + 6y, 12x - 6y> + <3x, 3y> - <3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y), 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y> = 0 147*649ef022SMatthew Knepley 148*649ef022SMatthew Knepley u \cdot \nabla T - \alpha\Delta T - Q = (x^3 + y^3) x + (2x^3 - 3x^2y) y - 2*\alpha - (x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2) = 0 149*649ef022SMatthew Knepley */ 150*649ef022SMatthew Knepley 151*649ef022SMatthew Knepley static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 152*649ef022SMatthew Knepley { 153*649ef022SMatthew Knepley u[0] = X[0]*X[0]*X[0] + X[1]*X[1]*X[1]; 154*649ef022SMatthew Knepley u[1] = 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1]; 155*649ef022SMatthew Knepley return 0; 156*649ef022SMatthew Knepley } 157*649ef022SMatthew Knepley 158*649ef022SMatthew Knepley static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 159*649ef022SMatthew Knepley { 160*649ef022SMatthew Knepley p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0; 161*649ef022SMatthew Knepley return 0; 162*649ef022SMatthew Knepley } 163*649ef022SMatthew Knepley 164*649ef022SMatthew Knepley static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 165*649ef022SMatthew Knepley { 166*649ef022SMatthew Knepley T[0] = X[0]*X[0]/2.0 + X[1]*X[1]/2.0; 167*649ef022SMatthew Knepley return 0; 168*649ef022SMatthew Knepley } 169*649ef022SMatthew Knepley 170*649ef022SMatthew Knepley 171*649ef022SMatthew Knepley static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 172*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 173*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 174*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 175*649ef022SMatthew Knepley { 176*649ef022SMatthew Knepley PetscInt c, d; 177*649ef022SMatthew Knepley PetscInt Nc = dim; 178*649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 179*649ef022SMatthew Knepley 180*649ef022SMatthew Knepley for (c=0; c<Nc; ++c) { 181*649ef022SMatthew Knepley for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 182*649ef022SMatthew Knepley } 183*649ef022SMatthew Knepley f0[0] -= (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]); 184*649ef022SMatthew Knepley f0[1] -= (6*X[0]*X[0]*X[1]*X[1]*X[1] + 3*X[0]*X[0]*X[0]*X[0]*X[1] - 6*X[0]*X[1]*X[1]*X[1]*X[1] - (12*X[0] - 6*X[1])*nu + 3*X[1]); 185*649ef022SMatthew Knepley } 186*649ef022SMatthew Knepley 187*649ef022SMatthew Knepley static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 188*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 189*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 190*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 191*649ef022SMatthew Knepley { 192*649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 193*649ef022SMatthew Knepley PetscInt d; 194*649ef022SMatthew Knepley 195*649ef022SMatthew Knepley for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 196*649ef022SMatthew Knepley f0[0] -= (X[0]*X[0]*X[0]*X[0] + X[0]*X[1]*X[1]*X[1] + 2.0*X[0]*X[0]*X[0]*X[1] - 3.0*X[0]*X[0]*X[1]*X[1] - 2.0*alpha); 197*649ef022SMatthew Knepley } 198*649ef022SMatthew Knepley 199*649ef022SMatthew Knepley static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 200*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 201*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 202*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 203*649ef022SMatthew Knepley { 204*649ef022SMatthew Knepley PetscInt d; 205*649ef022SMatthew Knepley for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; 206*649ef022SMatthew Knepley } 207*649ef022SMatthew Knepley 208*649ef022SMatthew Knepley static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 209*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 210*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 211*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 212*649ef022SMatthew Knepley { 213*649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 214*649ef022SMatthew Knepley const PetscInt Nc = dim; 215*649ef022SMatthew Knepley PetscInt c, d; 216*649ef022SMatthew Knepley 217*649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 218*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 219*649ef022SMatthew Knepley f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]); 220*649ef022SMatthew Knepley //f1[c*dim+d] = nu*u_x[c*dim+d]; 221*649ef022SMatthew Knepley } 222*649ef022SMatthew Knepley f1[c*dim+c] -= u[uOff[1]]; 223*649ef022SMatthew Knepley } 224*649ef022SMatthew Knepley } 225*649ef022SMatthew Knepley 226*649ef022SMatthew Knepley static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 227*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 228*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 229*649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 230*649ef022SMatthew Knepley { 231*649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 232*649ef022SMatthew Knepley PetscInt d; 233*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d]; 234*649ef022SMatthew Knepley } 235*649ef022SMatthew Knepley 236*649ef022SMatthew Knepley /*Jacobians*/ 237*649ef022SMatthew Knepley static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 238*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 239*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 240*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 241*649ef022SMatthew Knepley { 242*649ef022SMatthew Knepley PetscInt d; 243*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; 244*649ef022SMatthew Knepley } 245*649ef022SMatthew Knepley 246*649ef022SMatthew Knepley static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 247*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 248*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 249*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 250*649ef022SMatthew Knepley { 251*649ef022SMatthew Knepley const PetscInt Nc = dim; 252*649ef022SMatthew Knepley PetscInt c, d; 253*649ef022SMatthew Knepley 254*649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 255*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 256*649ef022SMatthew Knepley g0[c*Nc+d] = u_x[ c*Nc+d]; 257*649ef022SMatthew Knepley } 258*649ef022SMatthew Knepley } 259*649ef022SMatthew Knepley } 260*649ef022SMatthew Knepley 261*649ef022SMatthew Knepley static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 262*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 263*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 264*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 265*649ef022SMatthew Knepley { 266*649ef022SMatthew Knepley PetscInt NcI = dim; 267*649ef022SMatthew Knepley PetscInt NcJ = dim; 268*649ef022SMatthew Knepley PetscInt c, d, e; 269*649ef022SMatthew Knepley 270*649ef022SMatthew Knepley for (c = 0; c < NcI; ++c) { 271*649ef022SMatthew Knepley for (d = 0; d < NcJ; ++d) { 272*649ef022SMatthew Knepley for (e = 0; e < dim; ++e) { 273*649ef022SMatthew Knepley if (c == d) { 274*649ef022SMatthew Knepley g1[(c*NcJ+d)*dim+e] = u[e]; 275*649ef022SMatthew Knepley } 276*649ef022SMatthew Knepley } 277*649ef022SMatthew Knepley } 278*649ef022SMatthew Knepley } 279*649ef022SMatthew Knepley } 280*649ef022SMatthew Knepley 281*649ef022SMatthew Knepley static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 282*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 283*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 284*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 285*649ef022SMatthew Knepley { 286*649ef022SMatthew Knepley PetscInt d; 287*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; 288*649ef022SMatthew Knepley } 289*649ef022SMatthew Knepley 290*649ef022SMatthew Knepley static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 291*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 292*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 293*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 294*649ef022SMatthew Knepley { 295*649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 296*649ef022SMatthew Knepley const PetscInt Nc = dim; 297*649ef022SMatthew Knepley PetscInt c, d; 298*649ef022SMatthew Knepley 299*649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 300*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 301*649ef022SMatthew Knepley g3[((c*Nc+c)*dim+d)*dim+d] += nu; // gradU 302*649ef022SMatthew Knepley g3[((c*Nc+d)*dim+d)*dim+c] += nu; // gradU transpose 303*649ef022SMatthew Knepley } 304*649ef022SMatthew Knepley } 305*649ef022SMatthew Knepley } 306*649ef022SMatthew Knepley 307*649ef022SMatthew Knepley static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 308*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 309*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 310*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 311*649ef022SMatthew Knepley { 312*649ef022SMatthew Knepley PetscInt d; 313*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d]; 314*649ef022SMatthew Knepley } 315*649ef022SMatthew Knepley 316*649ef022SMatthew Knepley 317*649ef022SMatthew Knepley static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 318*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 319*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 320*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 321*649ef022SMatthew Knepley { 322*649ef022SMatthew Knepley PetscInt d; 323*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d]; 324*649ef022SMatthew Knepley } 325*649ef022SMatthew Knepley 326*649ef022SMatthew Knepley static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 327*649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 328*649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 329*649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 330*649ef022SMatthew Knepley { 331*649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 332*649ef022SMatthew Knepley PetscInt d; 333*649ef022SMatthew Knepley 334*649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha; 335*649ef022SMatthew Knepley } 336*649ef022SMatthew Knepley 337*649ef022SMatthew Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 338*649ef022SMatthew Knepley { 339*649ef022SMatthew Knepley PetscInt n = 3, sol; 340*649ef022SMatthew Knepley PetscErrorCode ierr; 341*649ef022SMatthew Knepley 342*649ef022SMatthew Knepley 343*649ef022SMatthew Knepley PetscFunctionBeginUser; 344*649ef022SMatthew Knepley options->dim = 2; 345*649ef022SMatthew Knepley options->simplex = PETSC_TRUE; 346*649ef022SMatthew Knepley options->cells[0] = 3; 347*649ef022SMatthew Knepley options->cells[1] = 3; 348*649ef022SMatthew Knepley options->cells[2] = 3; 349*649ef022SMatthew Knepley options->solType = SOL_QUADRATIC; 350*649ef022SMatthew Knepley options->showError= PETSC_FALSE; 351*649ef022SMatthew Knepley 352*649ef022SMatthew Knepley ierr = PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");CHKERRQ(ierr); 353*649ef022SMatthew Knepley ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex62.c", options->dim, &options->dim, NULL);CHKERRQ(ierr); 354*649ef022SMatthew Knepley ierr = PetscOptionsBool("-simplex", "Use simplices or tensor product cells", "ex62.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 355*649ef022SMatthew Knepley ierr = PetscOptionsIntArray("-cells", "The initial mesh division", "ex62.c", options->cells, &n, NULL);CHKERRQ(ierr); 356*649ef022SMatthew Knepley if (options->simplex) { 357*649ef022SMatthew Knepley options->cells[0] = 4 - options->dim; 358*649ef022SMatthew Knepley options->cells[1] = 4 - options->dim; 359*649ef022SMatthew Knepley options->cells[2] = 4 - options->dim; 360*649ef022SMatthew Knepley } 361*649ef022SMatthew Knepley sol = options->solType; 362*649ef022SMatthew Knepley ierr = PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr); 363*649ef022SMatthew Knepley options->solType = (SolType) sol; 364*649ef022SMatthew Knepley ierr = PetscOptionsBool("-show_error", "Output the error for verification", "ex62.c", options->showError, &options->showError, NULL);CHKERRQ(ierr); 365*649ef022SMatthew Knepley 366*649ef022SMatthew Knepley ierr = PetscOptionsEnd(); 367*649ef022SMatthew Knepley PetscFunctionReturn(0); 368*649ef022SMatthew Knepley } 369*649ef022SMatthew Knepley 370*649ef022SMatthew Knepley static PetscErrorCode SetupParameters(AppCtx *user) 371*649ef022SMatthew Knepley { 372*649ef022SMatthew Knepley PetscBag bag; 373*649ef022SMatthew Knepley Parameter *p; 374*649ef022SMatthew Knepley PetscErrorCode ierr; 375*649ef022SMatthew Knepley 376*649ef022SMatthew Knepley PetscFunctionBeginUser; 377*649ef022SMatthew Knepley /* setup PETSc parameter bag */ 378*649ef022SMatthew Knepley ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr); 379*649ef022SMatthew Knepley ierr = PetscBagSetName(user->bag, "par", "Poiseuille flow parameters");CHKERRQ(ierr); 380*649ef022SMatthew Knepley bag = user->bag; 381*649ef022SMatthew Knepley ierr = PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity");CHKERRQ(ierr); 382*649ef022SMatthew Knepley ierr = PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity");CHKERRQ(ierr); 383*649ef022SMatthew Knepley ierr = PetscBagRegisterReal(bag, &p->theta, 0.0, "theta", "Angle of pipe wall to x-axis");CHKERRQ(ierr); 384*649ef022SMatthew Knepley PetscFunctionReturn(0); 385*649ef022SMatthew Knepley } 386*649ef022SMatthew Knepley 387*649ef022SMatthew Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 388*649ef022SMatthew Knepley { 389*649ef022SMatthew Knepley PetscInt dim = user->dim; 390*649ef022SMatthew Knepley PetscErrorCode ierr; 391*649ef022SMatthew Knepley 392*649ef022SMatthew Knepley PetscFunctionBeginUser; 393*649ef022SMatthew Knepley ierr = DMPlexCreateBoxMesh(comm, dim, user->simplex, user->cells, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 394*649ef022SMatthew Knepley { 395*649ef022SMatthew Knepley Parameter *param; 396*649ef022SMatthew Knepley Vec coordinates; 397*649ef022SMatthew Knepley PetscScalar *coords; 398*649ef022SMatthew Knepley PetscReal theta; 399*649ef022SMatthew Knepley PetscInt cdim, N, bs, i; 400*649ef022SMatthew Knepley 401*649ef022SMatthew Knepley ierr = DMGetCoordinateDim(*dm, &cdim);CHKERRQ(ierr); 402*649ef022SMatthew Knepley ierr = DMGetCoordinates(*dm, &coordinates);CHKERRQ(ierr); 403*649ef022SMatthew Knepley ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr); 404*649ef022SMatthew Knepley ierr = VecGetBlockSize(coordinates, &bs);CHKERRQ(ierr); 405*649ef022SMatthew Knepley if (bs != cdim) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %D != embedding dimension %D", bs, cdim); 406*649ef022SMatthew Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 407*649ef022SMatthew Knepley ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 408*649ef022SMatthew Knepley theta = param->theta; 409*649ef022SMatthew Knepley for (i = 0; i < N; i += cdim) { 410*649ef022SMatthew Knepley PetscScalar x = coords[i+0]; 411*649ef022SMatthew Knepley PetscScalar y = coords[i+1]; 412*649ef022SMatthew Knepley 413*649ef022SMatthew Knepley coords[i+0] = PetscCosReal(theta)*x - PetscSinReal(theta)*y; 414*649ef022SMatthew Knepley coords[i+1] = PetscSinReal(theta)*x + PetscCosReal(theta)*y; 415*649ef022SMatthew Knepley } 416*649ef022SMatthew Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 417*649ef022SMatthew Knepley ierr = DMSetCoordinates(*dm, coordinates);CHKERRQ(ierr); 418*649ef022SMatthew Knepley } 419*649ef022SMatthew Knepley { 420*649ef022SMatthew Knepley DM pdm = NULL; 421*649ef022SMatthew Knepley PetscPartitioner part; 422*649ef022SMatthew Knepley 423*649ef022SMatthew Knepley ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 424*649ef022SMatthew Knepley ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 425*649ef022SMatthew Knepley ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr); 426*649ef022SMatthew Knepley if (pdm) { 427*649ef022SMatthew Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 428*649ef022SMatthew Knepley *dm = pdm; 429*649ef022SMatthew Knepley } 430*649ef022SMatthew Knepley } 431*649ef022SMatthew Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 432*649ef022SMatthew Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 433*649ef022SMatthew Knepley PetscFunctionReturn(0); 434*649ef022SMatthew Knepley } 435*649ef022SMatthew Knepley 436*649ef022SMatthew Knepley static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 437*649ef022SMatthew Knepley { 438*649ef022SMatthew Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 439*649ef022SMatthew Knepley PetscDS prob; 440*649ef022SMatthew Knepley Parameter *ctx; 441*649ef022SMatthew Knepley PetscInt id; 442*649ef022SMatthew Knepley PetscErrorCode ierr; 443*649ef022SMatthew Knepley 444*649ef022SMatthew Knepley PetscFunctionBeginUser; 445*649ef022SMatthew Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 446*649ef022SMatthew Knepley switch(user->solType){ 447*649ef022SMatthew Knepley case SOL_QUADRATIC: 448*649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v);CHKERRQ(ierr); 449*649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w);CHKERRQ(ierr); 450*649ef022SMatthew Knepley 451*649ef022SMatthew Knepley exactFuncs[0] = quadratic_u; 452*649ef022SMatthew Knepley exactFuncs[1] = linear_p; 453*649ef022SMatthew Knepley exactFuncs[2] = linear_T; 454*649ef022SMatthew Knepley break; 455*649ef022SMatthew Knepley case SOL_CUBIC: 456*649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v);CHKERRQ(ierr); 457*649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w);CHKERRQ(ierr); 458*649ef022SMatthew Knepley 459*649ef022SMatthew Knepley exactFuncs[0] = cubic_u; 460*649ef022SMatthew Knepley exactFuncs[1] = quadratic_p; 461*649ef022SMatthew Knepley exactFuncs[2] = quadratic_T; 462*649ef022SMatthew Knepley break; 463*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); 464*649ef022SMatthew Knepley } 465*649ef022SMatthew Knepley 466*649ef022SMatthew Knepley ierr = PetscDSSetResidual(prob, 1, f0_q, NULL);CHKERRQ(ierr); 467*649ef022SMatthew Knepley 468*649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu);CHKERRQ(ierr); 469*649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL);CHKERRQ(ierr); 470*649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL);CHKERRQ(ierr); 471*649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL);CHKERRQ(ierr); 472*649ef022SMatthew Knepley ierr = PetscDSSetJacobian(prob, 2, 2, NULL, g1_wT, NULL, g3_wT);CHKERRQ(ierr); 473*649ef022SMatthew Knepley /* Setup constants */ 474*649ef022SMatthew Knepley { 475*649ef022SMatthew Knepley Parameter *param; 476*649ef022SMatthew Knepley PetscScalar constants[3]; 477*649ef022SMatthew Knepley 478*649ef022SMatthew Knepley ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 479*649ef022SMatthew Knepley 480*649ef022SMatthew Knepley constants[0] = param->nu; 481*649ef022SMatthew Knepley constants[1] = param->alpha; 482*649ef022SMatthew Knepley constants[2] = param->theta; 483*649ef022SMatthew Knepley ierr = PetscDSSetConstants(prob, 3, constants);CHKERRQ(ierr); 484*649ef022SMatthew Knepley } 485*649ef022SMatthew Knepley /* Setup Boundary Conditions */ 486*649ef022SMatthew Knepley ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr); 487*649ef022SMatthew Knepley id = 3; 488*649ef022SMatthew Knepley ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity", "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, 1, &id, ctx);CHKERRQ(ierr); 489*649ef022SMatthew Knepley id = 1; 490*649ef022SMatthew Knepley ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, 1, &id, ctx);CHKERRQ(ierr); 491*649ef022SMatthew Knepley id = 2; 492*649ef022SMatthew Knepley ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity", "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, 1, &id, ctx);CHKERRQ(ierr); 493*649ef022SMatthew Knepley id = 4; 494*649ef022SMatthew Knepley ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity", "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, 1, &id, ctx);CHKERRQ(ierr); 495*649ef022SMatthew Knepley id = 3; 496*649ef022SMatthew Knepley ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp", "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, 1, &id, ctx);CHKERRQ(ierr); 497*649ef022SMatthew Knepley id = 1; 498*649ef022SMatthew Knepley ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, 1, &id, ctx);CHKERRQ(ierr); 499*649ef022SMatthew Knepley id = 2; 500*649ef022SMatthew Knepley ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp", "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, 1, &id, ctx);CHKERRQ(ierr); 501*649ef022SMatthew Knepley id = 4; 502*649ef022SMatthew Knepley ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp", "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, 1, &id, ctx);CHKERRQ(ierr); 503*649ef022SMatthew Knepley 504*649ef022SMatthew Knepley /*setup exact solution.*/ 505*649ef022SMatthew Knepley ierr = PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx);CHKERRQ(ierr); 506*649ef022SMatthew Knepley ierr = PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx);CHKERRQ(ierr); 507*649ef022SMatthew Knepley ierr = PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx);CHKERRQ(ierr); 508*649ef022SMatthew Knepley PetscFunctionReturn(0); 509*649ef022SMatthew Knepley } 510*649ef022SMatthew Knepley 511*649ef022SMatthew Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 512*649ef022SMatthew Knepley { 513*649ef022SMatthew Knepley DM cdm = dm; 514*649ef022SMatthew Knepley const PetscInt dim = user->dim; 515*649ef022SMatthew Knepley PetscFE fe[3]; 516*649ef022SMatthew Knepley Parameter *param; 517*649ef022SMatthew Knepley MPI_Comm comm; 518*649ef022SMatthew Knepley PetscErrorCode ierr; 519*649ef022SMatthew Knepley 520*649ef022SMatthew Knepley PetscFunctionBeginUser; 521*649ef022SMatthew Knepley /* Create finite element */ 522*649ef022SMatthew Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 523*649ef022SMatthew Knepley ierr = PetscFECreateDefault(comm, dim, dim, user->simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 524*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr); 525*649ef022SMatthew Knepley 526*649ef022SMatthew Knepley ierr = PetscFECreateDefault(comm, dim, 1, user->simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 527*649ef022SMatthew Knepley ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 528*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr); 529*649ef022SMatthew Knepley 530*649ef022SMatthew Knepley ierr = PetscFECreateDefault(comm, dim, 1, user->simplex, "temp_", PETSC_DEFAULT, &fe[2]);CHKERRQ(ierr); 531*649ef022SMatthew Knepley ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr); 532*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) fe[2], "temperature");CHKERRQ(ierr); 533*649ef022SMatthew Knepley 534*649ef022SMatthew Knepley /* Set discretization and boundary conditions for each mesh */ 535*649ef022SMatthew Knepley ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 536*649ef022SMatthew Knepley ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr); 537*649ef022SMatthew Knepley ierr = DMSetField(dm, 2, NULL, (PetscObject) fe[2]);CHKERRQ(ierr); 538*649ef022SMatthew Knepley ierr = DMCreateDS(dm);CHKERRQ(ierr); 539*649ef022SMatthew Knepley ierr = SetupProblem(dm, user);CHKERRQ(ierr); 540*649ef022SMatthew Knepley ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 541*649ef022SMatthew Knepley while (cdm) { 542*649ef022SMatthew Knepley ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 543*649ef022SMatthew Knepley ierr = DMPlexCreateBasisRotation(cdm, param->theta, 0.0, 0.0);CHKERRQ(ierr); 544*649ef022SMatthew Knepley ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 545*649ef022SMatthew Knepley } 546*649ef022SMatthew Knepley ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); 547*649ef022SMatthew Knepley ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); 548*649ef022SMatthew Knepley ierr = PetscFEDestroy(&fe[2]);CHKERRQ(ierr); 549*649ef022SMatthew Knepley PetscFunctionReturn(0); 550*649ef022SMatthew Knepley } 551*649ef022SMatthew Knepley 552*649ef022SMatthew Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) 553*649ef022SMatthew Knepley { 554*649ef022SMatthew Knepley Vec vec; 555*649ef022SMatthew Knepley PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 556*649ef022SMatthew Knepley PetscErrorCode ierr; 557*649ef022SMatthew Knepley 558*649ef022SMatthew Knepley PetscFunctionBeginUser; 559*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); 560*649ef022SMatthew Knepley funcs[nfield] = constant; 561*649ef022SMatthew Knepley ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr); 562*649ef022SMatthew Knepley ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr); 563*649ef022SMatthew Knepley ierr = VecNormalize(vec, NULL);CHKERRQ(ierr); 564*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr); 565*649ef022SMatthew Knepley ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr); 566*649ef022SMatthew Knepley ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace);CHKERRQ(ierr); 567*649ef022SMatthew Knepley ierr = VecDestroy(&vec);CHKERRQ(ierr); 568*649ef022SMatthew Knepley PetscFunctionReturn(0); 569*649ef022SMatthew Knepley } 570*649ef022SMatthew Knepley 571*649ef022SMatthew Knepley int main(int argc, char **argv) 572*649ef022SMatthew Knepley { 573*649ef022SMatthew Knepley SNES snes; /* nonlinear solver */ 574*649ef022SMatthew Knepley DM dm; /* problem definition */ 575*649ef022SMatthew Knepley Vec u, r; /* solution, residual vectors */ 576*649ef022SMatthew Knepley AppCtx user; /* user-defined work context */ 577*649ef022SMatthew Knepley PetscErrorCode ierr; 578*649ef022SMatthew Knepley 579*649ef022SMatthew Knepley ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 580*649ef022SMatthew Knepley ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 581*649ef022SMatthew Knepley ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr); 582*649ef022SMatthew Knepley ierr = SetupParameters(&user);CHKERRQ(ierr); 583*649ef022SMatthew Knepley ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 584*649ef022SMatthew Knepley ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 585*649ef022SMatthew Knepley ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 586*649ef022SMatthew Knepley ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr); 587*649ef022SMatthew Knepley /* Setup problem */ 588*649ef022SMatthew Knepley ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr); 589*649ef022SMatthew Knepley ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr); 590*649ef022SMatthew Knepley 591*649ef022SMatthew Knepley ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 592*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) u, "Solution");CHKERRQ(ierr); 593*649ef022SMatthew Knepley ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 594*649ef022SMatthew Knepley 595*649ef022SMatthew Knepley ierr = DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace);CHKERRQ(ierr); 596*649ef022SMatthew Knepley ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr); 597*649ef022SMatthew Knepley 598*649ef022SMatthew Knepley ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 599*649ef022SMatthew Knepley { 600*649ef022SMatthew Knepley PetscDS ds; 601*649ef022SMatthew Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 602*649ef022SMatthew Knepley void *ctxs[3]; 603*649ef022SMatthew Knepley 604*649ef022SMatthew Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 605*649ef022SMatthew Knepley ierr = PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0]);CHKERRQ(ierr); 606*649ef022SMatthew Knepley ierr = PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1]);CHKERRQ(ierr); 607*649ef022SMatthew Knepley ierr = PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2]);CHKERRQ(ierr); 608*649ef022SMatthew Knepley ierr = DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr); 609*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) u, "Exact Solution");CHKERRQ(ierr); 610*649ef022SMatthew Knepley ierr = VecViewFromOptions(u, NULL, "-exact_vec_view");CHKERRQ(ierr); 611*649ef022SMatthew Knepley } 612*649ef022SMatthew Knepley ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr); 613*649ef022SMatthew Knepley ierr = VecSet(u, 0.0);CHKERRQ(ierr); 614*649ef022SMatthew Knepley ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 615*649ef022SMatthew Knepley 616*649ef022SMatthew Knepley if (user.showError) { 617*649ef022SMatthew Knepley PetscDS ds; 618*649ef022SMatthew Knepley Vec r; 619*649ef022SMatthew Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 620*649ef022SMatthew Knepley void *ctxs[3]; 621*649ef022SMatthew Knepley 622*649ef022SMatthew Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 623*649ef022SMatthew Knepley ierr = PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0]);CHKERRQ(ierr); 624*649ef022SMatthew Knepley ierr = PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1]);CHKERRQ(ierr); 625*649ef022SMatthew Knepley ierr = PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2]);CHKERRQ(ierr); 626*649ef022SMatthew Knepley ierr = DMGetGlobalVector(dm, &r);CHKERRQ(ierr); 627*649ef022SMatthew Knepley ierr = DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, r);CHKERRQ(ierr); 628*649ef022SMatthew Knepley ierr = VecAXPY(r, -1.0, u);CHKERRQ(ierr); 629*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) r, "Solution Error");CHKERRQ(ierr); 630*649ef022SMatthew Knepley ierr = VecViewFromOptions(r, NULL, "-error_vec_view");CHKERRQ(ierr); 631*649ef022SMatthew Knepley ierr = DMRestoreGlobalVector(dm, &r);CHKERRQ(ierr); 632*649ef022SMatthew Knepley } 633*649ef022SMatthew Knepley ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr); 634*649ef022SMatthew Knepley ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 635*649ef022SMatthew Knepley 636*649ef022SMatthew Knepley ierr = VecDestroy(&u);CHKERRQ(ierr); 637*649ef022SMatthew Knepley ierr = VecDestroy(&r);CHKERRQ(ierr); 638*649ef022SMatthew Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 639*649ef022SMatthew Knepley ierr = SNESDestroy(&snes);CHKERRQ(ierr); 640*649ef022SMatthew Knepley ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr); 641*649ef022SMatthew Knepley ierr = PetscFinalize(); 642*649ef022SMatthew Knepley return ierr; 643*649ef022SMatthew Knepley } 644*649ef022SMatthew Knepley 645*649ef022SMatthew Knepley /*TEST 646*649ef022SMatthew Knepley 647*649ef022SMatthew Knepley test: 648*649ef022SMatthew Knepley suffix: 2d_tri_p2_p1_p1 649*649ef022SMatthew Knepley requires: triangle !single 650*649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \ 651*649ef022SMatthew Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 652*649ef022SMatthew Knepley -dmsnes_check .001 -snes_error_if_not_converged \ 653*649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 654*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 \ 655*649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 656*649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 657*649ef022SMatthew Knepley 658*649ef022SMatthew Knepley test: 659*649ef022SMatthew Knepley # Using -dm_refine 2 -convest_num_refine 3 gives L_2 convergence rate: [2.9, 2.3, 1.9] 660*649ef022SMatthew Knepley suffix: 2d_tri_p2_p1_p1_conv 661*649ef022SMatthew Knepley requires: triangle !single 662*649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 663*649ef022SMatthew Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 664*649ef022SMatthew Knepley -snes_error_if_not_converged -snes_convergence_test correct_pressure -snes_convergence_estimate -convest_num_refine 1 \ 665*649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 666*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 \ 667*649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 668*649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 669*649ef022SMatthew Knepley 670*649ef022SMatthew Knepley test: 671*649ef022SMatthew Knepley suffix: 2d_tri_p3_p2_p2 672*649ef022SMatthew Knepley requires: triangle !single 673*649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 674*649ef022SMatthew Knepley -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \ 675*649ef022SMatthew Knepley -dmsnes_check .001 -snes_error_if_not_converged \ 676*649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 677*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 \ 678*649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 679*649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 680*649ef022SMatthew Knepley 681*649ef022SMatthew Knepley TEST*/ 682