1c4762a1bSJed Brown static char help[] = "Heat Equation in 2d and 3d with finite elements.\n\ 2c4762a1bSJed Brown We solve the heat equation in a rectangular\n\ 3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4c4762a1bSJed Brown Contributed by: Julian Andrej <juan@tf.uni-kiel.de>\n\n\n"; 5c4762a1bSJed Brown 6c4762a1bSJed Brown #include <petscdmplex.h> 7c4762a1bSJed Brown #include <petscds.h> 8c4762a1bSJed Brown #include <petscts.h> 9c4762a1bSJed Brown 10c4762a1bSJed Brown /* 11c4762a1bSJed Brown Heat equation: 12c4762a1bSJed Brown 13a3d0cf85SMatthew G. Knepley du/dt - \Delta u + f = 0 14c4762a1bSJed Brown */ 15c4762a1bSJed Brown 16742ee2edSMatthew G. Knepley typedef enum {SOL_QUADRATIC_LINEAR, SOL_QUADRATIC_TRIG, SOL_TRIG_LINEAR, SOL_TRIG_TRIG, NUM_SOLUTION_TYPES} SolutionType; 17742ee2edSMatthew G. Knepley const char *solutionTypes[NUM_SOLUTION_TYPES+1] = {"quadratic_linear", "quadratic_trig", "trig_linear", "trig_trig", "unknown"}; 18a3d0cf85SMatthew G. Knepley 19c4762a1bSJed Brown typedef struct { 20a3d0cf85SMatthew G. Knepley SolutionType solType; /* Type of exact solution */ 21742ee2edSMatthew G. Knepley /* Solver setup */ 22742ee2edSMatthew G. Knepley PetscBool expTS; /* Flag for explicit timestepping */ 23742ee2edSMatthew G. Knepley PetscBool lumped; /* Lump the mass matrix */ 24c4762a1bSJed Brown } AppCtx; 25c4762a1bSJed Brown 26a3d0cf85SMatthew G. Knepley /* 27a3d0cf85SMatthew G. Knepley Exact 2D solution: 28a3d0cf85SMatthew G. Knepley u = 2t + x^2 + y^2 29742ee2edSMatthew G. Knepley u_t = 2 30742ee2edSMatthew G. Knepley \Delta u = 2 + 2 = 4 31742ee2edSMatthew G. Knepley f = 2 32a3d0cf85SMatthew G. Knepley F(u) = 2 - (2 + 2) + 2 = 0 33a3d0cf85SMatthew G. Knepley 34a3d0cf85SMatthew G. Knepley Exact 3D solution: 35a3d0cf85SMatthew G. Knepley u = 3t + x^2 + y^2 + z^2 36a3d0cf85SMatthew G. Knepley F(u) = 3 - (2 + 2 + 2) + 3 = 0 37a3d0cf85SMatthew G. Knepley */ 38a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_quad_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 39c4762a1bSJed Brown { 40c4762a1bSJed Brown PetscInt d; 41c4762a1bSJed Brown 42c4762a1bSJed Brown *u = dim*time; 43c4762a1bSJed Brown for (d = 0; d < dim; ++d) *u += x[d]*x[d]; 44c4762a1bSJed Brown return 0; 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_quad_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 48a3d0cf85SMatthew G. Knepley { 49a3d0cf85SMatthew G. Knepley *u = dim; 50a3d0cf85SMatthew G. Knepley return 0; 51a3d0cf85SMatthew G. Knepley } 52a3d0cf85SMatthew G. Knepley 53742ee2edSMatthew G. Knepley static void f0_quad_lin_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 54742ee2edSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 55742ee2edSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 56742ee2edSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 57742ee2edSMatthew G. Knepley { 58742ee2edSMatthew G. Knepley f0[0] = -(PetscScalar) dim; 59742ee2edSMatthew G. Knepley } 60a3d0cf85SMatthew G. Knepley static void f0_quad_lin(PetscInt dim, PetscInt Nf, PetscInt NfAux, 61c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 62c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 63c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 64c4762a1bSJed Brown { 65742ee2edSMatthew G. Knepley PetscScalar exp[1] = {0.}; 66742ee2edSMatthew G. Knepley f0_quad_lin_exp(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, exp); 67742ee2edSMatthew G. Knepley f0[0] = u_t[0] - exp[0]; 68c4762a1bSJed Brown } 69c4762a1bSJed Brown 70a3d0cf85SMatthew G. Knepley /* 71a3d0cf85SMatthew G. Knepley Exact 2D solution: 72a3d0cf85SMatthew G. Knepley u = 2*cos(t) + x^2 + y^2 73a3d0cf85SMatthew G. Knepley F(u) = -2*sint(t) - (2 + 2) + 2*sin(t) + 4 = 0 74a3d0cf85SMatthew G. Knepley 75a3d0cf85SMatthew G. Knepley Exact 3D solution: 76a3d0cf85SMatthew G. Knepley u = 3*cos(t) + x^2 + y^2 + z^2 77a3d0cf85SMatthew G. Knepley F(u) = -3*sin(t) - (2 + 2 + 2) + 3*sin(t) + 6 = 0 78a3d0cf85SMatthew G. Knepley */ 79a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_quad_trig(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 80a3d0cf85SMatthew G. Knepley { 81a3d0cf85SMatthew G. Knepley PetscInt d; 82a3d0cf85SMatthew G. Knepley 83a3d0cf85SMatthew G. Knepley *u = dim*PetscCosReal(time); 84a3d0cf85SMatthew G. Knepley for (d = 0; d < dim; ++d) *u += x[d]*x[d]; 85a3d0cf85SMatthew G. Knepley return 0; 86a3d0cf85SMatthew G. Knepley } 87a3d0cf85SMatthew G. Knepley 88a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_quad_trig_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 89a3d0cf85SMatthew G. Knepley { 90a3d0cf85SMatthew G. Knepley *u = -dim*PetscSinReal(time); 91a3d0cf85SMatthew G. Knepley return 0; 92a3d0cf85SMatthew G. Knepley } 93a3d0cf85SMatthew G. Knepley 94742ee2edSMatthew G. Knepley static void f0_quad_trig_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 95742ee2edSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 96742ee2edSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 97742ee2edSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 98742ee2edSMatthew G. Knepley { 99742ee2edSMatthew G. Knepley f0[0] = -dim*(PetscSinReal(t) + 2.0); 100742ee2edSMatthew G. Knepley } 101a3d0cf85SMatthew G. Knepley static void f0_quad_trig(PetscInt dim, PetscInt Nf, PetscInt NfAux, 102a3d0cf85SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 103a3d0cf85SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 104a3d0cf85SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 105a3d0cf85SMatthew G. Knepley { 106742ee2edSMatthew G. Knepley PetscScalar exp[1] = {0.}; 107742ee2edSMatthew G. Knepley f0_quad_trig_exp(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, exp); 108742ee2edSMatthew G. Knepley f0[0] = u_t[0] - exp[0]; 109a3d0cf85SMatthew G. Knepley } 110a3d0cf85SMatthew G. Knepley 111a3d0cf85SMatthew G. Knepley /* 112a3d0cf85SMatthew G. Knepley Exact 2D solution: 113a3d0cf85SMatthew G. Knepley u = 2\pi^2 t + cos(\pi x) + cos(\pi y) 114a3d0cf85SMatthew G. Knepley F(u) = 2\pi^2 - \pi^2 (cos(\pi x) + cos(\pi y)) + \pi^2 (cos(\pi x) + cos(\pi y)) - 2\pi^2 = 0 115a3d0cf85SMatthew G. Knepley 116a3d0cf85SMatthew G. Knepley Exact 3D solution: 117a3d0cf85SMatthew G. Knepley u = 3\pi^2 t + cos(\pi x) + cos(\pi y) + cos(\pi z) 118a3d0cf85SMatthew G. Knepley F(u) = 3\pi^2 - \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z)) + \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z)) - 3\pi^2 = 0 119a3d0cf85SMatthew G. Knepley */ 120a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_trig_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 121a3d0cf85SMatthew G. Knepley { 122a3d0cf85SMatthew G. Knepley PetscInt d; 123a3d0cf85SMatthew G. Knepley 124a3d0cf85SMatthew G. Knepley *u = dim*PetscSqr(PETSC_PI)*time; 125a3d0cf85SMatthew G. Knepley for (d = 0; d < dim; ++d) *u += PetscCosReal(PETSC_PI*x[d]); 126a3d0cf85SMatthew G. Knepley return 0; 127a3d0cf85SMatthew G. Knepley } 128a3d0cf85SMatthew G. Knepley 129a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_trig_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 130a3d0cf85SMatthew G. Knepley { 131a3d0cf85SMatthew G. Knepley *u = dim*PetscSqr(PETSC_PI); 132a3d0cf85SMatthew G. Knepley return 0; 133a3d0cf85SMatthew G. Knepley } 134a3d0cf85SMatthew G. Knepley 135a3d0cf85SMatthew G. Knepley static void f0_trig_lin(PetscInt dim, PetscInt Nf, PetscInt NfAux, 136a3d0cf85SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 137a3d0cf85SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 138a3d0cf85SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 139a3d0cf85SMatthew G. Knepley { 140a3d0cf85SMatthew G. Knepley PetscInt d; 141a3d0cf85SMatthew G. Knepley f0[0] = u_t[0]; 142a3d0cf85SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += PetscSqr(PETSC_PI)*(PetscCosReal(PETSC_PI*x[d]) - 1.0); 143a3d0cf85SMatthew G. Knepley } 144a3d0cf85SMatthew G. Knepley 145742ee2edSMatthew G. Knepley /* 146742ee2edSMatthew G. Knepley Exact 2D solution: 147742ee2edSMatthew G. Knepley u = pi^2 cos(t) + cos(\pi x) + cos(\pi y) 148742ee2edSMatthew G. Knepley u_t = -pi^2 sin(t) 149742ee2edSMatthew G. Knepley \Delta u = -\pi^2 (cos(\pi x) + cos(\pi y)) 150742ee2edSMatthew G. Knepley f = pi^2 sin(t) - \pi^2 (cos(\pi x) + cos(\pi y)) 151742ee2edSMatthew G. Knepley F(u) = -\pi^2 sin(t) + \pi^2 (cos(\pi x) + cos(\pi y)) - \pi^2 (cos(\pi x) + cos(\pi y)) + \pi^2 sin(t) = 0 152742ee2edSMatthew G. Knepley 153742ee2edSMatthew G. Knepley Exact 3D solution: 154742ee2edSMatthew G. Knepley u = pi^2 cos(t) + cos(\pi x) + cos(\pi y) + cos(\pi z) 155742ee2edSMatthew G. Knepley u_t = -pi^2 sin(t) 156742ee2edSMatthew G. Knepley \Delta u = -\pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z)) 157742ee2edSMatthew G. Knepley f = pi^2 sin(t) - \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z)) 158742ee2edSMatthew G. Knepley F(u) = -\pi^2 sin(t) + \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z)) - \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z)) + \pi^2 sin(t) = 0 159742ee2edSMatthew G. Knepley */ 160742ee2edSMatthew G. Knepley static PetscErrorCode mms_trig_trig(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 161742ee2edSMatthew G. Knepley { 162742ee2edSMatthew G. Knepley PetscInt d; 163742ee2edSMatthew G. Knepley 164742ee2edSMatthew G. Knepley *u = PetscSqr(PETSC_PI)*PetscCosReal(time); 165742ee2edSMatthew G. Knepley for (d = 0; d < dim; ++d) *u += PetscCosReal(PETSC_PI*x[d]); 166742ee2edSMatthew G. Knepley return 0; 167742ee2edSMatthew G. Knepley } 168742ee2edSMatthew G. Knepley 169742ee2edSMatthew G. Knepley static PetscErrorCode mms_trig_trig_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 170742ee2edSMatthew G. Knepley { 171742ee2edSMatthew G. Knepley *u = -PetscSqr(PETSC_PI)*PetscSinReal(time); 172742ee2edSMatthew G. Knepley return 0; 173742ee2edSMatthew G. Knepley } 174742ee2edSMatthew G. Knepley 175742ee2edSMatthew G. Knepley static void f0_trig_trig_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 176742ee2edSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 177742ee2edSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 178742ee2edSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 179742ee2edSMatthew G. Knepley { 180742ee2edSMatthew G. Knepley PetscInt d; 181742ee2edSMatthew G. Knepley f0[0] -= PetscSqr(PETSC_PI)*PetscSinReal(t); 182742ee2edSMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += PetscSqr(PETSC_PI)*PetscCosReal(PETSC_PI*x[d]); 183742ee2edSMatthew G. Knepley } 184742ee2edSMatthew G. Knepley static void f0_trig_trig(PetscInt dim, PetscInt Nf, PetscInt NfAux, 185742ee2edSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 186742ee2edSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 187742ee2edSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 188742ee2edSMatthew G. Knepley { 189742ee2edSMatthew G. Knepley PetscScalar exp[1] = {0.}; 190742ee2edSMatthew G. Knepley f0_trig_trig_exp(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, exp); 191742ee2edSMatthew G. Knepley f0[0] = u_t[0] - exp[0]; 192742ee2edSMatthew G. Knepley } 193742ee2edSMatthew G. Knepley 194742ee2edSMatthew G. Knepley static void f1_temp_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 195742ee2edSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 196742ee2edSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 197742ee2edSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 198742ee2edSMatthew G. Knepley { 199742ee2edSMatthew G. Knepley PetscInt d; 200742ee2edSMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = -u_x[d]; 201742ee2edSMatthew G. Knepley } 202c4762a1bSJed Brown static void f1_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 203c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 204c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 205c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 206c4762a1bSJed Brown { 207c4762a1bSJed Brown PetscInt d; 208a3d0cf85SMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 209c4762a1bSJed Brown } 210c4762a1bSJed Brown 211c4762a1bSJed Brown static void g3_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 212c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 213c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 214c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 215c4762a1bSJed Brown { 216c4762a1bSJed Brown PetscInt d; 217a3d0cf85SMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 220c4762a1bSJed Brown static void g0_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 221c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 222c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 223c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 224c4762a1bSJed Brown { 225c4762a1bSJed Brown g0[0] = u_tShift*1.0; 226c4762a1bSJed Brown } 227c4762a1bSJed Brown 228c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 229c4762a1bSJed Brown { 230a3d0cf85SMatthew G. Knepley PetscInt sol; 231c4762a1bSJed Brown PetscErrorCode ierr; 232c4762a1bSJed Brown 233c4762a1bSJed Brown PetscFunctionBeginUser; 234a3d0cf85SMatthew G. Knepley options->solType = SOL_QUADRATIC_LINEAR; 235742ee2edSMatthew G. Knepley options->expTS = PETSC_FALSE; 236742ee2edSMatthew G. Knepley options->lumped = PETSC_TRUE; 237c4762a1bSJed Brown 238c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Heat Equation Options", "DMPLEX");CHKERRQ(ierr); 239a3d0cf85SMatthew G. Knepley ierr = PetscOptionsEList("-sol_type", "Type of exact solution", "ex45.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL);CHKERRQ(ierr); 240a3d0cf85SMatthew G. Knepley options->solType = (SolutionType) sol; 241742ee2edSMatthew G. Knepley ierr = PetscOptionsBool("-explicit", "Use explicit timestepping", "ex45.c", options->expTS, &options->expTS, NULL);CHKERRQ(ierr); 242742ee2edSMatthew G. Knepley ierr = PetscOptionsBool("-lumped", "Lump the mass matrix", "ex45.c", options->lumped, &options->lumped, NULL);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 244c4762a1bSJed Brown PetscFunctionReturn(0); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown 247c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx) 248c4762a1bSJed Brown { 249c4762a1bSJed Brown PetscErrorCode ierr; 250c4762a1bSJed Brown 251c4762a1bSJed Brown PetscFunctionBeginUser; 25230602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 25330602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 254c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 255c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 256c4762a1bSJed Brown PetscFunctionReturn(0); 257c4762a1bSJed Brown } 258c4762a1bSJed Brown 259c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) 260c4762a1bSJed Brown { 261a3d0cf85SMatthew G. Knepley PetscDS ds; 26245480ffeSMatthew G. Knepley DMLabel label; 263c4762a1bSJed Brown const PetscInt id = 1; 264c4762a1bSJed Brown PetscErrorCode ierr; 265c4762a1bSJed Brown 266c4762a1bSJed Brown PetscFunctionBeginUser; 26745480ffeSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 268a3d0cf85SMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 269a3d0cf85SMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 0, 0, g0_temp, NULL, NULL, g3_temp);CHKERRQ(ierr); 270a3d0cf85SMatthew G. Knepley switch (ctx->solType) { 271a3d0cf85SMatthew G. Knepley case SOL_QUADRATIC_LINEAR: 272a3d0cf85SMatthew G. Knepley ierr = PetscDSSetResidual(ds, 0, f0_quad_lin, f1_temp);CHKERRQ(ierr); 273742ee2edSMatthew G. Knepley ierr = PetscDSSetRHSResidual(ds, 0, f0_quad_lin_exp, f1_temp_exp);CHKERRQ(ierr); 274a3d0cf85SMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, mms_quad_lin, ctx);CHKERRQ(ierr); 275a3d0cf85SMatthew G. Knepley ierr = PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_lin_t, ctx);CHKERRQ(ierr); 27645480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms_quad_lin, (void (*)(void)) mms_quad_lin_t, ctx, NULL);CHKERRQ(ierr); 277a3d0cf85SMatthew G. Knepley break; 278a3d0cf85SMatthew G. Knepley case SOL_QUADRATIC_TRIG: 279a3d0cf85SMatthew G. Knepley ierr = PetscDSSetResidual(ds, 0, f0_quad_trig, f1_temp);CHKERRQ(ierr); 280742ee2edSMatthew G. Knepley ierr = PetscDSSetRHSResidual(ds, 0, f0_quad_trig_exp, f1_temp_exp);CHKERRQ(ierr); 281a3d0cf85SMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, mms_quad_trig, ctx);CHKERRQ(ierr); 282a3d0cf85SMatthew G. Knepley ierr = PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_trig_t, ctx);CHKERRQ(ierr); 28345480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms_quad_trig, (void (*)(void)) mms_quad_trig_t, ctx, NULL);CHKERRQ(ierr); 284a3d0cf85SMatthew G. Knepley break; 285a3d0cf85SMatthew G. Knepley case SOL_TRIG_LINEAR: 286a3d0cf85SMatthew G. Knepley ierr = PetscDSSetResidual(ds, 0, f0_trig_lin, f1_temp);CHKERRQ(ierr); 287a3d0cf85SMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, mms_trig_lin, ctx);CHKERRQ(ierr); 288a3d0cf85SMatthew G. Knepley ierr = PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_lin_t, ctx);CHKERRQ(ierr); 28945480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms_trig_lin, (void (*)(void)) mms_trig_lin_t, ctx, NULL);CHKERRQ(ierr); 290a3d0cf85SMatthew G. Knepley break; 291742ee2edSMatthew G. Knepley case SOL_TRIG_TRIG: 292742ee2edSMatthew G. Knepley ierr = PetscDSSetResidual(ds, 0, f0_trig_trig, f1_temp);CHKERRQ(ierr); 293742ee2edSMatthew G. Knepley ierr = PetscDSSetRHSResidual(ds, 0, f0_trig_trig_exp, f1_temp_exp);CHKERRQ(ierr); 294742ee2edSMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, mms_trig_trig, ctx);CHKERRQ(ierr); 295742ee2edSMatthew G. Knepley ierr = PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_trig_t, ctx);CHKERRQ(ierr); 296742ee2edSMatthew G. Knepley break; 297*98921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%D)", solutionTypes[PetscMin(ctx->solType, NUM_SOLUTION_TYPES)], ctx->solType); 298a3d0cf85SMatthew G. Knepley } 299c4762a1bSJed Brown PetscFunctionReturn(0); 300c4762a1bSJed Brown } 301c4762a1bSJed Brown 302c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx) 303c4762a1bSJed Brown { 304c4762a1bSJed Brown DM cdm = dm; 305c4762a1bSJed Brown PetscFE fe; 306a3d0cf85SMatthew G. Knepley DMPolytopeType ct; 307a3d0cf85SMatthew G. Knepley PetscBool simplex; 308a3d0cf85SMatthew G. Knepley PetscInt dim, cStart; 309c4762a1bSJed Brown PetscErrorCode ierr; 310c4762a1bSJed Brown 311c4762a1bSJed Brown PetscFunctionBeginUser; 312a3d0cf85SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 313a3d0cf85SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 314a3d0cf85SMatthew G. Knepley ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 315a3d0cf85SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 316c4762a1bSJed Brown /* Create finite element */ 317a3d0cf85SMatthew G. Knepley ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "temp_", -1, &fe);CHKERRQ(ierr); 318c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe, "temperature");CHKERRQ(ierr); 319c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 320c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 321c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 322742ee2edSMatthew G. Knepley if (ctx->expTS) { 323742ee2edSMatthew G. Knepley PetscDS ds; 324742ee2edSMatthew G. Knepley 325742ee2edSMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 326742ee2edSMatthew G. Knepley ierr = PetscDSSetImplicit(ds, 0, PETSC_FALSE);CHKERRQ(ierr); 327742ee2edSMatthew G. Knepley } 328c4762a1bSJed Brown ierr = SetupProblem(dm, ctx);CHKERRQ(ierr); 329c4762a1bSJed Brown while (cdm) { 330408cafa0SMatthew G. Knepley ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 331c4762a1bSJed Brown ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 332c4762a1bSJed Brown } 333c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 334c4762a1bSJed Brown PetscFunctionReturn(0); 335c4762a1bSJed Brown } 336c4762a1bSJed Brown 337a3d0cf85SMatthew G. Knepley static PetscErrorCode SetInitialConditions(TS ts, Vec u) 338a3d0cf85SMatthew G. Knepley { 339a3d0cf85SMatthew G. Knepley DM dm; 340a3d0cf85SMatthew G. Knepley PetscReal t; 341a3d0cf85SMatthew G. Knepley PetscErrorCode ierr; 342a3d0cf85SMatthew G. Knepley 343a3d0cf85SMatthew G. Knepley PetscFunctionBegin; 344a3d0cf85SMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 345a3d0cf85SMatthew G. Knepley ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 346a3d0cf85SMatthew G. Knepley ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr); 347a3d0cf85SMatthew G. Knepley PetscFunctionReturn(0); 348a3d0cf85SMatthew G. Knepley } 349a3d0cf85SMatthew G. Knepley 350c4762a1bSJed Brown int main(int argc, char **argv) 351c4762a1bSJed Brown { 352c4762a1bSJed Brown DM dm; 353c4762a1bSJed Brown TS ts; 354a3d0cf85SMatthew G. Knepley Vec u; 355a3d0cf85SMatthew G. Knepley AppCtx ctx; 356c4762a1bSJed Brown PetscErrorCode ierr; 357c4762a1bSJed Brown 358c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 359c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr); 360c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &dm, &ctx);CHKERRQ(ierr); 361c4762a1bSJed Brown ierr = DMSetApplicationContext(dm, &ctx);CHKERRQ(ierr); 362c4762a1bSJed Brown ierr = SetupDiscretization(dm, &ctx);CHKERRQ(ierr); 363c4762a1bSJed Brown 364c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 365c4762a1bSJed Brown ierr = TSSetDM(ts, dm);CHKERRQ(ierr); 366c4762a1bSJed Brown ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx);CHKERRQ(ierr); 367742ee2edSMatthew G. Knepley if (ctx.expTS) { 368d64986aaSMatthew G. Knepley ierr = DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFEM, &ctx);CHKERRQ(ierr); 369742ee2edSMatthew G. Knepley if (ctx.lumped) {ierr = DMTSCreateRHSMassMatrixLumped(dm);CHKERRQ(ierr);} 370742ee2edSMatthew G. Knepley else {ierr = DMTSCreateRHSMassMatrix(dm);CHKERRQ(ierr);} 371742ee2edSMatthew G. Knepley } else { 372c4762a1bSJed Brown ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx);CHKERRQ(ierr); 373c4762a1bSJed Brown ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx);CHKERRQ(ierr); 374742ee2edSMatthew G. Knepley } 375a3d0cf85SMatthew G. Knepley ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 376c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 377a3d0cf85SMatthew G. Knepley ierr = TSSetComputeInitialCondition(ts, SetInitialConditions);CHKERRQ(ierr); 378c4762a1bSJed Brown 379a3d0cf85SMatthew G. Knepley ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 380a3d0cf85SMatthew G. Knepley ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 381a3d0cf85SMatthew G. Knepley ierr = SetInitialConditions(ts, u);CHKERRQ(ierr); 382a3d0cf85SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) u, "temperature");CHKERRQ(ierr); 383c4762a1bSJed Brown ierr = TSSolve(ts, u);CHKERRQ(ierr); 384a3d0cf85SMatthew G. Knepley ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 385742ee2edSMatthew G. Knepley if (ctx.expTS) {ierr = DMTSDestroyRHSMassMatrix(dm);CHKERRQ(ierr);} 386c4762a1bSJed Brown 387c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 388c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 389c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 390c4762a1bSJed Brown ierr = PetscFinalize(); 391c4762a1bSJed Brown return ierr; 392c4762a1bSJed Brown } 393c4762a1bSJed Brown 394c4762a1bSJed Brown /*TEST 395c4762a1bSJed Brown 396c4762a1bSJed Brown test: 397a3d0cf85SMatthew G. Knepley suffix: 2d_p1 398c4762a1bSJed Brown requires: triangle 399a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 400a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 401c4762a1bSJed Brown test: 402742ee2edSMatthew G. Knepley suffix: 2d_p1_exp 403742ee2edSMatthew G. Knepley requires: triangle 404742ee2edSMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -explicit \ 405742ee2edSMatthew G. Knepley -ts_type euler -ts_max_steps 4 -ts_dt 1e-3 -ts_monitor_error 406742ee2edSMatthew G. Knepley test: 407a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 408a3d0cf85SMatthew G. Knepley suffix: 2d_p1_sconv 409c4762a1bSJed Brown requires: triangle 410a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 411a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 412c4762a1bSJed Brown test: 413742ee2edSMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.1] 414742ee2edSMatthew G. Knepley suffix: 2d_p1_sconv_2 415742ee2edSMatthew G. Knepley requires: triangle 416742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 417742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 1e-6 -snes_error_if_not_converged -pc_type lu 418742ee2edSMatthew G. Knepley test: 419a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 420a3d0cf85SMatthew G. Knepley suffix: 2d_p1_tconv 421c4762a1bSJed Brown requires: triangle 422a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 423a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 424c4762a1bSJed Brown test: 425742ee2edSMatthew G. Knepley # -dm_refine 6 -convest_num_refine 3 get L_2 convergence rate: [1.0] 426742ee2edSMatthew G. Knepley suffix: 2d_p1_tconv_2 427742ee2edSMatthew G. Knepley requires: triangle 428742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 429742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 430742ee2edSMatthew G. Knepley test: 431742ee2edSMatthew G. Knepley # The L_2 convergence rate cannot be seen since stability of the explicit integrator requires that is be more accurate than the grid 432742ee2edSMatthew G. Knepley suffix: 2d_p1_exp_tconv_2 433742ee2edSMatthew G. Knepley requires: triangle 434742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \ 435742ee2edSMatthew G. Knepley -ts_type euler -ts_max_steps 4 -ts_dt 1e-4 -lumped 0 -mass_pc_type lu 436742ee2edSMatthew G. Knepley test: 437742ee2edSMatthew G. Knepley # The L_2 convergence rate cannot be seen since stability of the explicit integrator requires that is be more accurate than the grid 438742ee2edSMatthew G. Knepley suffix: 2d_p1_exp_tconv_2_lump 439742ee2edSMatthew G. Knepley requires: triangle 440742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \ 441742ee2edSMatthew G. Knepley -ts_type euler -ts_max_steps 4 -ts_dt 1e-4 442742ee2edSMatthew G. Knepley test: 443a3d0cf85SMatthew G. Knepley suffix: 2d_p2 444c4762a1bSJed Brown requires: triangle 445a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 446a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 447c4762a1bSJed Brown test: 448a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 449a3d0cf85SMatthew G. Knepley suffix: 2d_p2_sconv 450a3d0cf85SMatthew G. Knepley requires: triangle 451a3d0cf85SMatthew G. Knepley args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 452a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 453c4762a1bSJed Brown test: 454742ee2edSMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [3.1] 455742ee2edSMatthew G. Knepley suffix: 2d_p2_sconv_2 456742ee2edSMatthew G. Knepley requires: triangle 457742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 458742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 459742ee2edSMatthew G. Knepley test: 460a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 461a3d0cf85SMatthew G. Knepley suffix: 2d_p2_tconv 462a3d0cf85SMatthew G. Knepley requires: triangle 463a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 464a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 465c4762a1bSJed Brown test: 466742ee2edSMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 467742ee2edSMatthew G. Knepley suffix: 2d_p2_tconv_2 468742ee2edSMatthew G. Knepley requires: triangle 469742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 470742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 471742ee2edSMatthew G. Knepley test: 472a3d0cf85SMatthew G. Knepley suffix: 2d_q1 47330602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 474a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 475c4762a1bSJed Brown test: 476a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 477a3d0cf85SMatthew G. Knepley suffix: 2d_q1_sconv 47830602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 479a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 480c4762a1bSJed Brown test: 481a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 482a3d0cf85SMatthew G. Knepley suffix: 2d_q1_tconv 48330602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 484a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 485a3d0cf85SMatthew G. Knepley test: 486a3d0cf85SMatthew G. Knepley suffix: 2d_q2 48730602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 488a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 489a3d0cf85SMatthew G. Knepley test: 490a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 491a3d0cf85SMatthew G. Knepley suffix: 2d_q2_sconv 49230602db0SMatthew G. Knepley args: -sol_type trig_linear -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 493a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 494a3d0cf85SMatthew G. Knepley test: 495a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 496a3d0cf85SMatthew G. Knepley suffix: 2d_q2_tconv 49730602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 498a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 499a3d0cf85SMatthew G. Knepley 500a3d0cf85SMatthew G. Knepley test: 501a3d0cf85SMatthew G. Knepley suffix: 3d_p1 502c4762a1bSJed Brown requires: ctetgen 503a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 504a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 505c4762a1bSJed Brown test: 506a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 507a3d0cf85SMatthew G. Knepley suffix: 3d_p1_sconv 508c4762a1bSJed Brown requires: ctetgen 509a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 510a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 511c4762a1bSJed Brown test: 512a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 513a3d0cf85SMatthew G. Knepley suffix: 3d_p1_tconv 514c4762a1bSJed Brown requires: ctetgen 515a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 516a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 517c4762a1bSJed Brown test: 518a3d0cf85SMatthew G. Knepley suffix: 3d_p2 519c4762a1bSJed Brown requires: ctetgen 520a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 521a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 522c4762a1bSJed Brown test: 523a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 524a3d0cf85SMatthew G. Knepley suffix: 3d_p2_sconv 525a3d0cf85SMatthew G. Knepley requires: ctetgen 526a3d0cf85SMatthew G. Knepley args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 527a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 528c4762a1bSJed Brown test: 529a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 530a3d0cf85SMatthew G. Knepley suffix: 3d_p2_tconv 531a3d0cf85SMatthew G. Knepley requires: ctetgen 532a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 533a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 534c4762a1bSJed Brown test: 535a3d0cf85SMatthew G. Knepley suffix: 3d_q1 53630602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 537a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 538c4762a1bSJed Brown test: 539a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 540a3d0cf85SMatthew G. Knepley suffix: 3d_q1_sconv 54130602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 542a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 543a3d0cf85SMatthew G. Knepley test: 544a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 545a3d0cf85SMatthew G. Knepley suffix: 3d_q1_tconv 54630602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 547a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 548a3d0cf85SMatthew G. Knepley test: 549a3d0cf85SMatthew G. Knepley suffix: 3d_q2 55030602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 551a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 552a3d0cf85SMatthew G. Knepley test: 553a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 554a3d0cf85SMatthew G. Knepley suffix: 3d_q2_sconv 55530602db0SMatthew G. Knepley args: -sol_type trig_linear -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 556a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 557a3d0cf85SMatthew G. Knepley test: 558a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 559a3d0cf85SMatthew G. Knepley suffix: 3d_q2_tconv 56030602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 561a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 562a3d0cf85SMatthew G. Knepley 563a3d0cf85SMatthew G. Knepley test: 564a3d0cf85SMatthew G. Knepley # For a nice picture, -bd_dm_refine 2 -dm_refine 1 -dm_view hdf5:${PETSC_DIR}/sol.h5 -ts_monitor_solution hdf5:${PETSC_DIR}/sol.h5::append 565a3d0cf85SMatthew G. Knepley suffix: egads_sphere 566a3d0cf85SMatthew G. Knepley requires: egads ctetgen 56730602db0SMatthew G. Knepley args: -sol_type quadratic_linear \ 56830602db0SMatthew G. Knepley -dm_plex_boundary_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/unit_sphere.egadslite -dm_plex_boundary_label marker -bd_dm_plex_scale 40 \ 569a3d0cf85SMatthew G. Knepley -temp_petscspace_degree 2 -dmts_check .0001 \ 570a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 571c4762a1bSJed Brown 572c4762a1bSJed Brown TEST*/ 573