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 27649ef022SMatthew Knepley typedef enum {SOL_QUADRATIC, SOL_CUBIC, NUM_SOL_TYPES} SolType; 28649ef022SMatthew Knepley const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic", "unknown"}; 29649ef022SMatthew Knepley 30649ef022SMatthew Knepley typedef struct { 31649ef022SMatthew Knepley PetscReal nu; /* Kinematic viscosity */ 32649ef022SMatthew Knepley PetscReal theta; /* Angle of pipe wall to x-axis */ 33649ef022SMatthew Knepley PetscReal alpha; /* Thermal diffusivity */ 34649ef022SMatthew Knepley PetscReal T_in; /* Inlet temperature*/ 35649ef022SMatthew Knepley } Parameter; 36649ef022SMatthew Knepley 37649ef022SMatthew Knepley typedef struct { 3830602db0SMatthew G. Knepley PetscBool showError; 3930602db0SMatthew G. Knepley PetscBag bag; 40649ef022SMatthew Knepley SolType solType; 41649ef022SMatthew Knepley } AppCtx; 42649ef022SMatthew Knepley 43649ef022SMatthew Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 44649ef022SMatthew Knepley { 45649ef022SMatthew Knepley PetscInt d; 46649ef022SMatthew Knepley for (d = 0; d < Nc; ++d) u[d] = 0.0; 47649ef022SMatthew Knepley return 0; 48649ef022SMatthew Knepley } 49649ef022SMatthew Knepley 50649ef022SMatthew Knepley static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 51649ef022SMatthew Knepley { 52649ef022SMatthew Knepley PetscInt d; 53649ef022SMatthew Knepley for (d = 0; d < Nc; ++d) u[d] = 1.0; 54649ef022SMatthew Knepley return 0; 55649ef022SMatthew Knepley } 56649ef022SMatthew Knepley 57649ef022SMatthew Knepley /* 58649ef022SMatthew Knepley CASE: quadratic 59649ef022SMatthew Knepley In 2D we use exact solution: 60649ef022SMatthew Knepley 61649ef022SMatthew Knepley u = x^2 + y^2 62649ef022SMatthew Knepley v = 2x^2 - 2xy 63649ef022SMatthew Knepley p = x + y - 1 64649ef022SMatthew Knepley T = x + y 65649ef022SMatthew Knepley f = <2x^3 + 4x^2y - 2xy^2 -4\nu + 1, 4xy^2 + 2x^2y - 2y^3 -4\nu + 1> 66649ef022SMatthew Knepley Q = 3x^2 + y^2 - 2xy 67649ef022SMatthew Knepley 68649ef022SMatthew Knepley so that 69649ef022SMatthew Knepley 70649ef022SMatthew Knepley (1) \nabla \cdot u = 2x - 2x = 0 71649ef022SMatthew Knepley 72649ef022SMatthew Knepley (2) u \cdot \nabla u - \nu \Delta u + \nabla p - f 73649ef022SMatthew 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 74649ef022SMatthew Knepley 75649ef022SMatthew Knepley (3) u \cdot \nabla T - \alpha \Delta T - Q = 3x^2 + y^2 - 2xy - \alpha*0 - 3x^2 - y^2 + 2xy = 0 76649ef022SMatthew Knepley */ 77649ef022SMatthew Knepley 78649ef022SMatthew Knepley static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 79649ef022SMatthew Knepley { 80649ef022SMatthew Knepley u[0] = X[0]*X[0] + X[1]*X[1]; 81649ef022SMatthew Knepley u[1] = 2.0*X[0]*X[0] - 2.0*X[0]*X[1]; 82649ef022SMatthew Knepley return 0; 83649ef022SMatthew Knepley } 84649ef022SMatthew Knepley 85649ef022SMatthew Knepley static PetscErrorCode linear_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 86649ef022SMatthew Knepley { 87649ef022SMatthew Knepley p[0] = X[0] + X[1] - 1.0; 88649ef022SMatthew Knepley return 0; 89649ef022SMatthew Knepley } 90649ef022SMatthew Knepley 91649ef022SMatthew Knepley static PetscErrorCode linear_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 92649ef022SMatthew Knepley { 93649ef022SMatthew Knepley T[0] = X[0] + X[1]; 94649ef022SMatthew Knepley return 0; 95649ef022SMatthew Knepley } 96649ef022SMatthew Knepley 97649ef022SMatthew Knepley static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 98649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 99649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 100649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 101649ef022SMatthew Knepley { 102649ef022SMatthew Knepley PetscInt c, d; 103649ef022SMatthew Knepley PetscInt Nc = dim; 104649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 105649ef022SMatthew Knepley 106649ef022SMatthew Knepley for (c=0; c<Nc; ++c) { 107649ef022SMatthew Knepley for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 108649ef022SMatthew Knepley } 109649ef022SMatthew 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); 110649ef022SMatthew 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); 111649ef022SMatthew Knepley } 112649ef022SMatthew Knepley 113649ef022SMatthew Knepley static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 114649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 115649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 116649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 117649ef022SMatthew Knepley { 118649ef022SMatthew Knepley PetscInt d; 119649ef022SMatthew Knepley for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 120649ef022SMatthew Knepley f0[0] -= (3*X[0]*X[0] + X[1]*X[1] - 2*X[0]*X[1]); 121649ef022SMatthew Knepley } 122649ef022SMatthew Knepley 123649ef022SMatthew Knepley /* 124649ef022SMatthew Knepley CASE: cubic 125649ef022SMatthew Knepley In 2D we use exact solution: 126649ef022SMatthew Knepley 127649ef022SMatthew Knepley u = x^3 + y^3 128649ef022SMatthew Knepley v = 2x^3 - 3x^2y 129649ef022SMatthew Knepley p = 3/2 x^2 + 3/2 y^2 - 1 130649ef022SMatthew Knepley T = 1/2 x^2 + 1/2 y^2 131649ef022SMatthew Knepley f = <3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y), 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y> 132649ef022SMatthew Knepley Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2 133649ef022SMatthew Knepley 134649ef022SMatthew Knepley so that 135649ef022SMatthew Knepley 136649ef022SMatthew Knepley \nabla \cdot u = 3x^2 - 3x^2 = 0 137649ef022SMatthew Knepley 138649ef022SMatthew Knepley u \cdot \nabla u - \nu \Delta u + \nabla p - f 139649ef022SMatthew 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 140649ef022SMatthew Knepley 141649ef022SMatthew 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 142649ef022SMatthew Knepley */ 143649ef022SMatthew Knepley 144649ef022SMatthew Knepley static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 145649ef022SMatthew Knepley { 146649ef022SMatthew Knepley u[0] = X[0]*X[0]*X[0] + X[1]*X[1]*X[1]; 147649ef022SMatthew Knepley u[1] = 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1]; 148649ef022SMatthew Knepley return 0; 149649ef022SMatthew Knepley } 150649ef022SMatthew Knepley 151649ef022SMatthew Knepley static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 152649ef022SMatthew Knepley { 153649ef022SMatthew Knepley p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0; 154649ef022SMatthew Knepley return 0; 155649ef022SMatthew Knepley } 156649ef022SMatthew Knepley 157649ef022SMatthew Knepley static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 158649ef022SMatthew Knepley { 159649ef022SMatthew Knepley T[0] = X[0]*X[0]/2.0 + X[1]*X[1]/2.0; 160649ef022SMatthew Knepley return 0; 161649ef022SMatthew Knepley } 162649ef022SMatthew Knepley 163649ef022SMatthew Knepley static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 164649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 165649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 166649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 167649ef022SMatthew Knepley { 168649ef022SMatthew Knepley PetscInt c, d; 169649ef022SMatthew Knepley PetscInt Nc = dim; 170649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 171649ef022SMatthew Knepley 172649ef022SMatthew Knepley for (c=0; c<Nc; ++c) { 173649ef022SMatthew Knepley for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 174649ef022SMatthew Knepley } 175649ef022SMatthew 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]); 176649ef022SMatthew 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]); 177649ef022SMatthew Knepley } 178649ef022SMatthew Knepley 179649ef022SMatthew Knepley static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 180649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 181649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 182649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 183649ef022SMatthew Knepley { 184649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 185649ef022SMatthew Knepley PetscInt d; 186649ef022SMatthew Knepley 187649ef022SMatthew Knepley for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 188649ef022SMatthew 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); 189649ef022SMatthew Knepley } 190649ef022SMatthew Knepley 191649ef022SMatthew Knepley static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 192649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 193649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 194649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 195649ef022SMatthew Knepley { 196649ef022SMatthew Knepley PetscInt d; 197649ef022SMatthew Knepley for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; 198649ef022SMatthew Knepley } 199649ef022SMatthew Knepley 200649ef022SMatthew Knepley static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 201649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 202649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 203649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 204649ef022SMatthew Knepley { 205649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 206649ef022SMatthew Knepley const PetscInt Nc = dim; 207649ef022SMatthew Knepley PetscInt c, d; 208649ef022SMatthew Knepley 209649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 210649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 211649ef022SMatthew Knepley f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]); 212649ef022SMatthew Knepley //f1[c*dim+d] = nu*u_x[c*dim+d]; 213649ef022SMatthew Knepley } 214649ef022SMatthew Knepley f1[c*dim+c] -= u[uOff[1]]; 215649ef022SMatthew Knepley } 216649ef022SMatthew Knepley } 217649ef022SMatthew Knepley 218649ef022SMatthew Knepley static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 219649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 220649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 221649ef022SMatthew Knepley PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 222649ef022SMatthew Knepley { 223649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 224649ef022SMatthew Knepley PetscInt d; 225649ef022SMatthew Knepley for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d]; 226649ef022SMatthew Knepley } 227649ef022SMatthew Knepley 228649ef022SMatthew Knepley /*Jacobians*/ 229649ef022SMatthew Knepley static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 230649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 231649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 232649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 233649ef022SMatthew Knepley { 234649ef022SMatthew Knepley PetscInt d; 235649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; 236649ef022SMatthew Knepley } 237649ef022SMatthew Knepley 238649ef022SMatthew Knepley static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 239649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 240649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 241649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 242649ef022SMatthew Knepley { 243649ef022SMatthew Knepley const PetscInt Nc = dim; 244649ef022SMatthew Knepley PetscInt c, d; 245649ef022SMatthew Knepley 246649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 247649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 248649ef022SMatthew Knepley g0[c*Nc+d] = u_x[ c*Nc+d]; 249649ef022SMatthew Knepley } 250649ef022SMatthew Knepley } 251649ef022SMatthew Knepley } 252649ef022SMatthew Knepley 253649ef022SMatthew Knepley static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 254649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 255649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 256649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 257649ef022SMatthew Knepley { 258649ef022SMatthew Knepley PetscInt NcI = dim; 259649ef022SMatthew Knepley PetscInt NcJ = dim; 260649ef022SMatthew Knepley PetscInt c, d, e; 261649ef022SMatthew Knepley 262649ef022SMatthew Knepley for (c = 0; c < NcI; ++c) { 263649ef022SMatthew Knepley for (d = 0; d < NcJ; ++d) { 264649ef022SMatthew Knepley for (e = 0; e < dim; ++e) { 265649ef022SMatthew Knepley if (c == d) { 266649ef022SMatthew Knepley g1[(c*NcJ+d)*dim+e] = u[e]; 267649ef022SMatthew Knepley } 268649ef022SMatthew Knepley } 269649ef022SMatthew Knepley } 270649ef022SMatthew Knepley } 271649ef022SMatthew Knepley } 272649ef022SMatthew Knepley 273649ef022SMatthew Knepley static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 274649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 275649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 276649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 277649ef022SMatthew Knepley { 278649ef022SMatthew Knepley PetscInt d; 279649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; 280649ef022SMatthew Knepley } 281649ef022SMatthew Knepley 282649ef022SMatthew Knepley static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 283649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 284649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 285649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 286649ef022SMatthew Knepley { 287649ef022SMatthew Knepley const PetscReal nu = PetscRealPart(constants[0]); 288649ef022SMatthew Knepley const PetscInt Nc = dim; 289649ef022SMatthew Knepley PetscInt c, d; 290649ef022SMatthew Knepley 291649ef022SMatthew Knepley for (c = 0; c < Nc; ++c) { 292649ef022SMatthew Knepley for (d = 0; d < dim; ++d) { 293649ef022SMatthew Knepley g3[((c*Nc+c)*dim+d)*dim+d] += nu; // gradU 294649ef022SMatthew Knepley g3[((c*Nc+d)*dim+d)*dim+c] += nu; // gradU transpose 295649ef022SMatthew Knepley } 296649ef022SMatthew Knepley } 297649ef022SMatthew Knepley } 298649ef022SMatthew Knepley 299649ef022SMatthew Knepley static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 300649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 301649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 302649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 303649ef022SMatthew Knepley { 304649ef022SMatthew Knepley PetscInt d; 305649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d]; 306649ef022SMatthew Knepley } 307649ef022SMatthew Knepley 308649ef022SMatthew Knepley static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 309649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 310649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 311649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 312649ef022SMatthew Knepley { 313649ef022SMatthew Knepley PetscInt d; 314649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d]; 315649ef022SMatthew Knepley } 316649ef022SMatthew Knepley 317649ef022SMatthew Knepley static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 318649ef022SMatthew Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 319649ef022SMatthew Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 320649ef022SMatthew Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 321649ef022SMatthew Knepley { 322649ef022SMatthew Knepley const PetscReal alpha = PetscRealPart(constants[1]); 323649ef022SMatthew Knepley PetscInt d; 324649ef022SMatthew Knepley 325649ef022SMatthew Knepley for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha; 326649ef022SMatthew Knepley } 327649ef022SMatthew Knepley 328649ef022SMatthew Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 329649ef022SMatthew Knepley { 33030602db0SMatthew G. Knepley PetscInt sol; 331649ef022SMatthew Knepley PetscErrorCode ierr; 332649ef022SMatthew Knepley 333649ef022SMatthew Knepley PetscFunctionBeginUser; 334649ef022SMatthew Knepley options->solType = SOL_QUADRATIC; 335649ef022SMatthew Knepley options->showError = PETSC_FALSE; 336649ef022SMatthew Knepley 337649ef022SMatthew Knepley ierr = PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");CHKERRQ(ierr); 338649ef022SMatthew Knepley sol = options->solType; 339*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL)); 340649ef022SMatthew Knepley options->solType = (SolType) sol; 341*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-show_error", "Output the error for verification", "ex62.c", options->showError, &options->showError, NULL)); 3421e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 343649ef022SMatthew Knepley PetscFunctionReturn(0); 344649ef022SMatthew Knepley } 345649ef022SMatthew Knepley 346649ef022SMatthew Knepley static PetscErrorCode SetupParameters(AppCtx *user) 347649ef022SMatthew Knepley { 348649ef022SMatthew Knepley PetscBag bag; 349649ef022SMatthew Knepley Parameter *p; 350649ef022SMatthew Knepley 351649ef022SMatthew Knepley PetscFunctionBeginUser; 352649ef022SMatthew Knepley /* setup PETSc parameter bag */ 353*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(user->bag, (void **) &p)); 354*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagSetName(user->bag, "par", "Poiseuille flow parameters")); 355649ef022SMatthew Knepley bag = user->bag; 356*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity")); 357*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity")); 358*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagRegisterReal(bag, &p->theta, 0.0, "theta", "Angle of pipe wall to x-axis")); 359649ef022SMatthew Knepley PetscFunctionReturn(0); 360649ef022SMatthew Knepley } 361649ef022SMatthew Knepley 362649ef022SMatthew Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 363649ef022SMatthew Knepley { 364649ef022SMatthew Knepley PetscFunctionBeginUser; 365*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 366*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 367*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 368649ef022SMatthew Knepley { 369649ef022SMatthew Knepley Parameter *param; 370649ef022SMatthew Knepley Vec coordinates; 371649ef022SMatthew Knepley PetscScalar *coords; 372649ef022SMatthew Knepley PetscReal theta; 373649ef022SMatthew Knepley PetscInt cdim, N, bs, i; 374649ef022SMatthew Knepley 375*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDim(*dm, &cdim)); 376*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(*dm, &coordinates)); 377*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(coordinates, &N)); 378*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetBlockSize(coordinates, &bs)); 3792c71b3e2SJacob Faibussowitsch PetscCheckFalse(bs != cdim,comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %D != embedding dimension %D", bs, cdim); 380*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(coordinates, &coords)); 381*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(user->bag, (void **) ¶m)); 382649ef022SMatthew Knepley theta = param->theta; 383649ef022SMatthew Knepley for (i = 0; i < N; i += cdim) { 384649ef022SMatthew Knepley PetscScalar x = coords[i+0]; 385649ef022SMatthew Knepley PetscScalar y = coords[i+1]; 386649ef022SMatthew Knepley 387649ef022SMatthew Knepley coords[i+0] = PetscCosReal(theta)*x - PetscSinReal(theta)*y; 388649ef022SMatthew Knepley coords[i+1] = PetscSinReal(theta)*x + PetscCosReal(theta)*y; 389649ef022SMatthew Knepley } 390*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(coordinates, &coords)); 391*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetCoordinates(*dm, coordinates)); 392649ef022SMatthew Knepley } 393*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 394649ef022SMatthew Knepley PetscFunctionReturn(0); 395649ef022SMatthew Knepley } 396649ef022SMatthew Knepley 397649ef022SMatthew Knepley static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 398649ef022SMatthew Knepley { 399649ef022SMatthew Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 400649ef022SMatthew Knepley PetscDS prob; 40145480ffeSMatthew G. Knepley DMLabel label; 402649ef022SMatthew Knepley Parameter *ctx; 403649ef022SMatthew Knepley PetscInt id; 404649ef022SMatthew Knepley 405649ef022SMatthew Knepley PetscFunctionBeginUser; 406*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 407*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &prob)); 408649ef022SMatthew Knepley switch(user->solType) { 409649ef022SMatthew Knepley case SOL_QUADRATIC: 410*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v)); 411*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w)); 412649ef022SMatthew Knepley 413649ef022SMatthew Knepley exactFuncs[0] = quadratic_u; 414649ef022SMatthew Knepley exactFuncs[1] = linear_p; 415649ef022SMatthew Knepley exactFuncs[2] = linear_T; 416649ef022SMatthew Knepley break; 417649ef022SMatthew Knepley case SOL_CUBIC: 418*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v)); 419*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w)); 420649ef022SMatthew Knepley 421649ef022SMatthew Knepley exactFuncs[0] = cubic_u; 422649ef022SMatthew Knepley exactFuncs[1] = quadratic_p; 423649ef022SMatthew Knepley exactFuncs[2] = quadratic_T; 424649ef022SMatthew Knepley break; 42598921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 426649ef022SMatthew Knepley } 427649ef022SMatthew Knepley 428*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 1, f0_q, NULL)); 429649ef022SMatthew Knepley 430*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu)); 431*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL)); 432*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL)); 433*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL)); 434*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 2, 2, NULL, g1_wT, NULL, g3_wT)); 435649ef022SMatthew Knepley /* Setup constants */ 436649ef022SMatthew Knepley { 437649ef022SMatthew Knepley Parameter *param; 438649ef022SMatthew Knepley PetscScalar constants[3]; 439649ef022SMatthew Knepley 440*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(user->bag, (void **) ¶m)); 441649ef022SMatthew Knepley 442649ef022SMatthew Knepley constants[0] = param->nu; 443649ef022SMatthew Knepley constants[1] = param->alpha; 444649ef022SMatthew Knepley constants[2] = param->theta; 445*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetConstants(prob, 3, constants)); 446649ef022SMatthew Knepley } 447649ef022SMatthew Knepley /* Setup Boundary Conditions */ 448*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(user->bag, (void **) &ctx)); 449649ef022SMatthew Knepley id = 3; 450*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL)); 451649ef022SMatthew Knepley id = 1; 452*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL)); 453649ef022SMatthew Knepley id = 2; 454*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL)); 455649ef022SMatthew Knepley id = 4; 456*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL)); 457649ef022SMatthew Knepley id = 3; 458*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL)); 459649ef022SMatthew Knepley id = 1; 460*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL)); 461649ef022SMatthew Knepley id = 2; 462*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL)); 463649ef022SMatthew Knepley id = 4; 464*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL)); 465649ef022SMatthew Knepley 466649ef022SMatthew Knepley /*setup exact solution.*/ 467*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx)); 468*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx)); 469*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx)); 470649ef022SMatthew Knepley PetscFunctionReturn(0); 471649ef022SMatthew Knepley } 472649ef022SMatthew Knepley 473649ef022SMatthew Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 474649ef022SMatthew Knepley { 475649ef022SMatthew Knepley DM cdm = dm; 476649ef022SMatthew Knepley PetscFE fe[3]; 477649ef022SMatthew Knepley Parameter *param; 478649ef022SMatthew Knepley MPI_Comm comm; 47930602db0SMatthew G. Knepley PetscInt dim; 48030602db0SMatthew G. Knepley PetscBool simplex; 481649ef022SMatthew Knepley 482649ef022SMatthew Knepley PetscFunctionBeginUser; 483*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 484*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 485649ef022SMatthew Knepley /* Create finite element */ 486*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 487*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0])); 488*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe[0], "velocity")); 489649ef022SMatthew Knepley 490*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1])); 491*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECopyQuadrature(fe[0], fe[1])); 492*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe[1], "pressure")); 493649ef022SMatthew Knepley 494*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2])); 495*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECopyQuadrature(fe[0], fe[2])); 496*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe[2], "temperature")); 497649ef022SMatthew Knepley 498649ef022SMatthew Knepley /* Set discretization and boundary conditions for each mesh */ 499*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe[0])); 500*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 1, NULL, (PetscObject) fe[1])); 501*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 2, NULL, (PetscObject) fe[2])); 502*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 503*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupProblem(dm, user)); 504*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(user->bag, (void **) ¶m)); 505649ef022SMatthew Knepley while (cdm) { 506*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCopyDisc(dm, cdm)); 507*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateBasisRotation(cdm, param->theta, 0.0, 0.0)); 508*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 509649ef022SMatthew Knepley } 510*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe[0])); 511*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe[1])); 512*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe[2])); 513649ef022SMatthew Knepley PetscFunctionReturn(0); 514649ef022SMatthew Knepley } 515649ef022SMatthew Knepley 516649ef022SMatthew Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) 517649ef022SMatthew Knepley { 518649ef022SMatthew Knepley Vec vec; 519649ef022SMatthew Knepley PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 520649ef022SMatthew Knepley 521649ef022SMatthew Knepley PetscFunctionBeginUser; 5222c71b3e2SJacob Faibussowitsch PetscCheckFalse(ofield != 1,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield); 523649ef022SMatthew Knepley funcs[nfield] = constant; 524*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &vec)); 525*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec)); 526*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNormalize(vec, NULL)); 527*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) vec, "Pressure Null Space")); 528*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(vec, NULL, "-pressure_nullspace_view")); 529*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace)); 530*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&vec)); 531649ef022SMatthew Knepley PetscFunctionReturn(0); 532649ef022SMatthew Knepley } 533649ef022SMatthew Knepley 534649ef022SMatthew Knepley int main(int argc, char **argv) 535649ef022SMatthew Knepley { 536649ef022SMatthew Knepley SNES snes; /* nonlinear solver */ 537649ef022SMatthew Knepley DM dm; /* problem definition */ 538649ef022SMatthew Knepley Vec u, r; /* solution, residual vectors */ 539649ef022SMatthew Knepley AppCtx user; /* user-defined work context */ 540649ef022SMatthew Knepley PetscErrorCode ierr; 541649ef022SMatthew Knepley 542649ef022SMatthew Knepley ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 543*5f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 544*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag)); 545*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupParameters(&user)); 546*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD, &snes)); 547*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 548*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetDM(snes, dm)); 549*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(dm, &user)); 550649ef022SMatthew Knepley /* Setup problem */ 551*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(dm, &user)); 552*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateClosureIndex(dm, NULL)); 553649ef022SMatthew Knepley 554*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &u)); 555*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) u, "Solution")); 556*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u, &r)); 557649ef022SMatthew Knepley 558*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace)); 559*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetSNESLocalFEM(dm,&user,&user,&user)); 560649ef022SMatthew Knepley 561*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 562649ef022SMatthew Knepley { 563649ef022SMatthew Knepley PetscDS ds; 564649ef022SMatthew Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 565649ef022SMatthew Knepley void *ctxs[3]; 566649ef022SMatthew Knepley 567*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 568*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0])); 569*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1])); 570*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2])); 571*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, u)); 572*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) u, "Exact Solution")); 573*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-exact_vec_view")); 574649ef022SMatthew Knepley } 575*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSNESCheckFromOptions(snes, u)); 576*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(u, 0.0)); 577*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes, NULL, u)); 578649ef022SMatthew Knepley 579649ef022SMatthew Knepley if (user.showError) { 580649ef022SMatthew Knepley PetscDS ds; 581649ef022SMatthew Knepley Vec r; 582649ef022SMatthew Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 583649ef022SMatthew Knepley void *ctxs[3]; 584649ef022SMatthew Knepley 585*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 586*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0])); 587*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1])); 588*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2])); 589*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dm, &r)); 590*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, r)); 591*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(r, -1.0, u)); 592*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) r, "Solution Error")); 593*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(r, NULL, "-error_vec_view")); 594*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dm, &r)); 595649ef022SMatthew Knepley } 596*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) u, "Numerical Solution")); 597*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-sol_vec_view")); 598649ef022SMatthew Knepley 599*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 600*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 601*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 602*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 603*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagDestroy(&user.bag)); 604649ef022SMatthew Knepley ierr = PetscFinalize(); 605649ef022SMatthew Knepley return ierr; 606649ef022SMatthew Knepley } 607649ef022SMatthew Knepley 608649ef022SMatthew Knepley /*TEST 609649ef022SMatthew Knepley 610649ef022SMatthew Knepley test: 611649ef022SMatthew Knepley suffix: 2d_tri_p2_p1_p1 612649ef022SMatthew Knepley requires: triangle !single 613649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \ 614649ef022SMatthew Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 615649ef022SMatthew Knepley -dmsnes_check .001 -snes_error_if_not_converged \ 616649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 617649ef022SMatthew 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 \ 618649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 619649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 620649ef022SMatthew Knepley 621649ef022SMatthew Knepley test: 622649ef022SMatthew Knepley # Using -dm_refine 2 -convest_num_refine 3 gives L_2 convergence rate: [2.9, 2.3, 1.9] 623649ef022SMatthew Knepley suffix: 2d_tri_p2_p1_p1_conv 624649ef022SMatthew Knepley requires: triangle !single 625649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 626649ef022SMatthew Knepley -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 627649ef022SMatthew Knepley -snes_error_if_not_converged -snes_convergence_test correct_pressure -snes_convergence_estimate -convest_num_refine 1 \ 628649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 629649ef022SMatthew 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 \ 630649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 631649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 632649ef022SMatthew Knepley 633649ef022SMatthew Knepley test: 634649ef022SMatthew Knepley suffix: 2d_tri_p3_p2_p2 635649ef022SMatthew Knepley requires: triangle !single 636649ef022SMatthew Knepley args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 637649ef022SMatthew Knepley -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \ 638649ef022SMatthew Knepley -dmsnes_check .001 -snes_error_if_not_converged \ 639649ef022SMatthew Knepley -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 640649ef022SMatthew 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 \ 641649ef022SMatthew Knepley -fieldsplit_0_pc_type lu \ 642649ef022SMatthew Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 643649ef022SMatthew Knepley 644649ef022SMatthew Knepley TEST*/ 645