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 169371c9d4SSatish Balay typedef enum { 179371c9d4SSatish Balay SOL_QUADRATIC_LINEAR, 189371c9d4SSatish Balay SOL_QUADRATIC_TRIG, 199371c9d4SSatish Balay SOL_TRIG_LINEAR, 209371c9d4SSatish Balay SOL_TRIG_TRIG, 219371c9d4SSatish Balay NUM_SOLUTION_TYPES 229371c9d4SSatish Balay } SolutionType; 23742ee2edSMatthew G. Knepley const char *solutionTypes[NUM_SOLUTION_TYPES + 1] = {"quadratic_linear", "quadratic_trig", "trig_linear", "trig_trig", "unknown"}; 24a3d0cf85SMatthew G. Knepley 25c4762a1bSJed Brown typedef struct { 26a3d0cf85SMatthew G. Knepley SolutionType solType; /* Type of exact solution */ 27742ee2edSMatthew G. Knepley /* Solver setup */ 28742ee2edSMatthew G. Knepley PetscBool expTS; /* Flag for explicit timestepping */ 29742ee2edSMatthew G. Knepley PetscBool lumped; /* Lump the mass matrix */ 30c4762a1bSJed Brown } AppCtx; 31c4762a1bSJed Brown 32a3d0cf85SMatthew G. Knepley /* 33a3d0cf85SMatthew G. Knepley Exact 2D solution: 34a3d0cf85SMatthew G. Knepley u = 2t + x^2 + y^2 35742ee2edSMatthew G. Knepley u_t = 2 36742ee2edSMatthew G. Knepley \Delta u = 2 + 2 = 4 37742ee2edSMatthew G. Knepley f = 2 38a3d0cf85SMatthew G. Knepley F(u) = 2 - (2 + 2) + 2 = 0 39a3d0cf85SMatthew G. Knepley 40a3d0cf85SMatthew G. Knepley Exact 3D solution: 41a3d0cf85SMatthew G. Knepley u = 3t + x^2 + y^2 + z^2 42a3d0cf85SMatthew G. Knepley F(u) = 3 - (2 + 2 + 2) + 3 = 0 43a3d0cf85SMatthew G. Knepley */ 44d71ae5a4SJacob Faibussowitsch static PetscErrorCode mms_quad_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 45d71ae5a4SJacob Faibussowitsch { 46c4762a1bSJed Brown PetscInt d; 47c4762a1bSJed Brown 48c4762a1bSJed Brown *u = dim * time; 49c4762a1bSJed Brown for (d = 0; d < dim; ++d) *u += x[d] * x[d]; 50*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 51c4762a1bSJed Brown } 52c4762a1bSJed Brown 53d71ae5a4SJacob Faibussowitsch static PetscErrorCode mms_quad_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 54d71ae5a4SJacob Faibussowitsch { 55a3d0cf85SMatthew G. Knepley *u = dim; 56*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 57a3d0cf85SMatthew G. Knepley } 58a3d0cf85SMatthew G. Knepley 59d71ae5a4SJacob Faibussowitsch static void f0_quad_lin_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 60d71ae5a4SJacob Faibussowitsch { 61742ee2edSMatthew G. Knepley f0[0] = -(PetscScalar)dim; 62742ee2edSMatthew G. Knepley } 63d71ae5a4SJacob Faibussowitsch static void f0_quad_lin(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 64d71ae5a4SJacob Faibussowitsch { 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 */ 79d71ae5a4SJacob Faibussowitsch static PetscErrorCode mms_quad_trig(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 80d71ae5a4SJacob Faibussowitsch { 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]; 85*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 86a3d0cf85SMatthew G. Knepley } 87a3d0cf85SMatthew G. Knepley 88d71ae5a4SJacob Faibussowitsch static PetscErrorCode mms_quad_trig_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 89d71ae5a4SJacob Faibussowitsch { 90a3d0cf85SMatthew G. Knepley *u = -dim * PetscSinReal(time); 91*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 92a3d0cf85SMatthew G. Knepley } 93a3d0cf85SMatthew G. Knepley 94d71ae5a4SJacob Faibussowitsch static void f0_quad_trig_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 95d71ae5a4SJacob Faibussowitsch { 96742ee2edSMatthew G. Knepley f0[0] = -dim * (PetscSinReal(t) + 2.0); 97742ee2edSMatthew G. Knepley } 98d71ae5a4SJacob Faibussowitsch static void f0_quad_trig(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 99d71ae5a4SJacob Faibussowitsch { 100742ee2edSMatthew G. Knepley PetscScalar exp[1] = {0.}; 101742ee2edSMatthew 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); 102742ee2edSMatthew G. Knepley f0[0] = u_t[0] - exp[0]; 103a3d0cf85SMatthew G. Knepley } 104a3d0cf85SMatthew G. Knepley 105a3d0cf85SMatthew G. Knepley /* 106a3d0cf85SMatthew G. Knepley Exact 2D solution: 107a3d0cf85SMatthew G. Knepley u = 2\pi^2 t + cos(\pi x) + cos(\pi y) 108a3d0cf85SMatthew 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 109a3d0cf85SMatthew G. Knepley 110a3d0cf85SMatthew G. Knepley Exact 3D solution: 111a3d0cf85SMatthew G. Knepley u = 3\pi^2 t + cos(\pi x) + cos(\pi y) + cos(\pi z) 112a3d0cf85SMatthew 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 113a3d0cf85SMatthew G. Knepley */ 114d71ae5a4SJacob Faibussowitsch static PetscErrorCode mms_trig_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 115d71ae5a4SJacob Faibussowitsch { 116a3d0cf85SMatthew G. Knepley PetscInt d; 117a3d0cf85SMatthew G. Knepley 118a3d0cf85SMatthew G. Knepley *u = dim * PetscSqr(PETSC_PI) * time; 119a3d0cf85SMatthew G. Knepley for (d = 0; d < dim; ++d) *u += PetscCosReal(PETSC_PI * x[d]); 120*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 121a3d0cf85SMatthew G. Knepley } 122a3d0cf85SMatthew G. Knepley 123d71ae5a4SJacob Faibussowitsch static PetscErrorCode mms_trig_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 124d71ae5a4SJacob Faibussowitsch { 125a3d0cf85SMatthew G. Knepley *u = dim * PetscSqr(PETSC_PI); 126*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 127a3d0cf85SMatthew G. Knepley } 128a3d0cf85SMatthew G. Knepley 129d71ae5a4SJacob Faibussowitsch static void f0_trig_lin(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 130d71ae5a4SJacob Faibussowitsch { 131a3d0cf85SMatthew G. Knepley PetscInt d; 132a3d0cf85SMatthew G. Knepley f0[0] = u_t[0]; 133a3d0cf85SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += PetscSqr(PETSC_PI) * (PetscCosReal(PETSC_PI * x[d]) - 1.0); 134a3d0cf85SMatthew G. Knepley } 135a3d0cf85SMatthew G. Knepley 136742ee2edSMatthew G. Knepley /* 137742ee2edSMatthew G. Knepley Exact 2D solution: 138742ee2edSMatthew G. Knepley u = pi^2 cos(t) + cos(\pi x) + cos(\pi y) 139742ee2edSMatthew G. Knepley u_t = -pi^2 sin(t) 140742ee2edSMatthew G. Knepley \Delta u = -\pi^2 (cos(\pi x) + cos(\pi y)) 141742ee2edSMatthew G. Knepley f = pi^2 sin(t) - \pi^2 (cos(\pi x) + cos(\pi y)) 142742ee2edSMatthew 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 143742ee2edSMatthew G. Knepley 144742ee2edSMatthew G. Knepley Exact 3D solution: 145742ee2edSMatthew G. Knepley u = pi^2 cos(t) + cos(\pi x) + cos(\pi y) + cos(\pi z) 146742ee2edSMatthew G. Knepley u_t = -pi^2 sin(t) 147742ee2edSMatthew G. Knepley \Delta u = -\pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z)) 148742ee2edSMatthew G. Knepley f = pi^2 sin(t) - \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z)) 149742ee2edSMatthew 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 150742ee2edSMatthew G. Knepley */ 151d71ae5a4SJacob Faibussowitsch static PetscErrorCode mms_trig_trig(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 152d71ae5a4SJacob Faibussowitsch { 153742ee2edSMatthew G. Knepley PetscInt d; 154742ee2edSMatthew G. Knepley 155742ee2edSMatthew G. Knepley *u = PetscSqr(PETSC_PI) * PetscCosReal(time); 156742ee2edSMatthew G. Knepley for (d = 0; d < dim; ++d) *u += PetscCosReal(PETSC_PI * x[d]); 157*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 158742ee2edSMatthew G. Knepley } 159742ee2edSMatthew G. Knepley 160d71ae5a4SJacob Faibussowitsch static PetscErrorCode mms_trig_trig_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 161d71ae5a4SJacob Faibussowitsch { 162742ee2edSMatthew G. Knepley *u = -PetscSqr(PETSC_PI) * PetscSinReal(time); 163*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 164742ee2edSMatthew G. Knepley } 165742ee2edSMatthew G. Knepley 166d71ae5a4SJacob Faibussowitsch static void f0_trig_trig_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 167d71ae5a4SJacob Faibussowitsch { 168742ee2edSMatthew G. Knepley PetscInt d; 169742ee2edSMatthew G. Knepley f0[0] -= PetscSqr(PETSC_PI) * PetscSinReal(t); 170742ee2edSMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += PetscSqr(PETSC_PI) * PetscCosReal(PETSC_PI * x[d]); 171742ee2edSMatthew G. Knepley } 172d71ae5a4SJacob Faibussowitsch static void f0_trig_trig(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 173d71ae5a4SJacob Faibussowitsch { 174742ee2edSMatthew G. Knepley PetscScalar exp[1] = {0.}; 175742ee2edSMatthew 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); 176742ee2edSMatthew G. Knepley f0[0] = u_t[0] - exp[0]; 177742ee2edSMatthew G. Knepley } 178742ee2edSMatthew G. Knepley 179d71ae5a4SJacob Faibussowitsch static void f1_temp_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 180d71ae5a4SJacob Faibussowitsch { 181742ee2edSMatthew G. Knepley PetscInt d; 182742ee2edSMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = -u_x[d]; 183742ee2edSMatthew G. Knepley } 184d71ae5a4SJacob Faibussowitsch static void f1_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 185d71ae5a4SJacob Faibussowitsch { 186c4762a1bSJed Brown PetscInt d; 187a3d0cf85SMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 188c4762a1bSJed Brown } 189c4762a1bSJed Brown 190d71ae5a4SJacob Faibussowitsch static void g3_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 191d71ae5a4SJacob Faibussowitsch { 192c4762a1bSJed Brown PetscInt d; 193a3d0cf85SMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196d71ae5a4SJacob Faibussowitsch static void g0_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 197d71ae5a4SJacob Faibussowitsch { 198c4762a1bSJed Brown g0[0] = u_tShift * 1.0; 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 202d71ae5a4SJacob Faibussowitsch { 203a3d0cf85SMatthew G. Knepley PetscInt sol; 204c4762a1bSJed Brown 205c4762a1bSJed Brown PetscFunctionBeginUser; 206a3d0cf85SMatthew G. Knepley options->solType = SOL_QUADRATIC_LINEAR; 207742ee2edSMatthew G. Knepley options->expTS = PETSC_FALSE; 208742ee2edSMatthew G. Knepley options->lumped = PETSC_TRUE; 209c4762a1bSJed Brown 210d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Heat Equation Options", "DMPLEX"); 2119566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-sol_type", "Type of exact solution", "ex45.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL)); 212a3d0cf85SMatthew G. Knepley options->solType = (SolutionType)sol; 2139566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-explicit", "Use explicit timestepping", "ex45.c", options->expTS, &options->expTS, NULL)); 2149566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-lumped", "Lump the mass matrix", "ex45.c", options->lumped, &options->lumped, NULL)); 215d0609cedSBarry Smith PetscOptionsEnd(); 216*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 219d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx) 220d71ae5a4SJacob Faibussowitsch { 221c4762a1bSJed Brown PetscFunctionBeginUser; 2229566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 2239566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 2249566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2259566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 226*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown 229d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) 230d71ae5a4SJacob Faibussowitsch { 231a3d0cf85SMatthew G. Knepley PetscDS ds; 23245480ffeSMatthew G. Knepley DMLabel label; 233c4762a1bSJed Brown const PetscInt id = 1; 234c4762a1bSJed Brown 235c4762a1bSJed Brown PetscFunctionBeginUser; 2369566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 2379566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2389566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_temp, NULL, NULL, g3_temp)); 239a3d0cf85SMatthew G. Knepley switch (ctx->solType) { 240a3d0cf85SMatthew G. Knepley case SOL_QUADRATIC_LINEAR: 2419566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_quad_lin, f1_temp)); 2429566063dSJacob Faibussowitsch PetscCall(PetscDSSetRHSResidual(ds, 0, f0_quad_lin_exp, f1_temp_exp)); 2439566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms_quad_lin, ctx)); 2449566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_lin_t, ctx)); 2459566063dSJacob 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)); 246a3d0cf85SMatthew G. Knepley break; 247a3d0cf85SMatthew G. Knepley case SOL_QUADRATIC_TRIG: 2489566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_quad_trig, f1_temp)); 2499566063dSJacob Faibussowitsch PetscCall(PetscDSSetRHSResidual(ds, 0, f0_quad_trig_exp, f1_temp_exp)); 2509566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms_quad_trig, ctx)); 2519566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_trig_t, ctx)); 2529566063dSJacob 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)); 253a3d0cf85SMatthew G. Knepley break; 254a3d0cf85SMatthew G. Knepley case SOL_TRIG_LINEAR: 2559566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_lin, f1_temp)); 2569566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms_trig_lin, ctx)); 2579566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_lin_t, ctx)); 2589566063dSJacob 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)); 259a3d0cf85SMatthew G. Knepley break; 260742ee2edSMatthew G. Knepley case SOL_TRIG_TRIG: 2619566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_trig, f1_temp)); 2629566063dSJacob Faibussowitsch PetscCall(PetscDSSetRHSResidual(ds, 0, f0_trig_trig_exp, f1_temp_exp)); 2639566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms_trig_trig, ctx)); 2649566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_trig_t, ctx)); 265742ee2edSMatthew G. Knepley break; 266d71ae5a4SJacob Faibussowitsch default: 267d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%d)", solutionTypes[PetscMin(ctx->solType, NUM_SOLUTION_TYPES)], ctx->solType); 268a3d0cf85SMatthew G. Knepley } 269*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 272d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx) 273d71ae5a4SJacob Faibussowitsch { 274c4762a1bSJed Brown DM cdm = dm; 275c4762a1bSJed Brown PetscFE fe; 276a3d0cf85SMatthew G. Knepley DMPolytopeType ct; 277a3d0cf85SMatthew G. Knepley PetscBool simplex; 278a3d0cf85SMatthew G. Knepley PetscInt dim, cStart; 279c4762a1bSJed Brown 280c4762a1bSJed Brown PetscFunctionBeginUser; 2819566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2829566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 2839566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 284a3d0cf85SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 285c4762a1bSJed Brown /* Create finite element */ 2869566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "temp_", -1, &fe)); 2879566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "temperature")); 288c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 2899566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 2909566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 291742ee2edSMatthew G. Knepley if (ctx->expTS) { 292742ee2edSMatthew G. Knepley PetscDS ds; 293742ee2edSMatthew G. Knepley 2949566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2959566063dSJacob Faibussowitsch PetscCall(PetscDSSetImplicit(ds, 0, PETSC_FALSE)); 296742ee2edSMatthew G. Knepley } 2979566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, ctx)); 298c4762a1bSJed Brown while (cdm) { 2999566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 3009566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 301c4762a1bSJed Brown } 3029566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 303*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 304c4762a1bSJed Brown } 305c4762a1bSJed Brown 306d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetInitialConditions(TS ts, Vec u) 307d71ae5a4SJacob Faibussowitsch { 308a3d0cf85SMatthew G. Knepley DM dm; 309a3d0cf85SMatthew G. Knepley PetscReal t; 310a3d0cf85SMatthew G. Knepley 3117510d9b0SBarry Smith PetscFunctionBeginUser; 3129566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 3139566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 3149566063dSJacob Faibussowitsch PetscCall(DMComputeExactSolution(dm, t, u, NULL)); 315*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 316a3d0cf85SMatthew G. Knepley } 317a3d0cf85SMatthew G. Knepley 318d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 319d71ae5a4SJacob Faibussowitsch { 320c4762a1bSJed Brown DM dm; 321c4762a1bSJed Brown TS ts; 322a3d0cf85SMatthew G. Knepley Vec u; 323a3d0cf85SMatthew G. Knepley AppCtx ctx; 324c4762a1bSJed Brown 325327415f7SBarry Smith PetscFunctionBeginUser; 3269566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3279566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 3289566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx)); 3299566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &ctx)); 3309566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &ctx)); 331c4762a1bSJed Brown 3329566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 3339566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm)); 3349566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx)); 335742ee2edSMatthew G. Knepley if (ctx.expTS) { 3369566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFEM, &ctx)); 3379566063dSJacob Faibussowitsch if (ctx.lumped) PetscCall(DMTSCreateRHSMassMatrixLumped(dm)); 3389566063dSJacob Faibussowitsch else PetscCall(DMTSCreateRHSMassMatrix(dm)); 339742ee2edSMatthew G. Knepley } else { 3409566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx)); 3419566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx)); 342742ee2edSMatthew G. Knepley } 3439566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 3449566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 3459566063dSJacob Faibussowitsch PetscCall(TSSetComputeInitialCondition(ts, SetInitialConditions)); 346c4762a1bSJed Brown 3479566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 3489566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u)); 3499566063dSJacob Faibussowitsch PetscCall(SetInitialConditions(ts, u)); 3509566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "temperature")); 3519566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 3529566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u)); 3539566063dSJacob Faibussowitsch if (ctx.expTS) PetscCall(DMTSDestroyRHSMassMatrix(dm)); 354c4762a1bSJed Brown 3559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 3569566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 3579566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3589566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 359b122ec5aSJacob Faibussowitsch return 0; 360c4762a1bSJed Brown } 361c4762a1bSJed Brown 362c4762a1bSJed Brown /*TEST 363c4762a1bSJed Brown 364c4762a1bSJed Brown test: 365a3d0cf85SMatthew G. Knepley suffix: 2d_p1 366c4762a1bSJed Brown requires: triangle 367a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 368a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 369c4762a1bSJed Brown test: 370742ee2edSMatthew G. Knepley suffix: 2d_p1_exp 371742ee2edSMatthew G. Knepley requires: triangle 372742ee2edSMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -explicit \ 373742ee2edSMatthew G. Knepley -ts_type euler -ts_max_steps 4 -ts_dt 1e-3 -ts_monitor_error 374742ee2edSMatthew G. Knepley test: 375a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 376a3d0cf85SMatthew G. Knepley suffix: 2d_p1_sconv 377c4762a1bSJed Brown requires: triangle 378a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 379a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 380c4762a1bSJed Brown test: 381742ee2edSMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.1] 382742ee2edSMatthew G. Knepley suffix: 2d_p1_sconv_2 383742ee2edSMatthew G. Knepley requires: triangle 384742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 385742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 1e-6 -snes_error_if_not_converged -pc_type lu 386742ee2edSMatthew G. Knepley test: 387a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 388a3d0cf85SMatthew G. Knepley suffix: 2d_p1_tconv 389c4762a1bSJed Brown requires: triangle 390a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 391a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 392c4762a1bSJed Brown test: 393742ee2edSMatthew G. Knepley # -dm_refine 6 -convest_num_refine 3 get L_2 convergence rate: [1.0] 394742ee2edSMatthew G. Knepley suffix: 2d_p1_tconv_2 395742ee2edSMatthew G. Knepley requires: triangle 396742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 397742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 398742ee2edSMatthew G. Knepley test: 399742ee2edSMatthew 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 400742ee2edSMatthew G. Knepley suffix: 2d_p1_exp_tconv_2 401742ee2edSMatthew G. Knepley requires: triangle 402742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \ 403742ee2edSMatthew G. Knepley -ts_type euler -ts_max_steps 4 -ts_dt 1e-4 -lumped 0 -mass_pc_type lu 404742ee2edSMatthew G. Knepley test: 405742ee2edSMatthew 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 406742ee2edSMatthew G. Knepley suffix: 2d_p1_exp_tconv_2_lump 407742ee2edSMatthew G. Knepley requires: triangle 408742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \ 409742ee2edSMatthew G. Knepley -ts_type euler -ts_max_steps 4 -ts_dt 1e-4 410742ee2edSMatthew G. Knepley test: 411a3d0cf85SMatthew G. Knepley suffix: 2d_p2 412c4762a1bSJed Brown requires: triangle 413a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 414a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 415c4762a1bSJed Brown test: 416a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 417a3d0cf85SMatthew G. Knepley suffix: 2d_p2_sconv 418a3d0cf85SMatthew G. Knepley requires: triangle 419a3d0cf85SMatthew G. Knepley args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 420a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 421c4762a1bSJed Brown test: 422742ee2edSMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [3.1] 423742ee2edSMatthew G. Knepley suffix: 2d_p2_sconv_2 424742ee2edSMatthew G. Knepley requires: triangle 425742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 426742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 427742ee2edSMatthew G. Knepley test: 428a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 429a3d0cf85SMatthew G. Knepley suffix: 2d_p2_tconv 430a3d0cf85SMatthew G. Knepley requires: triangle 431a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 432a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 433c4762a1bSJed Brown test: 434742ee2edSMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 435742ee2edSMatthew G. Knepley suffix: 2d_p2_tconv_2 436742ee2edSMatthew G. Knepley requires: triangle 437742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 438742ee2edSMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 439742ee2edSMatthew G. Knepley test: 440a3d0cf85SMatthew G. Knepley suffix: 2d_q1 44130602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 442a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 443c4762a1bSJed Brown test: 444a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 445a3d0cf85SMatthew G. Knepley suffix: 2d_q1_sconv 44630602db0SMatthew 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 \ 447a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 448c4762a1bSJed Brown test: 449a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 450a3d0cf85SMatthew G. Knepley suffix: 2d_q1_tconv 45130602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 452a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 453a3d0cf85SMatthew G. Knepley test: 454a3d0cf85SMatthew G. Knepley suffix: 2d_q2 45530602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 456a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 457a3d0cf85SMatthew G. Knepley test: 458a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 459a3d0cf85SMatthew G. Knepley suffix: 2d_q2_sconv 46030602db0SMatthew 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 \ 461a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 462a3d0cf85SMatthew G. Knepley test: 463a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 464a3d0cf85SMatthew G. Knepley suffix: 2d_q2_tconv 46530602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 466a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 467a3d0cf85SMatthew G. Knepley 468a3d0cf85SMatthew G. Knepley test: 469a3d0cf85SMatthew G. Knepley suffix: 3d_p1 470c4762a1bSJed Brown requires: ctetgen 471a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 472a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 473c4762a1bSJed Brown test: 474a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 475a3d0cf85SMatthew G. Knepley suffix: 3d_p1_sconv 476c4762a1bSJed Brown requires: ctetgen 477a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 478a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 479c4762a1bSJed Brown test: 480a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 481a3d0cf85SMatthew G. Knepley suffix: 3d_p1_tconv 482c4762a1bSJed Brown requires: ctetgen 483a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -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 485c4762a1bSJed Brown test: 486a3d0cf85SMatthew G. Knepley suffix: 3d_p2 487c4762a1bSJed Brown requires: ctetgen 488a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 489a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 490c4762a1bSJed Brown test: 491a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 492a3d0cf85SMatthew G. Knepley suffix: 3d_p2_sconv 493a3d0cf85SMatthew G. Knepley requires: ctetgen 494a3d0cf85SMatthew G. Knepley args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 495a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 496c4762a1bSJed Brown test: 497a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 498a3d0cf85SMatthew G. Knepley suffix: 3d_p2_tconv 499a3d0cf85SMatthew G. Knepley requires: ctetgen 500a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 501a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 502c4762a1bSJed Brown test: 503a3d0cf85SMatthew G. Knepley suffix: 3d_q1 50430602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \ 505a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 506c4762a1bSJed Brown test: 507a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9] 508a3d0cf85SMatthew G. Knepley suffix: 3d_q1_sconv 50930602db0SMatthew 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 \ 510a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu 511a3d0cf85SMatthew G. Knepley test: 512a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2] 513a3d0cf85SMatthew G. Knepley suffix: 3d_q1_tconv 51430602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \ 515a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 516a3d0cf85SMatthew G. Knepley test: 517a3d0cf85SMatthew G. Knepley suffix: 3d_q2 51830602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \ 519a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 520a3d0cf85SMatthew G. Knepley test: 521a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9] 522a3d0cf85SMatthew G. Knepley suffix: 3d_q2_sconv 52330602db0SMatthew 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 \ 524a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu 525a3d0cf85SMatthew G. Knepley test: 526a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0] 527a3d0cf85SMatthew G. Knepley suffix: 3d_q2_tconv 52830602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \ 529a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 530a3d0cf85SMatthew G. Knepley 531a3d0cf85SMatthew G. Knepley test: 532a3d0cf85SMatthew 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 533a3d0cf85SMatthew G. Knepley suffix: egads_sphere 534a3d0cf85SMatthew G. Knepley requires: egads ctetgen 53530602db0SMatthew G. Knepley args: -sol_type quadratic_linear \ 53630602db0SMatthew 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 \ 537a3d0cf85SMatthew G. Knepley -temp_petscspace_degree 2 -dmts_check .0001 \ 538a3d0cf85SMatthew G. Knepley -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu 539c4762a1bSJed Brown 540c4762a1bSJed Brown TEST*/ 541