1649ef022SMatthew Knepley static char help[] = "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 a steady-state 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, 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 <petscds.h> 25649ef022SMatthew Knepley #include <petscbag.h> 26649ef022SMatthew Knepley 279371c9d4SSatish Balay typedef enum { 289371c9d4SSatish Balay SOL_QUADRATIC, 299371c9d4SSatish Balay SOL_CUBIC, 309371c9d4SSatish Balay NUM_SOL_TYPES 319371c9d4SSatish Balay } SolType; 32649ef022SMatthew Knepley const char *solTypes[NUM_SOL_TYPES + 1] = {"quadratic", "cubic", "unknown"}; 33649ef022SMatthew Knepley 34649ef022SMatthew Knepley typedef struct { 35649ef022SMatthew Knepley PetscReal nu; /* Kinematic viscosity */ 36649ef022SMatthew Knepley PetscReal theta; /* Angle of pipe wall to x-axis */ 37649ef022SMatthew Knepley PetscReal alpha; /* Thermal diffusivity */ 38649ef022SMatthew Knepley PetscReal T_in; /* Inlet temperature*/ 39649ef022SMatthew Knepley } Parameter; 40649ef022SMatthew Knepley 41649ef022SMatthew Knepley typedef struct { 4230602db0SMatthew G. Knepley PetscBool showError; 4330602db0SMatthew G. Knepley PetscBag bag; 44649ef022SMatthew Knepley SolType solType; 45649ef022SMatthew Knepley } AppCtx; 46649ef022SMatthew Knepley 47d71ae5a4SJacob Faibussowitsch static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 48d71ae5a4SJacob Faibussowitsch { 49649ef022SMatthew Knepley PetscInt d; 50649ef022SMatthew Knepley for (d = 0; d < Nc; ++d) u[d] = 0.0; 51*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 52649ef022SMatthew Knepley } 53649ef022SMatthew Knepley 54d71ae5a4SJacob Faibussowitsch static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 55d71ae5a4SJacob Faibussowitsch { 56649ef022SMatthew Knepley PetscInt d; 57649ef022SMatthew Knepley for (d = 0; d < Nc; ++d) u[d] = 1.0; 58*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 59649ef022SMatthew Knepley } 60649ef022SMatthew Knepley 61649ef022SMatthew Knepley /* 62649ef022SMatthew Knepley CASE: quadratic 63649ef022SMatthew Knepley In 2D we use exact solution: 64649ef022SMatthew Knepley 65649ef022SMatthew Knepley u = x^2 + y^2 66649ef022SMatthew Knepley v = 2x^2 - 2xy 67649ef022SMatthew Knepley p = x + y - 1 68649ef022SMatthew Knepley T = x + y 69649ef022SMatthew Knepley f = <2x^3 + 4x^2y - 2xy^2 -4\nu + 1, 4xy^2 + 2x^2y - 2y^3 -4\nu + 1> 70649ef022SMatthew Knepley Q = 3x^2 + y^2 - 2xy 71649ef022SMatthew Knepley 72649ef022SMatthew Knepley so that 73649ef022SMatthew Knepley 74649ef022SMatthew Knepley (1) \nabla \cdot u = 2x - 2x = 0 75649ef022SMatthew Knepley 76649ef022SMatthew Knepley (2) u \cdot \nabla u - \nu \Delta u + \nabla p - f 77649ef022SMatthew 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 78649ef022SMatthew Knepley 79649ef022SMatthew Knepley (3) u \cdot \nabla T - \alpha \Delta T - Q = 3x^2 + y^2 - 2xy - \alpha*0 - 3x^2 - y^2 + 2xy = 0 80649ef022SMatthew Knepley */ 81649ef022SMatthew Knepley 82d71ae5a4SJacob Faibussowitsch static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 83d71ae5a4SJacob Faibussowitsch { 84649ef022SMatthew Knepley u[0] = X[0] * X[0] + X[1] * X[1]; 85649ef022SMatthew Knepley u[1] = 2.0 * X[0] * X[0] - 2.0 * X[0] * X[1]; 86*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 87649ef022SMatthew Knepley } 88649ef022SMatthew Knepley 89d71ae5a4SJacob Faibussowitsch static PetscErrorCode linear_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 90d71ae5a4SJacob Faibussowitsch { 91649ef022SMatthew Knepley p[0] = X[0] + X[1] - 1.0; 92*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 93649ef022SMatthew Knepley } 94649ef022SMatthew Knepley 95d71ae5a4SJacob Faibussowitsch static PetscErrorCode linear_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 96d71ae5a4SJacob Faibussowitsch { 97649ef022SMatthew Knepley T[0] = X[0] + X[1]; 98*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 99649ef022SMatthew Knepley } 100649ef022SMatthew Knepley 101d71ae5a4SJacob Faibussowitsch static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 102d71ae5a4SJacob Faibussowitsch { 103649ef022SMatthew Knepley PetscInt c, d; 104649ef022SMatthew Knepley PetscInt Nc = dim; 105649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 106649ef022SMatthew Knepley 107649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 108649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f0[c] += u[d] * u_x[c * dim + d]; 109649ef022SMatthew Knepley } 110649ef022SMatthew 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); 111649ef022SMatthew 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); 112649ef022SMatthew Knepley } 113649ef022SMatthew Knepley 114d71ae5a4SJacob Faibussowitsch static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 115d71ae5a4SJacob Faibussowitsch { 116649ef022SMatthew Knepley PetscInt d; 117649ef022SMatthew Knepley for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0] + d] * u_x[uOff_x[2] + d]; 118649ef022SMatthew Knepley f0[0] -= (3 * X[0] * X[0] + X[1] * X[1] - 2 * X[0] * X[1]); 119649ef022SMatthew Knepley } 120649ef022SMatthew Knepley 121649ef022SMatthew Knepley /* 122649ef022SMatthew Knepley CASE: cubic 123649ef022SMatthew Knepley In 2D we use exact solution: 124649ef022SMatthew Knepley 125649ef022SMatthew Knepley u = x^3 + y^3 126649ef022SMatthew Knepley v = 2x^3 - 3x^2y 127649ef022SMatthew Knepley p = 3/2 x^2 + 3/2 y^2 - 1 128649ef022SMatthew Knepley T = 1/2 x^2 + 1/2 y^2 129649ef022SMatthew Knepley f = <3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y), 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y> 130649ef022SMatthew Knepley Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2 131649ef022SMatthew Knepley 132649ef022SMatthew Knepley so that 133649ef022SMatthew Knepley 134649ef022SMatthew Knepley \nabla \cdot u = 3x^2 - 3x^2 = 0 135649ef022SMatthew Knepley 136649ef022SMatthew Knepley u \cdot \nabla u - \nu \Delta u + \nabla p - f 137649ef022SMatthew 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 138649ef022SMatthew Knepley 139649ef022SMatthew 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 140649ef022SMatthew Knepley */ 141649ef022SMatthew Knepley 142d71ae5a4SJacob Faibussowitsch static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 143d71ae5a4SJacob Faibussowitsch { 144649ef022SMatthew Knepley u[0] = X[0] * X[0] * X[0] + X[1] * X[1] * X[1]; 145649ef022SMatthew Knepley u[1] = 2.0 * X[0] * X[0] * X[0] - 3.0 * X[0] * X[0] * X[1]; 146*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 147649ef022SMatthew Knepley } 148649ef022SMatthew Knepley 149d71ae5a4SJacob Faibussowitsch static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 150d71ae5a4SJacob Faibussowitsch { 151649ef022SMatthew Knepley p[0] = 3.0 * X[0] * X[0] / 2.0 + 3.0 * X[1] * X[1] / 2.0 - 1.0; 152*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 153649ef022SMatthew Knepley } 154649ef022SMatthew Knepley 155d71ae5a4SJacob Faibussowitsch static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 156d71ae5a4SJacob Faibussowitsch { 157649ef022SMatthew Knepley T[0] = X[0] * X[0] / 2.0 + X[1] * X[1] / 2.0; 158*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 159649ef022SMatthew Knepley } 160649ef022SMatthew Knepley 161d71ae5a4SJacob Faibussowitsch static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 162d71ae5a4SJacob Faibussowitsch { 163649ef022SMatthew Knepley PetscInt c, d; 164649ef022SMatthew Knepley PetscInt Nc = dim; 165649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 166649ef022SMatthew Knepley 167649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 168649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f0[c] += u[d] * u_x[c * dim + d]; 169649ef022SMatthew Knepley } 170649ef022SMatthew 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]); 171649ef022SMatthew 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]); 172649ef022SMatthew Knepley } 173649ef022SMatthew Knepley 174d71ae5a4SJacob Faibussowitsch static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 175d71ae5a4SJacob Faibussowitsch { 176649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 177649ef022SMatthew Knepley PetscInt d; 178649ef022SMatthew Knepley 179649ef022SMatthew Knepley for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0] + d] * u_x[uOff_x[2] + d]; 180649ef022SMatthew 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); 181649ef022SMatthew Knepley } 182649ef022SMatthew Knepley 183d71ae5a4SJacob Faibussowitsch static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 184d71ae5a4SJacob Faibussowitsch { 185649ef022SMatthew Knepley PetscInt d; 186649ef022SMatthew Knepley for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d * dim + d]; 187649ef022SMatthew Knepley } 188649ef022SMatthew Knepley 189d71ae5a4SJacob Faibussowitsch static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 190d71ae5a4SJacob Faibussowitsch { 191649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 192649ef022SMatthew Knepley const PetscInt Nc = dim; 193649ef022SMatthew Knepley PetscInt c, d; 194649ef022SMatthew Knepley 195649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 196649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 197649ef022SMatthew Knepley f1[c * dim + d] = nu * (u_x[c * dim + d] + u_x[d * dim + c]); 198649ef022SMatthew Knepley //f1[c*dim+d] = nu*u_x[c*dim+d]; 199649ef022SMatthew Knepley } 200649ef022SMatthew Knepley f1[c * dim + c] -= u[uOff[1]]; 201649ef022SMatthew Knepley } 202649ef022SMatthew Knepley } 203649ef022SMatthew Knepley 204d71ae5a4SJacob Faibussowitsch static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 205d71ae5a4SJacob Faibussowitsch { 206649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 207649ef022SMatthew Knepley PetscInt d; 208649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f1[d] = alpha * u_x[uOff_x[2] + d]; 209649ef022SMatthew Knepley } 210649ef022SMatthew Knepley 211649ef022SMatthew Knepley /* Jacobians */ 212d71ae5a4SJacob Faibussowitsch static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 213d71ae5a4SJacob Faibussowitsch { 214649ef022SMatthew Knepley PetscInt d; 215649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; 216649ef022SMatthew Knepley } 217649ef022SMatthew Knepley 218d71ae5a4SJacob Faibussowitsch static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 219d71ae5a4SJacob Faibussowitsch { 220649ef022SMatthew Knepley const PetscInt Nc = dim; 221649ef022SMatthew Knepley PetscInt c, d; 222649ef022SMatthew Knepley 223649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 224ad540459SPierre Jolivet for (d = 0; d < dim; ++d) g0[c * Nc + d] = u_x[c * Nc + d]; 225649ef022SMatthew Knepley } 226649ef022SMatthew Knepley } 227649ef022SMatthew Knepley 228d71ae5a4SJacob Faibussowitsch static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 229d71ae5a4SJacob Faibussowitsch { 230649ef022SMatthew Knepley PetscInt NcI = dim; 231649ef022SMatthew Knepley PetscInt NcJ = dim; 232649ef022SMatthew Knepley PetscInt c, d, e; 233649ef022SMatthew Knepley 234649ef022SMatthew Knepley for (c = 0; c < NcI; ++c) { 235649ef022SMatthew Knepley for (d = 0; d < NcJ; ++d) { 236649ef022SMatthew Knepley for (e = 0; e < dim; ++e) { 237ad540459SPierre Jolivet if (c == d) g1[(c * NcJ + d) * dim + e] = u[e]; 238649ef022SMatthew Knepley } 239649ef022SMatthew Knepley } 240649ef022SMatthew Knepley } 241649ef022SMatthew Knepley } 242649ef022SMatthew Knepley 243d71ae5a4SJacob Faibussowitsch static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 244d71ae5a4SJacob Faibussowitsch { 245649ef022SMatthew Knepley PetscInt d; 246649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0; 247649ef022SMatthew Knepley } 248649ef022SMatthew Knepley 249d71ae5a4SJacob Faibussowitsch static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 250d71ae5a4SJacob Faibussowitsch { 251649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 252649ef022SMatthew Knepley const PetscInt Nc = dim; 253649ef022SMatthew Knepley PetscInt c, d; 254649ef022SMatthew Knepley 255649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 256649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 257649ef022SMatthew Knepley g3[((c * Nc + c) * dim + d) * dim + d] += nu; // gradU 258649ef022SMatthew Knepley g3[((c * Nc + d) * dim + d) * dim + c] += nu; // gradU transpose 259649ef022SMatthew Knepley } 260649ef022SMatthew Knepley } 261649ef022SMatthew Knepley } 262649ef022SMatthew Knepley 263d71ae5a4SJacob Faibussowitsch static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 264d71ae5a4SJacob Faibussowitsch { 265649ef022SMatthew Knepley PetscInt d; 266649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2] + d]; 267649ef022SMatthew Knepley } 268649ef022SMatthew Knepley 269d71ae5a4SJacob Faibussowitsch static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 270d71ae5a4SJacob Faibussowitsch { 271649ef022SMatthew Knepley PetscInt d; 272649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g1[d] = u[uOff[0] + d]; 273649ef022SMatthew Knepley } 274649ef022SMatthew Knepley 275d71ae5a4SJacob Faibussowitsch static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 276d71ae5a4SJacob Faibussowitsch { 277649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 278649ef022SMatthew Knepley PetscInt d; 279649ef022SMatthew Knepley 280649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g3[d * dim + d] = alpha; 281649ef022SMatthew Knepley } 282649ef022SMatthew Knepley 283d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 284d71ae5a4SJacob Faibussowitsch { 28530602db0SMatthew G. Knepley PetscInt sol; 286649ef022SMatthew Knepley 287649ef022SMatthew Knepley PetscFunctionBeginUser; 288649ef022SMatthew Knepley options->solType = SOL_QUADRATIC; 289649ef022SMatthew Knepley options->showError = PETSC_FALSE; 290d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX"); 291649ef022SMatthew Knepley sol = options->solType; 2929566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL)); 293649ef022SMatthew Knepley options->solType = (SolType)sol; 2949566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-show_error", "Output the error for verification", "ex62.c", options->showError, &options->showError, NULL)); 295d0609cedSBarry Smith PetscOptionsEnd(); 296*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 297649ef022SMatthew Knepley } 298649ef022SMatthew Knepley 299d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupParameters(AppCtx *user) 300d71ae5a4SJacob Faibussowitsch { 301649ef022SMatthew Knepley PetscBag bag; 302649ef022SMatthew Knepley Parameter *p; 303649ef022SMatthew Knepley 304649ef022SMatthew Knepley PetscFunctionBeginUser; 305649ef022SMatthew Knepley /* setup PETSc parameter bag */ 3069566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)&p)); 3079566063dSJacob Faibussowitsch PetscCall(PetscBagSetName(user->bag, "par", "Poiseuille flow parameters")); 308649ef022SMatthew Knepley bag = user->bag; 3099566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity")); 3109566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity")); 3119566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->theta, 0.0, "theta", "Angle of pipe wall to x-axis")); 312*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 313649ef022SMatthew Knepley } 314649ef022SMatthew Knepley 315d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 316d71ae5a4SJacob Faibussowitsch { 317649ef022SMatthew Knepley PetscFunctionBeginUser; 3189566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 3199566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 3209566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 321649ef022SMatthew Knepley { 322649ef022SMatthew Knepley Parameter *param; 323649ef022SMatthew Knepley Vec coordinates; 324649ef022SMatthew Knepley PetscScalar *coords; 325649ef022SMatthew Knepley PetscReal theta; 326649ef022SMatthew Knepley PetscInt cdim, N, bs, i; 327649ef022SMatthew Knepley 3289566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(*dm, &cdim)); 3299566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(*dm, &coordinates)); 3309566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coordinates, &N)); 3319566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coordinates, &bs)); 33263a3b9bcSJacob Faibussowitsch PetscCheck(bs == cdim, comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %" PetscInt_FMT " != embedding dimension %" PetscInt_FMT, bs, cdim); 3339566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 3349566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 335649ef022SMatthew Knepley theta = param->theta; 336649ef022SMatthew Knepley for (i = 0; i < N; i += cdim) { 337649ef022SMatthew Knepley PetscScalar x = coords[i + 0]; 338649ef022SMatthew Knepley PetscScalar y = coords[i + 1]; 339649ef022SMatthew Knepley 340649ef022SMatthew Knepley coords[i + 0] = PetscCosReal(theta) * x - PetscSinReal(theta) * y; 341649ef022SMatthew Knepley coords[i + 1] = PetscSinReal(theta) * x + PetscCosReal(theta) * y; 342649ef022SMatthew Knepley } 3439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 3449566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(*dm, coordinates)); 345649ef022SMatthew Knepley } 3469566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 347*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 348649ef022SMatthew Knepley } 349649ef022SMatthew Knepley 350d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 351d71ae5a4SJacob Faibussowitsch { 352649ef022SMatthew Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 353649ef022SMatthew Knepley PetscDS prob; 35445480ffeSMatthew G. Knepley DMLabel label; 355649ef022SMatthew Knepley Parameter *ctx; 356649ef022SMatthew Knepley PetscInt id; 357649ef022SMatthew Knepley 358649ef022SMatthew Knepley PetscFunctionBeginUser; 3599566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 3609566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 361649ef022SMatthew Knepley switch (user->solType) { 362649ef022SMatthew Knepley case SOL_QUADRATIC: 3639566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v)); 3649566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w)); 365649ef022SMatthew Knepley 366649ef022SMatthew Knepley exactFuncs[0] = quadratic_u; 367649ef022SMatthew Knepley exactFuncs[1] = linear_p; 368649ef022SMatthew Knepley exactFuncs[2] = linear_T; 369649ef022SMatthew Knepley break; 370649ef022SMatthew Knepley case SOL_CUBIC: 3719566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v)); 3729566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w)); 373649ef022SMatthew Knepley 374649ef022SMatthew Knepley exactFuncs[0] = cubic_u; 375649ef022SMatthew Knepley exactFuncs[1] = quadratic_p; 376649ef022SMatthew Knepley exactFuncs[2] = quadratic_T; 377649ef022SMatthew Knepley break; 378d71ae5a4SJacob Faibussowitsch default: 379d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 380649ef022SMatthew Knepley } 381649ef022SMatthew Knepley 3829566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 1, f0_q, NULL)); 383649ef022SMatthew Knepley 3849566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu)); 3859566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL)); 3869566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL)); 3879566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL)); 3889566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 2, 2, NULL, g1_wT, NULL, g3_wT)); 389649ef022SMatthew Knepley /* Setup constants */ 390649ef022SMatthew Knepley { 391649ef022SMatthew Knepley Parameter *param; 392649ef022SMatthew Knepley PetscScalar constants[3]; 393649ef022SMatthew Knepley 3949566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 395649ef022SMatthew Knepley 396649ef022SMatthew Knepley constants[0] = param->nu; 397649ef022SMatthew Knepley constants[1] = param->alpha; 398649ef022SMatthew Knepley constants[2] = param->theta; 3999566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 3, constants)); 400649ef022SMatthew Knepley } 401649ef022SMatthew Knepley /* Setup Boundary Conditions */ 4029566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)&ctx)); 403649ef022SMatthew Knepley id = 3; 4049566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], NULL, ctx, NULL)); 405649ef022SMatthew Knepley id = 1; 4069566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], NULL, ctx, NULL)); 407649ef022SMatthew Knepley id = 2; 4089566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], NULL, ctx, NULL)); 409649ef022SMatthew Knepley id = 4; 4109566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], NULL, ctx, NULL)); 411649ef022SMatthew Knepley id = 3; 4129566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void))exactFuncs[2], NULL, ctx, NULL)); 413649ef022SMatthew Knepley id = 1; 4149566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void))exactFuncs[2], NULL, ctx, NULL)); 415649ef022SMatthew Knepley id = 2; 4169566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void))exactFuncs[2], NULL, ctx, NULL)); 417649ef022SMatthew Knepley id = 4; 4189566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void))exactFuncs[2], NULL, ctx, NULL)); 419649ef022SMatthew Knepley 420649ef022SMatthew Knepley /*setup exact solution.*/ 4219566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx)); 4229566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx)); 4239566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx)); 424*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 425649ef022SMatthew Knepley } 426649ef022SMatthew Knepley 427d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 428d71ae5a4SJacob Faibussowitsch { 429649ef022SMatthew Knepley DM cdm = dm; 430649ef022SMatthew Knepley PetscFE fe[3]; 431649ef022SMatthew Knepley Parameter *param; 432649ef022SMatthew Knepley MPI_Comm comm; 43330602db0SMatthew G. Knepley PetscInt dim; 43430602db0SMatthew G. Knepley PetscBool simplex; 435649ef022SMatthew Knepley 436649ef022SMatthew Knepley PetscFunctionBeginUser; 4379566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4389566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 439649ef022SMatthew Knepley /* Create finite element */ 4409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 4419566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0])); 4429566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity")); 443649ef022SMatthew Knepley 4449566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1])); 4459566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 4469566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure")); 447649ef022SMatthew Knepley 4489566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2])); 4499566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[2])); 4509566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[2], "temperature")); 451649ef022SMatthew Knepley 452649ef022SMatthew Knepley /* Set discretization and boundary conditions for each mesh */ 4539566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0])); 4549566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1])); 4559566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 2, NULL, (PetscObject)fe[2])); 4569566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 4579566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user)); 4589566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 459649ef022SMatthew Knepley while (cdm) { 4609566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 4619566063dSJacob Faibussowitsch PetscCall(DMPlexCreateBasisRotation(cdm, param->theta, 0.0, 0.0)); 4629566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 463649ef022SMatthew Knepley } 4649566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0])); 4659566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1])); 4669566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[2])); 467*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 468649ef022SMatthew Knepley } 469649ef022SMatthew Knepley 470d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) 471d71ae5a4SJacob Faibussowitsch { 472649ef022SMatthew Knepley Vec vec; 473649ef022SMatthew Knepley PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 474649ef022SMatthew Knepley 475649ef022SMatthew Knepley PetscFunctionBeginUser; 47663a3b9bcSJacob Faibussowitsch PetscCheck(ofield == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %" PetscInt_FMT, ofield); 477649ef022SMatthew Knepley funcs[nfield] = constant; 4789566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &vec)); 4799566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec)); 4809566063dSJacob Faibussowitsch PetscCall(VecNormalize(vec, NULL)); 4819566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)vec, "Pressure Null Space")); 4829566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(vec, NULL, "-pressure_nullspace_view")); 4839566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullSpace)); 4849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec)); 485*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 486649ef022SMatthew Knepley } 487649ef022SMatthew Knepley 488d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 489d71ae5a4SJacob Faibussowitsch { 490649ef022SMatthew Knepley SNES snes; /* nonlinear solver */ 491649ef022SMatthew Knepley DM dm; /* problem definition */ 492649ef022SMatthew Knepley Vec u, r; /* solution, residual vectors */ 493649ef022SMatthew Knepley AppCtx user; /* user-defined work context */ 494649ef022SMatthew Knepley 495327415f7SBarry Smith PetscFunctionBeginUser; 4969566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 4979566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 4989566063dSJacob Faibussowitsch PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag)); 4999566063dSJacob Faibussowitsch PetscCall(SetupParameters(&user)); 5009566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 5019566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 5029566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 5039566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &user)); 504649ef022SMatthew Knepley /* Setup problem */ 5059566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &user)); 5069566063dSJacob Faibussowitsch PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 507649ef022SMatthew Knepley 5089566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 5099566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "Solution")); 5109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 511649ef022SMatthew Knepley 5129566063dSJacob Faibussowitsch PetscCall(DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace)); 5139566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 514649ef022SMatthew Knepley 5159566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 516649ef022SMatthew Knepley { 517649ef022SMatthew Knepley PetscDS ds; 518649ef022SMatthew Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 519649ef022SMatthew Knepley void *ctxs[3]; 520649ef022SMatthew Knepley 5219566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 5229566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0])); 5239566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1])); 5249566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2])); 5259566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, u)); 5269566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution")); 5279566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-exact_vec_view")); 528649ef022SMatthew Knepley } 5299566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 5309566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 5319566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 532649ef022SMatthew Knepley 533649ef022SMatthew Knepley if (user.showError) { 534649ef022SMatthew Knepley PetscDS ds; 535649ef022SMatthew Knepley Vec r; 536649ef022SMatthew Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 537649ef022SMatthew Knepley void *ctxs[3]; 538649ef022SMatthew Knepley 5399566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 5409566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0])); 5419566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1])); 5429566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2])); 5439566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &r)); 5449566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, r)); 5459566063dSJacob Faibussowitsch PetscCall(VecAXPY(r, -1.0, u)); 5469566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r, "Solution Error")); 5479566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(r, NULL, "-error_vec_view")); 5489566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &r)); 549649ef022SMatthew Knepley } 5509566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "Numerical Solution")); 5519566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 552649ef022SMatthew Knepley 5539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 5549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 5559566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5569566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 5579566063dSJacob Faibussowitsch PetscCall(PetscBagDestroy(&user.bag)); 5589566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 559b122ec5aSJacob Faibussowitsch return 0; 560649ef022SMatthew Knepley } 561649ef022SMatthew Knepley 562649ef022SMatthew Knepley /*TEST 563649ef022SMatthew Knepley 564649ef022SMatthew Knepley test: 565649ef022SMatthew Knepley suffix: 2d_tri_p2_p1_p1 566649ef022SMatthew Knepley requires: triangle !single 567649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \ 568649ef022SMatthew Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 569649ef022SMatthew Knepley -dmsnes_check .001 -snes_error_if_not_converged \ 570649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 571649ef022SMatthew 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 \ 572649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 573649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 574649ef022SMatthew Knepley 575649ef022SMatthew Knepley test: 576649ef022SMatthew Knepley # Using -dm_refine 2 -convest_num_refine 3 gives L_2 convergence rate: [2.9, 2.3, 1.9] 577649ef022SMatthew Knepley suffix: 2d_tri_p2_p1_p1_conv 578649ef022SMatthew Knepley requires: triangle !single 579649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 580649ef022SMatthew Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 581649ef022SMatthew Knepley -snes_error_if_not_converged -snes_convergence_test correct_pressure -snes_convergence_estimate -convest_num_refine 1 \ 582649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 583649ef022SMatthew 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 \ 584649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 585649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 586649ef022SMatthew Knepley 587649ef022SMatthew Knepley test: 588649ef022SMatthew Knepley suffix: 2d_tri_p3_p2_p2 589649ef022SMatthew Knepley requires: triangle !single 590649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 591649ef022SMatthew Knepley -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \ 592649ef022SMatthew Knepley -dmsnes_check .001 -snes_error_if_not_converged \ 593649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 594649ef022SMatthew 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 \ 595649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 596649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 597649ef022SMatthew Knepley 598649ef022SMatthew Knepley TEST*/ 599