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 232c4762a1bSJed Brown PetscFunctionBeginUser; 233a3d0cf85SMatthew G. Knepley options->solType = SOL_QUADRATIC_LINEAR; 234742ee2edSMatthew G. Knepley options->expTS = PETSC_FALSE; 235742ee2edSMatthew G. Knepley options->lumped = PETSC_TRUE; 236c4762a1bSJed Brown 237*d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Heat Equation Options", "DMPLEX"); 2389566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-sol_type", "Type of exact solution", "ex45.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL)); 239a3d0cf85SMatthew G. Knepley options->solType = (SolutionType) sol; 2409566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-explicit", "Use explicit timestepping", "ex45.c", options->expTS, &options->expTS, NULL)); 2419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-lumped", "Lump the mass matrix", "ex45.c", options->lumped, &options->lumped, NULL)); 242*d0609cedSBarry Smith PetscOptionsEnd(); 243c4762a1bSJed Brown PetscFunctionReturn(0); 244c4762a1bSJed Brown } 245c4762a1bSJed Brown 246c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx) 247c4762a1bSJed Brown { 248c4762a1bSJed Brown PetscFunctionBeginUser; 2499566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 2509566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 2519566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2529566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 253c4762a1bSJed Brown PetscFunctionReturn(0); 254c4762a1bSJed Brown } 255c4762a1bSJed Brown 256c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) 257c4762a1bSJed Brown { 258a3d0cf85SMatthew G. Knepley PetscDS ds; 25945480ffeSMatthew G. Knepley DMLabel label; 260c4762a1bSJed Brown const PetscInt id = 1; 261c4762a1bSJed Brown 262c4762a1bSJed Brown PetscFunctionBeginUser; 2639566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 2649566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2659566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_temp, NULL, NULL, g3_temp)); 266a3d0cf85SMatthew G. Knepley switch (ctx->solType) { 267a3d0cf85SMatthew G. Knepley case SOL_QUADRATIC_LINEAR: 2689566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_quad_lin, f1_temp)); 2699566063dSJacob Faibussowitsch PetscCall(PetscDSSetRHSResidual(ds, 0, f0_quad_lin_exp, f1_temp_exp)); 2709566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms_quad_lin, ctx)); 2719566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_lin_t, ctx)); 2729566063dSJacob Faibussowitsch PetscCall(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)); 273a3d0cf85SMatthew G. Knepley break; 274a3d0cf85SMatthew G. Knepley case SOL_QUADRATIC_TRIG: 2759566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_quad_trig, f1_temp)); 2769566063dSJacob Faibussowitsch PetscCall(PetscDSSetRHSResidual(ds, 0, f0_quad_trig_exp, f1_temp_exp)); 2779566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms_quad_trig, ctx)); 2789566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_trig_t, ctx)); 2799566063dSJacob Faibussowitsch PetscCall(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)); 280a3d0cf85SMatthew G. Knepley break; 281a3d0cf85SMatthew G. Knepley case SOL_TRIG_LINEAR: 2829566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_lin, f1_temp)); 2839566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms_trig_lin, ctx)); 2849566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_lin_t, ctx)); 2859566063dSJacob Faibussowitsch PetscCall(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)); 286a3d0cf85SMatthew G. Knepley break; 287742ee2edSMatthew G. Knepley case SOL_TRIG_TRIG: 2889566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_trig, f1_temp)); 2899566063dSJacob Faibussowitsch PetscCall(PetscDSSetRHSResidual(ds, 0, f0_trig_trig_exp, f1_temp_exp)); 2909566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms_trig_trig, ctx)); 2919566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_trig_t, ctx)); 292742ee2edSMatthew G. Knepley break; 29398921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%D)", solutionTypes[PetscMin(ctx->solType, NUM_SOLUTION_TYPES)], ctx->solType); 294a3d0cf85SMatthew G. Knepley } 295c4762a1bSJed Brown PetscFunctionReturn(0); 296c4762a1bSJed Brown } 297c4762a1bSJed Brown 298c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx) 299c4762a1bSJed Brown { 300c4762a1bSJed Brown DM cdm = dm; 301c4762a1bSJed Brown PetscFE fe; 302a3d0cf85SMatthew G. Knepley DMPolytopeType ct; 303a3d0cf85SMatthew G. Knepley PetscBool simplex; 304a3d0cf85SMatthew G. Knepley PetscInt dim, cStart; 305c4762a1bSJed Brown 306c4762a1bSJed Brown PetscFunctionBeginUser; 3079566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 3089566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 3099566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 310a3d0cf85SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 311c4762a1bSJed Brown /* Create finite element */ 3129566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "temp_", -1, &fe)); 3139566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe, "temperature")); 314c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 3159566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 3169566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 317742ee2edSMatthew G. Knepley if (ctx->expTS) { 318742ee2edSMatthew G. Knepley PetscDS ds; 319742ee2edSMatthew G. Knepley 3209566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3219566063dSJacob Faibussowitsch PetscCall(PetscDSSetImplicit(ds, 0, PETSC_FALSE)); 322742ee2edSMatthew G. Knepley } 3239566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, ctx)); 324c4762a1bSJed Brown while (cdm) { 3259566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 3269566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 327c4762a1bSJed Brown } 3289566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 329c4762a1bSJed Brown PetscFunctionReturn(0); 330c4762a1bSJed Brown } 331c4762a1bSJed Brown 332a3d0cf85SMatthew G. Knepley static PetscErrorCode SetInitialConditions(TS ts, Vec u) 333a3d0cf85SMatthew G. Knepley { 334a3d0cf85SMatthew G. Knepley DM dm; 335a3d0cf85SMatthew G. Knepley PetscReal t; 336a3d0cf85SMatthew G. Knepley 337a3d0cf85SMatthew G. Knepley PetscFunctionBegin; 3389566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 3399566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 3409566063dSJacob Faibussowitsch PetscCall(DMComputeExactSolution(dm, t, u, NULL)); 341a3d0cf85SMatthew G. Knepley PetscFunctionReturn(0); 342a3d0cf85SMatthew G. Knepley } 343a3d0cf85SMatthew G. Knepley 344c4762a1bSJed Brown int main(int argc, char **argv) 345c4762a1bSJed Brown { 346c4762a1bSJed Brown DM dm; 347c4762a1bSJed Brown TS ts; 348a3d0cf85SMatthew G. Knepley Vec u; 349a3d0cf85SMatthew G. Knepley AppCtx ctx; 350c4762a1bSJed Brown 3519566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3529566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 3539566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx)); 3549566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &ctx)); 3559566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &ctx)); 356c4762a1bSJed Brown 3579566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 3589566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm)); 3599566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx)); 360742ee2edSMatthew G. Knepley if (ctx.expTS) { 3619566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFEM, &ctx)); 3629566063dSJacob Faibussowitsch if (ctx.lumped) PetscCall(DMTSCreateRHSMassMatrixLumped(dm)); 3639566063dSJacob Faibussowitsch else PetscCall(DMTSCreateRHSMassMatrix(dm)); 364742ee2edSMatthew G. Knepley } else { 3659566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx)); 3669566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx)); 367742ee2edSMatthew G. Knepley } 3689566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 3699566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 3709566063dSJacob Faibussowitsch PetscCall(TSSetComputeInitialCondition(ts, SetInitialConditions)); 371c4762a1bSJed Brown 3729566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 3739566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u)); 3749566063dSJacob Faibussowitsch PetscCall(SetInitialConditions(ts, u)); 3759566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "temperature")); 3769566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 3779566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u)); 3789566063dSJacob Faibussowitsch if (ctx.expTS) PetscCall(DMTSDestroyRHSMassMatrix(dm)); 379c4762a1bSJed Brown 3809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 3819566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 3829566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3839566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 384b122ec5aSJacob Faibussowitsch return 0; 385c4762a1bSJed Brown } 386c4762a1bSJed Brown 387c4762a1bSJed Brown /*TEST 388c4762a1bSJed Brown 389c4762a1bSJed Brown test: 390a3d0cf85SMatthew G. Knepley suffix: 2d_p1 391c4762a1bSJed Brown requires: triangle 392a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 393a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 394c4762a1bSJed Brown test: 395742ee2edSMatthew G. Knepley suffix: 2d_p1_exp 396742ee2edSMatthew G. Knepley requires: triangle 397742ee2edSMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -explicit \ 398742ee2edSMatthew G. Knepley -ts_type euler -ts_max_steps 4 -ts_dt 1e-3 -ts_monitor_error 399742ee2edSMatthew G. Knepley test: 400a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 401a3d0cf85SMatthew G. Knepley suffix: 2d_p1_sconv 402c4762a1bSJed Brown requires: triangle 403a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 404a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 405c4762a1bSJed Brown test: 406742ee2edSMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.1] 407742ee2edSMatthew G. Knepley suffix: 2d_p1_sconv_2 408742ee2edSMatthew G. Knepley requires: triangle 409742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 410742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 1e-6 -snes_error_if_not_converged -pc_type lu 411742ee2edSMatthew G. Knepley test: 412a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 413a3d0cf85SMatthew G. Knepley suffix: 2d_p1_tconv 414c4762a1bSJed Brown requires: triangle 415a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 416a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 417c4762a1bSJed Brown test: 418742ee2edSMatthew G. Knepley # -dm_refine 6 -convest_num_refine 3 get L_2 convergence rate: [1.0] 419742ee2edSMatthew G. Knepley suffix: 2d_p1_tconv_2 420742ee2edSMatthew G. Knepley requires: triangle 421742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 422742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 423742ee2edSMatthew G. Knepley test: 424742ee2edSMatthew 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 425742ee2edSMatthew G. Knepley suffix: 2d_p1_exp_tconv_2 426742ee2edSMatthew G. Knepley requires: triangle 427742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \ 428742ee2edSMatthew G. Knepley -ts_type euler -ts_max_steps 4 -ts_dt 1e-4 -lumped 0 -mass_pc_type lu 429742ee2edSMatthew G. Knepley test: 430742ee2edSMatthew 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 431742ee2edSMatthew G. Knepley suffix: 2d_p1_exp_tconv_2_lump 432742ee2edSMatthew G. Knepley requires: triangle 433742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \ 434742ee2edSMatthew G. Knepley -ts_type euler -ts_max_steps 4 -ts_dt 1e-4 435742ee2edSMatthew G. Knepley test: 436a3d0cf85SMatthew G. Knepley suffix: 2d_p2 437c4762a1bSJed Brown requires: triangle 438a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 439a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 440c4762a1bSJed Brown test: 441a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 442a3d0cf85SMatthew G. Knepley suffix: 2d_p2_sconv 443a3d0cf85SMatthew G. Knepley requires: triangle 444a3d0cf85SMatthew G. Knepley args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 445a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 446c4762a1bSJed Brown test: 447742ee2edSMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [3.1] 448742ee2edSMatthew G. Knepley suffix: 2d_p2_sconv_2 449742ee2edSMatthew G. Knepley requires: triangle 450742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 451742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 452742ee2edSMatthew G. Knepley test: 453a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 454a3d0cf85SMatthew G. Knepley suffix: 2d_p2_tconv 455a3d0cf85SMatthew G. Knepley requires: triangle 456a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 457a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 458c4762a1bSJed Brown test: 459742ee2edSMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 460742ee2edSMatthew G. Knepley suffix: 2d_p2_tconv_2 461742ee2edSMatthew G. Knepley requires: triangle 462742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 463742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 464742ee2edSMatthew G. Knepley test: 465a3d0cf85SMatthew G. Knepley suffix: 2d_q1 46630602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 467a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 468c4762a1bSJed Brown test: 469a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 470a3d0cf85SMatthew G. Knepley suffix: 2d_q1_sconv 47130602db0SMatthew 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 \ 472a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 473c4762a1bSJed Brown test: 474a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 475a3d0cf85SMatthew G. Knepley suffix: 2d_q1_tconv 47630602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 477a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 478a3d0cf85SMatthew G. Knepley test: 479a3d0cf85SMatthew G. Knepley suffix: 2d_q2 48030602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 481a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 482a3d0cf85SMatthew G. Knepley test: 483a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 484a3d0cf85SMatthew G. Knepley suffix: 2d_q2_sconv 48530602db0SMatthew 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 \ 486a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 487a3d0cf85SMatthew G. Knepley test: 488a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 489a3d0cf85SMatthew G. Knepley suffix: 2d_q2_tconv 49030602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 491a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 492a3d0cf85SMatthew G. Knepley 493a3d0cf85SMatthew G. Knepley test: 494a3d0cf85SMatthew G. Knepley suffix: 3d_p1 495c4762a1bSJed Brown requires: ctetgen 496a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 497a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 498c4762a1bSJed Brown test: 499a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 500a3d0cf85SMatthew G. Knepley suffix: 3d_p1_sconv 501c4762a1bSJed Brown requires: ctetgen 502a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 503a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 504c4762a1bSJed Brown test: 505a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 506a3d0cf85SMatthew G. Knepley suffix: 3d_p1_tconv 507c4762a1bSJed Brown requires: ctetgen 508a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 509a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 510c4762a1bSJed Brown test: 511a3d0cf85SMatthew G. Knepley suffix: 3d_p2 512c4762a1bSJed Brown requires: ctetgen 513a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 514a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 515c4762a1bSJed Brown test: 516a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 517a3d0cf85SMatthew G. Knepley suffix: 3d_p2_sconv 518a3d0cf85SMatthew G. Knepley requires: ctetgen 519a3d0cf85SMatthew G. Knepley args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 520a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 521c4762a1bSJed Brown test: 522a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 523a3d0cf85SMatthew G. Knepley suffix: 3d_p2_tconv 524a3d0cf85SMatthew G. Knepley requires: ctetgen 525a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 526a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 527c4762a1bSJed Brown test: 528a3d0cf85SMatthew G. Knepley suffix: 3d_q1 52930602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 530a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 531c4762a1bSJed Brown test: 532a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 533a3d0cf85SMatthew G. Knepley suffix: 3d_q1_sconv 53430602db0SMatthew 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 \ 535a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 536a3d0cf85SMatthew G. Knepley test: 537a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 538a3d0cf85SMatthew G. Knepley suffix: 3d_q1_tconv 53930602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 540a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 541a3d0cf85SMatthew G. Knepley test: 542a3d0cf85SMatthew G. Knepley suffix: 3d_q2 54330602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 544a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 545a3d0cf85SMatthew G. Knepley test: 546a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 547a3d0cf85SMatthew G. Knepley suffix: 3d_q2_sconv 54830602db0SMatthew 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 \ 549a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 550a3d0cf85SMatthew G. Knepley test: 551a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 552a3d0cf85SMatthew G. Knepley suffix: 3d_q2_tconv 55330602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 554a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 555a3d0cf85SMatthew G. Knepley 556a3d0cf85SMatthew G. Knepley test: 557a3d0cf85SMatthew 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 558a3d0cf85SMatthew G. Knepley suffix: egads_sphere 559a3d0cf85SMatthew G. Knepley requires: egads ctetgen 56030602db0SMatthew G. Knepley args: -sol_type quadratic_linear \ 56130602db0SMatthew 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 \ 562a3d0cf85SMatthew G. Knepley -temp_petscspace_degree 2 -dmts_check .0001 \ 563a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 564c4762a1bSJed Brown 565c4762a1bSJed Brown TEST*/ 566