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); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-sol_type", "Type of exact solution", "ex45.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL)); 240a3d0cf85SMatthew G. Knepley options->solType = (SolutionType) sol; 2415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-explicit", "Use explicit timestepping", "ex45.c", options->expTS, &options->expTS, NULL)); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-lumped", "Lump the mass matrix", "ex45.c", options->lumped, &options->lumped, NULL)); 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 PetscFunctionBeginUser; 2505f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 2515f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 2525f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 254c4762a1bSJed Brown PetscFunctionReturn(0); 255c4762a1bSJed Brown } 256c4762a1bSJed Brown 257c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) 258c4762a1bSJed Brown { 259a3d0cf85SMatthew G. Knepley PetscDS ds; 26045480ffeSMatthew G. Knepley DMLabel label; 261c4762a1bSJed Brown const PetscInt id = 1; 262c4762a1bSJed Brown 263c4762a1bSJed Brown PetscFunctionBeginUser; 2645f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 2655f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, g0_temp, NULL, NULL, g3_temp)); 267a3d0cf85SMatthew G. Knepley switch (ctx->solType) { 268a3d0cf85SMatthew G. Knepley case SOL_QUADRATIC_LINEAR: 2695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_quad_lin, f1_temp)); 2705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetRHSResidual(ds, 0, f0_quad_lin_exp, f1_temp_exp)); 2715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(ds, 0, mms_quad_lin, ctx)); 2725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_lin_t, ctx)); 2735f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 274a3d0cf85SMatthew G. Knepley break; 275a3d0cf85SMatthew G. Knepley case SOL_QUADRATIC_TRIG: 2765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_quad_trig, f1_temp)); 2775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetRHSResidual(ds, 0, f0_quad_trig_exp, f1_temp_exp)); 2785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(ds, 0, mms_quad_trig, ctx)); 2795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_trig_t, ctx)); 2805f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 281a3d0cf85SMatthew G. Knepley break; 282a3d0cf85SMatthew G. Knepley case SOL_TRIG_LINEAR: 2835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_trig_lin, f1_temp)); 2845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(ds, 0, mms_trig_lin, ctx)); 2855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_lin_t, ctx)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 287a3d0cf85SMatthew G. Knepley break; 288742ee2edSMatthew G. Knepley case SOL_TRIG_TRIG: 2895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_trig_trig, f1_temp)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetRHSResidual(ds, 0, f0_trig_trig_exp, f1_temp_exp)); 2915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(ds, 0, mms_trig_trig, ctx)); 2925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_trig_t, ctx)); 293742ee2edSMatthew G. Knepley break; 29498921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%D)", solutionTypes[PetscMin(ctx->solType, NUM_SOLUTION_TYPES)], ctx->solType); 295a3d0cf85SMatthew G. Knepley } 296c4762a1bSJed Brown PetscFunctionReturn(0); 297c4762a1bSJed Brown } 298c4762a1bSJed Brown 299c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx) 300c4762a1bSJed Brown { 301c4762a1bSJed Brown DM cdm = dm; 302c4762a1bSJed Brown PetscFE fe; 303a3d0cf85SMatthew G. Knepley DMPolytopeType ct; 304a3d0cf85SMatthew G. Knepley PetscBool simplex; 305a3d0cf85SMatthew G. Knepley PetscInt dim, cStart; 306c4762a1bSJed Brown 307c4762a1bSJed Brown PetscFunctionBeginUser; 3085f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 3095f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 3105f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dm, cStart, &ct)); 311a3d0cf85SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 312c4762a1bSJed Brown /* Create finite element */ 3135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "temp_", -1, &fe)); 3145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe, "temperature")); 315c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 3165f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe)); 3175f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 318742ee2edSMatthew G. Knepley if (ctx->expTS) { 319742ee2edSMatthew G. Knepley PetscDS ds; 320742ee2edSMatthew G. Knepley 3215f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 3225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetImplicit(ds, 0, PETSC_FALSE)); 323742ee2edSMatthew G. Knepley } 3245f80ce2aSJacob Faibussowitsch CHKERRQ(SetupProblem(dm, ctx)); 325c4762a1bSJed Brown while (cdm) { 3265f80ce2aSJacob Faibussowitsch CHKERRQ(DMCopyDisc(dm, cdm)); 3275f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 328c4762a1bSJed Brown } 3295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 330c4762a1bSJed Brown PetscFunctionReturn(0); 331c4762a1bSJed Brown } 332c4762a1bSJed Brown 333a3d0cf85SMatthew G. Knepley static PetscErrorCode SetInitialConditions(TS ts, Vec u) 334a3d0cf85SMatthew G. Knepley { 335a3d0cf85SMatthew G. Knepley DM dm; 336a3d0cf85SMatthew G. Knepley PetscReal t; 337a3d0cf85SMatthew G. Knepley 338a3d0cf85SMatthew G. Knepley PetscFunctionBegin; 3395f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dm)); 3405f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTime(ts, &t)); 3415f80ce2aSJacob Faibussowitsch CHKERRQ(DMComputeExactSolution(dm, t, u, NULL)); 342a3d0cf85SMatthew G. Knepley PetscFunctionReturn(0); 343a3d0cf85SMatthew G. Knepley } 344a3d0cf85SMatthew G. Knepley 345c4762a1bSJed Brown int main(int argc, char **argv) 346c4762a1bSJed Brown { 347c4762a1bSJed Brown DM dm; 348c4762a1bSJed Brown TS ts; 349a3d0cf85SMatthew G. Knepley Vec u; 350a3d0cf85SMatthew G. Knepley AppCtx ctx; 351c4762a1bSJed Brown 352*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 3535f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 3545f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx)); 3555f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(dm, &ctx)); 3565f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(dm, &ctx)); 357c4762a1bSJed Brown 3585f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD, &ts)); 3595f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts, dm)); 3605f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx)); 361742ee2edSMatthew G. Knepley if (ctx.expTS) { 3625f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFEM, &ctx)); 3635f80ce2aSJacob Faibussowitsch if (ctx.lumped) CHKERRQ(DMTSCreateRHSMassMatrixLumped(dm)); 3645f80ce2aSJacob Faibussowitsch else CHKERRQ(DMTSCreateRHSMassMatrix(dm)); 365742ee2edSMatthew G. Knepley } else { 3665f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx)); 3675f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx)); 368742ee2edSMatthew G. Knepley } 3695f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 3705f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 3715f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetComputeInitialCondition(ts, SetInitialConditions)); 372c4762a1bSJed Brown 3735f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &u)); 3745f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSCheckFromOptions(ts, u)); 3755f80ce2aSJacob Faibussowitsch CHKERRQ(SetInitialConditions(ts, u)); 3765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) u, "temperature")); 3775f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts, u)); 3785f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSCheckFromOptions(ts, u)); 3795f80ce2aSJacob Faibussowitsch if (ctx.expTS) CHKERRQ(DMTSDestroyRHSMassMatrix(dm)); 380c4762a1bSJed Brown 3815f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 3825f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 3835f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 384*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 385*b122ec5aSJacob Faibussowitsch return 0; 386c4762a1bSJed Brown } 387c4762a1bSJed Brown 388c4762a1bSJed Brown /*TEST 389c4762a1bSJed Brown 390c4762a1bSJed Brown test: 391a3d0cf85SMatthew G. Knepley suffix: 2d_p1 392c4762a1bSJed Brown requires: triangle 393a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 394a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 395c4762a1bSJed Brown test: 396742ee2edSMatthew G. Knepley suffix: 2d_p1_exp 397742ee2edSMatthew G. Knepley requires: triangle 398742ee2edSMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -explicit \ 399742ee2edSMatthew G. Knepley -ts_type euler -ts_max_steps 4 -ts_dt 1e-3 -ts_monitor_error 400742ee2edSMatthew G. Knepley test: 401a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 402a3d0cf85SMatthew G. Knepley suffix: 2d_p1_sconv 403c4762a1bSJed Brown requires: triangle 404a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 405a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 406c4762a1bSJed Brown test: 407742ee2edSMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.1] 408742ee2edSMatthew G. Knepley suffix: 2d_p1_sconv_2 409742ee2edSMatthew G. Knepley requires: triangle 410742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 411742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 1e-6 -snes_error_if_not_converged -pc_type lu 412742ee2edSMatthew G. Knepley test: 413a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 414a3d0cf85SMatthew G. Knepley suffix: 2d_p1_tconv 415c4762a1bSJed Brown requires: triangle 416a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 417a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 418c4762a1bSJed Brown test: 419742ee2edSMatthew G. Knepley # -dm_refine 6 -convest_num_refine 3 get L_2 convergence rate: [1.0] 420742ee2edSMatthew G. Knepley suffix: 2d_p1_tconv_2 421742ee2edSMatthew G. Knepley requires: triangle 422742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 423742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 424742ee2edSMatthew G. Knepley test: 425742ee2edSMatthew 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 426742ee2edSMatthew G. Knepley suffix: 2d_p1_exp_tconv_2 427742ee2edSMatthew G. Knepley requires: triangle 428742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \ 429742ee2edSMatthew G. Knepley -ts_type euler -ts_max_steps 4 -ts_dt 1e-4 -lumped 0 -mass_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_lump 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 436742ee2edSMatthew G. Knepley test: 437a3d0cf85SMatthew G. Knepley suffix: 2d_p2 438c4762a1bSJed Brown requires: triangle 439a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 440a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 441c4762a1bSJed Brown test: 442a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 443a3d0cf85SMatthew G. Knepley suffix: 2d_p2_sconv 444a3d0cf85SMatthew G. Knepley requires: triangle 445a3d0cf85SMatthew G. Knepley args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 446a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 447c4762a1bSJed Brown test: 448742ee2edSMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [3.1] 449742ee2edSMatthew G. Knepley suffix: 2d_p2_sconv_2 450742ee2edSMatthew G. Knepley requires: triangle 451742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 452742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 453742ee2edSMatthew G. Knepley test: 454a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 455a3d0cf85SMatthew G. Knepley suffix: 2d_p2_tconv 456a3d0cf85SMatthew G. Knepley requires: triangle 457a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 458a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 459c4762a1bSJed Brown test: 460742ee2edSMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 461742ee2edSMatthew G. Knepley suffix: 2d_p2_tconv_2 462742ee2edSMatthew G. Knepley requires: triangle 463742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 464742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 465742ee2edSMatthew G. Knepley test: 466a3d0cf85SMatthew G. Knepley suffix: 2d_q1 46730602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 468a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 469c4762a1bSJed Brown test: 470a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 471a3d0cf85SMatthew G. Knepley suffix: 2d_q1_sconv 47230602db0SMatthew 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 \ 473a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 474c4762a1bSJed Brown test: 475a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 476a3d0cf85SMatthew G. Knepley suffix: 2d_q1_tconv 47730602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 478a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 479a3d0cf85SMatthew G. Knepley test: 480a3d0cf85SMatthew G. Knepley suffix: 2d_q2 48130602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 482a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 483a3d0cf85SMatthew G. Knepley test: 484a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 485a3d0cf85SMatthew G. Knepley suffix: 2d_q2_sconv 48630602db0SMatthew 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 \ 487a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 488a3d0cf85SMatthew G. Knepley test: 489a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 490a3d0cf85SMatthew G. Knepley suffix: 2d_q2_tconv 49130602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 492a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 493a3d0cf85SMatthew G. Knepley 494a3d0cf85SMatthew G. Knepley test: 495a3d0cf85SMatthew G. Knepley suffix: 3d_p1 496c4762a1bSJed Brown requires: ctetgen 497a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 498a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 499c4762a1bSJed Brown test: 500a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 501a3d0cf85SMatthew G. Knepley suffix: 3d_p1_sconv 502c4762a1bSJed Brown requires: ctetgen 503a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 504a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 505c4762a1bSJed Brown test: 506a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 507a3d0cf85SMatthew G. Knepley suffix: 3d_p1_tconv 508c4762a1bSJed Brown requires: ctetgen 509a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 510a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 511c4762a1bSJed Brown test: 512a3d0cf85SMatthew G. Knepley suffix: 3d_p2 513c4762a1bSJed Brown requires: ctetgen 514a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 515a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 516c4762a1bSJed Brown test: 517a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 518a3d0cf85SMatthew G. Knepley suffix: 3d_p2_sconv 519a3d0cf85SMatthew G. Knepley requires: ctetgen 520a3d0cf85SMatthew G. Knepley args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 521a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 522c4762a1bSJed Brown test: 523a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 524a3d0cf85SMatthew G. Knepley suffix: 3d_p2_tconv 525a3d0cf85SMatthew G. Knepley requires: ctetgen 526a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 527a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 528c4762a1bSJed Brown test: 529a3d0cf85SMatthew G. Knepley suffix: 3d_q1 53030602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 531a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 532c4762a1bSJed Brown test: 533a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 534a3d0cf85SMatthew G. Knepley suffix: 3d_q1_sconv 53530602db0SMatthew 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 \ 536a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 537a3d0cf85SMatthew G. Knepley test: 538a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 539a3d0cf85SMatthew G. Knepley suffix: 3d_q1_tconv 54030602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 541a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 542a3d0cf85SMatthew G. Knepley test: 543a3d0cf85SMatthew G. Knepley suffix: 3d_q2 54430602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 545a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 546a3d0cf85SMatthew G. Knepley test: 547a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 548a3d0cf85SMatthew G. Knepley suffix: 3d_q2_sconv 54930602db0SMatthew 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 \ 550a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 551a3d0cf85SMatthew G. Knepley test: 552a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 553a3d0cf85SMatthew G. Knepley suffix: 3d_q2_tconv 55430602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 555a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 556a3d0cf85SMatthew G. Knepley 557a3d0cf85SMatthew G. Knepley test: 558a3d0cf85SMatthew 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 559a3d0cf85SMatthew G. Knepley suffix: egads_sphere 560a3d0cf85SMatthew G. Knepley requires: egads ctetgen 56130602db0SMatthew G. Knepley args: -sol_type quadratic_linear \ 56230602db0SMatthew 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 \ 563a3d0cf85SMatthew G. Knepley -temp_petscspace_degree 2 -dmts_check .0001 \ 564a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 565c4762a1bSJed Brown 566c4762a1bSJed Brown TEST*/ 567