130602db0SMatthew G. Knepley const char help[] = "A test of H-div conforming discretizations on different cell types.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown #include <petscds.h> 5c4762a1bSJed Brown #include <petscsnes.h> 6c4762a1bSJed Brown #include <petscconvest.h> 7c4762a1bSJed Brown #include <petscfe.h> 8c4762a1bSJed Brown #include <petsc/private/petscfeimpl.h> 9c4762a1bSJed Brown 10c4762a1bSJed Brown /* 110e3d61c9SBarry Smith We are using the system 120e3d61c9SBarry Smith 130e3d61c9SBarry Smith \vec{u} = \vec{\hat{u}} 140e3d61c9SBarry Smith p = \div{\vec{u}} in low degree approximation space 150e3d61c9SBarry Smith d = \div{\vec{u}} - p == 0 in higher degree approximation space 160e3d61c9SBarry Smith 170e3d61c9SBarry Smith That is, we are using the field d to compute the error between \div{\vec{u}} 180e3d61c9SBarry Smith computed in a space 1 degree higher than p and the value of p which is 190e3d61c9SBarry Smith \div{u} computed in the low degree space. If H-div 200e3d61c9SBarry Smith elements are implemented correctly then this should be identically zero since 210e3d61c9SBarry Smith the divergence of a function in H(div) should be exactly representable in L_2 220e3d61c9SBarry Smith by definition. 23c4762a1bSJed Brown */ 249371c9d4SSatish Balay static PetscErrorCode zero_func(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 25c4762a1bSJed Brown PetscInt c; 26c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = 0; 27c4762a1bSJed Brown return 0; 28c4762a1bSJed Brown } 29c4762a1bSJed Brown /* Linear Exact Functions 30c4762a1bSJed Brown \vec{u} = \vec{x}; 31c4762a1bSJed Brown p = dim; 32c4762a1bSJed Brown */ 339371c9d4SSatish Balay static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 34c4762a1bSJed Brown PetscInt c; 35c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = x[c]; 36c4762a1bSJed Brown return 0; 37c4762a1bSJed Brown } 389371c9d4SSatish Balay static PetscErrorCode linear_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 39c4762a1bSJed Brown u[0] = dim; 40c4762a1bSJed Brown return 0; 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Sinusoidal Exact Functions 44c4762a1bSJed Brown * u_i = \sin{2*\pi*x_i} 45c4762a1bSJed Brown * p = \Sum_{i=1}^{dim} 2*\pi*cos{2*\pi*x_i} 46c4762a1bSJed Brown * */ 47c4762a1bSJed Brown 489371c9d4SSatish Balay static PetscErrorCode sinusoid_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 49c4762a1bSJed Brown PetscInt c; 50c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = PetscSinReal(2 * PETSC_PI * x[c]); 51c4762a1bSJed Brown return 0; 52c4762a1bSJed Brown } 539371c9d4SSatish Balay static PetscErrorCode sinusoid_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 54c4762a1bSJed Brown PetscInt d; 55c4762a1bSJed Brown u[0] = 0; 56c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[0] += 2 * PETSC_PI * PetscCosReal(2 * PETSC_PI * x[d]); 57c4762a1bSJed Brown return 0; 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Pointwise residual for u = u*. Need one of these for each possible u* */ 619371c9d4SSatish Balay static void f0_v_linear(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[]) { 62c4762a1bSJed Brown PetscInt i; 63c4762a1bSJed Brown PetscScalar *u_rhs; 64c4762a1bSJed Brown 65c4762a1bSJed Brown PetscCalloc1(dim, &u_rhs); 66c4762a1bSJed Brown (void)linear_u(dim, t, x, dim, u_rhs, NULL); 67c4762a1bSJed Brown for (i = 0; i < dim; ++i) f0[i] = u[uOff[0] + i] - u_rhs[i]; 68c4762a1bSJed Brown PetscFree(u_rhs); 69c4762a1bSJed Brown } 70c4762a1bSJed Brown 719371c9d4SSatish Balay static void f0_v_sinusoid(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[]) { 72c4762a1bSJed Brown PetscInt i; 73c4762a1bSJed Brown PetscScalar *u_rhs; 74c4762a1bSJed Brown 75c4762a1bSJed Brown PetscCalloc1(dim, &u_rhs); 76c4762a1bSJed Brown (void)sinusoid_u(dim, t, x, dim, u_rhs, NULL); 77c4762a1bSJed Brown for (i = 0; i < dim; ++i) f0[i] = u[uOff[0] + i] - u_rhs[i]; 78c4762a1bSJed Brown PetscFree(u_rhs); 79c4762a1bSJed Brown } 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Residual function for enforcing p = \div{u}. */ 829371c9d4SSatish Balay static void f0_q(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[]) { 83c4762a1bSJed Brown PetscInt i; 84c4762a1bSJed Brown PetscScalar divu; 85c4762a1bSJed Brown 86c4762a1bSJed Brown divu = 0.; 87c4762a1bSJed Brown for (i = 0; i < dim; ++i) divu += u_x[uOff_x[0] + i * dim + i]; 88c4762a1bSJed Brown f0[0] = u[uOff[1]] - divu; 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Residual function for p_err = \div{u} - p. */ 929371c9d4SSatish Balay static void f0_w(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[]) { 93c4762a1bSJed Brown PetscInt i; 94c4762a1bSJed Brown PetscScalar divu; 95c4762a1bSJed Brown 96c4762a1bSJed Brown divu = 0.; 97c4762a1bSJed Brown for (i = 0; i < dim; ++i) divu += u_x[uOff_x[0] + i * dim + i]; 98c4762a1bSJed Brown f0[0] = u[uOff[2]] - u[uOff[1]] + divu; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* Boundary residual for the embedding system. Need one for each form of 102c4762a1bSJed Brown * solution. These enforce u = \hat{u} at the boundary. */ 1039371c9d4SSatish Balay static void f0_bd_u_sinusoid(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[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 104c4762a1bSJed Brown PetscInt d; 105c4762a1bSJed Brown PetscScalar *u_rhs; 106c4762a1bSJed Brown PetscCalloc1(dim, &u_rhs); 107c4762a1bSJed Brown (void)sinusoid_u(dim, t, x, dim, u_rhs, NULL); 108c4762a1bSJed Brown 109c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] = u_rhs[d]; 110c4762a1bSJed Brown PetscFree(u_rhs); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 1139371c9d4SSatish Balay static void f0_bd_u_linear(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[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 114c4762a1bSJed Brown PetscInt d; 115c4762a1bSJed Brown PetscScalar *u_rhs; 116c4762a1bSJed Brown PetscCalloc1(dim, &u_rhs); 117c4762a1bSJed Brown (void)linear_u(dim, t, x, dim, u_rhs, NULL); 118c4762a1bSJed Brown 119c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] = u_rhs[d]; 120c4762a1bSJed Brown PetscFree(u_rhs); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown /* Jacobian functions. For the following, v is the test function associated with 123c4762a1bSJed Brown * u, q the test function associated with p, and w the test function associated 124c4762a1bSJed Brown * with d. */ 125c4762a1bSJed Brown /* <v, u> */ 1269371c9d4SSatish Balay static void g0_vu(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[]) { 127c4762a1bSJed Brown PetscInt c; 128c4762a1bSJed Brown for (c = 0; c < dim; ++c) g0[c * dim + c] = 1.0; 129c4762a1bSJed Brown } 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* <q, p> */ 1329371c9d4SSatish Balay static void g0_qp(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[]) { 133c4762a1bSJed Brown PetscInt d; 134c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[d * dim + d] = 1.0; 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* -<q, \div{u}> For the embedded system. This is different from the method of 138c4762a1bSJed Brown * manufactured solution because instead of computing <q,\div{u}> - <q,f> we 139c4762a1bSJed Brown * need <q,p> - <q,\div{u}.*/ 1409371c9d4SSatish Balay static void g1_qu(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 g1[]) { 141c4762a1bSJed Brown PetscInt d; 142c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d * dim + d] = -1.0; 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* <w, p> This is only used by the embedded system. Where we need to compute 146c4762a1bSJed Brown * <w,d> - <w,p> + <w, \div{u}>*/ 1479371c9d4SSatish Balay static void g0_wp(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[]) { 148c4762a1bSJed Brown PetscInt d; 149c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[d * dim + d] = -1.0; 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* <w, d> */ 1539371c9d4SSatish Balay static void g0_wd(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[]) { 154c4762a1bSJed Brown PetscInt c; 155c4762a1bSJed Brown for (c = 0; c < dim; ++c) g0[c * dim + c] = 1.0; 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* <w, \div{u}> for the embedded system. */ 1599371c9d4SSatish Balay static void g1_wu(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 g1[]) { 160c4762a1bSJed Brown PetscInt d; 161c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; 162c4762a1bSJed Brown } 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* Enum and string array for selecting mesh perturbation options */ 1659371c9d4SSatish Balay typedef enum { 1669371c9d4SSatish Balay NONE = 0, 1679371c9d4SSatish Balay PERTURB = 1, 1689371c9d4SSatish Balay SKEW = 2, 1699371c9d4SSatish Balay SKEW_PERTURB = 3 1709371c9d4SSatish Balay } Transform; 171c4762a1bSJed Brown const char *const TransformTypes[] = {"none", "perturb", "skew", "skew_perturb", "Perturbation", "", NULL}; 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* Enum and string array for selecting the form of the exact solution*/ 1749371c9d4SSatish Balay typedef enum { 1759371c9d4SSatish Balay LINEAR = 0, 1769371c9d4SSatish Balay SINUSOIDAL = 1 1779371c9d4SSatish Balay } Solution; 178c4762a1bSJed Brown const char *const SolutionTypes[] = {"linear", "sinusoidal", "Solution", "", NULL}; 179c4762a1bSJed Brown 1809371c9d4SSatish Balay typedef struct { 181c4762a1bSJed Brown Transform mesh_transform; 182c4762a1bSJed Brown Solution sol_form; 183c4762a1bSJed Brown } UserCtx; 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* Process command line options and initialize the UserCtx struct */ 1869371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, UserCtx *user) { 187c4762a1bSJed Brown PetscFunctionBegin; 188c4762a1bSJed Brown /* Default to 2D, unperturbed triangle mesh and Linear solution.*/ 189c4762a1bSJed Brown user->mesh_transform = NONE; 190c4762a1bSJed Brown user->sol_form = LINEAR; 191c4762a1bSJed Brown 192d0609cedSBarry Smith PetscOptionsBegin(comm, "", "H-div Test Options", "DMPLEX"); 1939566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-mesh_transform", "Method used to perturb the mesh vertices. Options are skew, perturb, skew_perturb,or none", "ex39.c", TransformTypes, (PetscEnum)user->mesh_transform, (PetscEnum *)&user->mesh_transform, NULL)); 1949566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-sol_form", "Form of the exact solution. Options are Linear or Sinusoidal", "ex39.c", SolutionTypes, (PetscEnum)user->sol_form, (PetscEnum *)&user->sol_form, NULL)); 195d0609cedSBarry Smith PetscOptionsEnd(); 196c4762a1bSJed Brown PetscFunctionReturn(0); 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199c4762a1bSJed Brown /* Perturb the position of each mesh vertex by a small amount.*/ 2009371c9d4SSatish Balay static PetscErrorCode PerturbMesh(DM *mesh, PetscScalar *coordVals, PetscInt npoints, PetscInt dim) { 201c4762a1bSJed Brown PetscInt i, j, k; 202d092c84bSBrandon Whitchurch PetscReal minCoords[3], maxCoords[3], maxPert[3], randVal, amp; 203c4762a1bSJed Brown PetscRandom ran; 204c4762a1bSJed Brown 205c4762a1bSJed Brown PetscFunctionBegin; 2069566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(*mesh, &dim)); 2079566063dSJacob Faibussowitsch PetscCall(DMGetLocalBoundingBox(*mesh, minCoords, maxCoords)); 2089566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &ran)); 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* Compute something approximately equal to half an edge length. This is the 211a5b23f4aSJose E. Roman * most we can perturb points and guarantee that there won't be any topology 212c4762a1bSJed Brown * issues. */ 213d092c84bSBrandon Whitchurch for (k = 0; k < dim; ++k) maxPert[k] = 0.025 * (maxCoords[k] - minCoords[k]) / (PetscPowReal(npoints, 1. / dim) - 1); 214c4762a1bSJed Brown /* For each mesh vertex */ 215c4762a1bSJed Brown for (i = 0; i < npoints; ++i) { 216c4762a1bSJed Brown /* For each coordinate of the vertex */ 217c4762a1bSJed Brown for (j = 0; j < dim; ++j) { 218c4762a1bSJed Brown /* Generate a random amplitude in [-0.5*maxPert, 0.5*maxPert] */ 2199566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(ran, &randVal)); 220c4762a1bSJed Brown amp = maxPert[j] * (randVal - 0.5); 221c4762a1bSJed Brown /* Add the perturbation to the vertex*/ 222d092c84bSBrandon Whitchurch coordVals[dim * i + j] += amp; 223c4762a1bSJed Brown } 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown PetscRandomDestroy(&ran); 227c4762a1bSJed Brown PetscFunctionReturn(0); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown 230c4762a1bSJed Brown /* Apply a global skew transformation to the mesh. */ 2319371c9d4SSatish Balay static PetscErrorCode SkewMesh(DM *mesh, PetscScalar *coordVals, PetscInt npoints, PetscInt dim) { 232c4762a1bSJed Brown PetscInt i, j, k, l; 233c4762a1bSJed Brown PetscScalar *transMat; 234c4762a1bSJed Brown PetscScalar tmpcoord[3]; 235c4762a1bSJed Brown PetscRandom ran; 236c4762a1bSJed Brown PetscReal randVal; 237c4762a1bSJed Brown 238c4762a1bSJed Brown PetscFunctionBegin; 2399566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(dim * dim, &transMat)); 2409566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &ran)); 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* Make a matrix representing a skew transformation */ 243c4762a1bSJed Brown for (i = 0; i < dim; ++i) { 244c4762a1bSJed Brown for (j = 0; j < dim; ++j) { 2459566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(ran, &randVal)); 246d092c84bSBrandon Whitchurch if (i == j) transMat[i * dim + j] = 1.; 247c4762a1bSJed Brown else if (j < i) transMat[i * dim + j] = 2 * (j + i) * randVal; 248c4762a1bSJed Brown else transMat[i * dim + j] = 0; 249c4762a1bSJed Brown } 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252c4762a1bSJed Brown /* Multiply each coordinate vector by our tranformation.*/ 253c4762a1bSJed Brown for (i = 0; i < npoints; ++i) { 254c4762a1bSJed Brown for (j = 0; j < dim; ++j) { 255c4762a1bSJed Brown tmpcoord[j] = 0; 256c4762a1bSJed Brown for (k = 0; k < dim; ++k) tmpcoord[j] += coordVals[dim * i + k] * transMat[dim * k + j]; 257c4762a1bSJed Brown } 258c4762a1bSJed Brown for (l = 0; l < dim; ++l) coordVals[dim * i + l] = tmpcoord[l]; 259c4762a1bSJed Brown } 2609566063dSJacob Faibussowitsch PetscCall(PetscFree(transMat)); 2619566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&ran)); 262c4762a1bSJed Brown PetscFunctionReturn(0); 263c4762a1bSJed Brown } 264c4762a1bSJed Brown 265c4762a1bSJed Brown /* Accesses the mesh coordinate array and performs the transformation operations 266c4762a1bSJed Brown * specified by the user options */ 2679371c9d4SSatish Balay static PetscErrorCode TransformMesh(UserCtx *user, DM *mesh) { 268c4762a1bSJed Brown PetscInt dim, npoints; 269c4762a1bSJed Brown PetscScalar *coordVals; 270c4762a1bSJed Brown Vec coords; 271c4762a1bSJed Brown 272c4762a1bSJed Brown PetscFunctionBegin; 2739566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(*mesh, &coords)); 2749566063dSJacob Faibussowitsch PetscCall(VecGetArray(coords, &coordVals)); 2759566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coords, &npoints)); 2769566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(*mesh, &dim)); 277c4762a1bSJed Brown npoints = npoints / dim; 278c4762a1bSJed Brown 279c4762a1bSJed Brown switch (user->mesh_transform) { 2809371c9d4SSatish Balay case PERTURB: PetscCall(PerturbMesh(mesh, coordVals, npoints, dim)); break; 2819371c9d4SSatish Balay case SKEW: PetscCall(SkewMesh(mesh, coordVals, npoints, dim)); break; 282c4762a1bSJed Brown case SKEW_PERTURB: 2839566063dSJacob Faibussowitsch PetscCall(SkewMesh(mesh, coordVals, npoints, dim)); 2849566063dSJacob Faibussowitsch PetscCall(PerturbMesh(mesh, coordVals, npoints, dim)); 285c4762a1bSJed Brown break; 2869371c9d4SSatish Balay default: PetscFunctionReturn(-1); 287c4762a1bSJed Brown } 2889566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coords, &coordVals)); 2899566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(*mesh, coords)); 290c4762a1bSJed Brown PetscFunctionReturn(0); 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 2939371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, UserCtx *user, DM *mesh) { 294c4762a1bSJed Brown PetscFunctionBegin; 2959566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, mesh)); 2969566063dSJacob Faibussowitsch PetscCall(DMSetType(*mesh, DMPLEX)); 2979566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*mesh)); 298c4762a1bSJed Brown 299c4762a1bSJed Brown /* Perform any mesh transformations if specified by user */ 300*48a46eb9SPierre Jolivet if (user->mesh_transform != NONE) PetscCall(TransformMesh(user, mesh)); 301c4762a1bSJed Brown 302c4762a1bSJed Brown /* Get any other mesh options from the command line */ 3039566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*mesh, user)); 3049566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*mesh, NULL, "-dm_view")); 305c4762a1bSJed Brown PetscFunctionReturn(0); 306c4762a1bSJed Brown } 307c4762a1bSJed Brown 308c4762a1bSJed Brown /* Setup the system of equations that we wish to solve */ 3099371c9d4SSatish Balay static PetscErrorCode SetupProblem(DM dm, UserCtx *user) { 310c4762a1bSJed Brown PetscDS prob; 31145480ffeSMatthew G. Knepley DMLabel label; 312c4762a1bSJed Brown const PetscInt id = 1; 313c4762a1bSJed Brown 314c4762a1bSJed Brown PetscFunctionBegin; 3159566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 316c4762a1bSJed Brown /* All of these are independent of the user's choice of solution */ 3179566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 1, f0_q, NULL)); 3189566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 2, f0_w, NULL)); 3199566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_vu, NULL, NULL, NULL)); 3209566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL)); 3219566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_qp, NULL, NULL, NULL)); 3229566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 2, 0, NULL, g1_wu, NULL, NULL)); 3239566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 2, 1, g0_wp, NULL, NULL, NULL)); 3249566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 2, 2, g0_wd, NULL, NULL, NULL)); 325c4762a1bSJed Brown 326c4762a1bSJed Brown /* Field 2 is the error between \div{u} and pressure in a higher dimensional 327c4762a1bSJed Brown * space. If all is right this should be machine zero. */ 3289566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 2, zero_func, NULL)); 329c4762a1bSJed Brown 330c4762a1bSJed Brown switch (user->sol_form) { 331c4762a1bSJed Brown case LINEAR: 3329566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_v_linear, NULL)); 3339566063dSJacob Faibussowitsch PetscCall(PetscDSSetBdResidual(prob, 0, f0_bd_u_linear, NULL)); 3349566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 0, linear_u, NULL)); 3359566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 1, linear_p, NULL)); 336c4762a1bSJed Brown break; 337c4762a1bSJed Brown case SINUSOIDAL: 3389566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_v_sinusoid, NULL)); 3399566063dSJacob Faibussowitsch PetscCall(PetscDSSetBdResidual(prob, 0, f0_bd_u_sinusoid, NULL)); 3409566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 0, sinusoid_u, NULL)); 3419566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 1, sinusoid_p, NULL)); 342c4762a1bSJed Brown break; 3439371c9d4SSatish Balay default: PetscFunctionReturn(-1); 344c4762a1bSJed Brown } 345c4762a1bSJed Brown 3469566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 3479566063dSJacob Faibussowitsch PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL, "Boundary Integral", label, 1, &id, 0, 0, NULL, (void (*)(void))NULL, NULL, user, NULL)); 348c4762a1bSJed Brown PetscFunctionReturn(0); 349c4762a1bSJed Brown } 350c4762a1bSJed Brown 351c4762a1bSJed Brown /* Create the finite element spaces we will use for this system */ 3529371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM mesh, PetscErrorCode (*setup)(DM, UserCtx *), UserCtx *user) { 353c4762a1bSJed Brown DM cdm = mesh; 354c4762a1bSJed Brown PetscFE fevel, fepres, fedivErr; 35530602db0SMatthew G. Knepley PetscInt dim; 35630602db0SMatthew G. Knepley PetscBool simplex; 357c4762a1bSJed Brown 358c4762a1bSJed Brown PetscFunctionBegin; 3599566063dSJacob Faibussowitsch PetscCall(DMGetDimension(mesh, &dim)); 3609566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(mesh, &simplex)); 361c4762a1bSJed Brown /* Create FE objects and give them names so that options can be set from 362c4762a1bSJed Brown * command line */ 3639566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh), dim, dim, simplex, "velocity_", -1, &fevel)); 3649566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fevel, "velocity")); 365c4762a1bSJed Brown 3669566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh), dim, 1, simplex, "pressure_", -1, &fepres)); 3679566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fepres, "pressure")); 368c4762a1bSJed Brown 369d0609cedSBarry Smith PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh), dim, 1, simplex, "divErr_", -1, &fedivErr)); 3709566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fedivErr, "divErr")); 371c4762a1bSJed Brown 3729566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fevel, fepres)); 3739566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fevel, fedivErr)); 374c4762a1bSJed Brown 375c4762a1bSJed Brown /* Associate the FE objects with the mesh and setup the system */ 3769566063dSJacob Faibussowitsch PetscCall(DMSetField(mesh, 0, NULL, (PetscObject)fevel)); 3779566063dSJacob Faibussowitsch PetscCall(DMSetField(mesh, 1, NULL, (PetscObject)fepres)); 3789566063dSJacob Faibussowitsch PetscCall(DMSetField(mesh, 2, NULL, (PetscObject)fedivErr)); 3799566063dSJacob Faibussowitsch PetscCall(DMCreateDS(mesh)); 3809566063dSJacob Faibussowitsch PetscCall((*setup)(mesh, user)); 381c4762a1bSJed Brown 382c4762a1bSJed Brown while (cdm) { 3839566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(mesh, cdm)); 3849566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 385c4762a1bSJed Brown } 386c4762a1bSJed Brown 387c4762a1bSJed Brown /* The Mesh now owns the fields, so we can destroy the FEs created in this 388c4762a1bSJed Brown * function */ 3899566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fevel)); 3909566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fepres)); 3919566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fedivErr)); 3929566063dSJacob Faibussowitsch PetscCall(DMDestroy(&cdm)); 393c4762a1bSJed Brown PetscFunctionReturn(0); 394c4762a1bSJed Brown } 395c4762a1bSJed Brown 3969371c9d4SSatish Balay int main(int argc, char **argv) { 397c4762a1bSJed Brown PetscInt i; 398c4762a1bSJed Brown UserCtx user; 399c4762a1bSJed Brown DM mesh; 400c4762a1bSJed Brown SNES snes; 401c4762a1bSJed Brown Vec computed, divErr; 402c4762a1bSJed Brown PetscReal divErrNorm; 403c4762a1bSJed Brown IS *fieldIS; 404c4762a1bSJed Brown PetscBool exampleSuccess = PETSC_FALSE; 405d092c84bSBrandon Whitchurch const PetscReal errTol = 10. * PETSC_SMALL; 406c4762a1bSJed Brown 407c4762a1bSJed Brown char stdFormat[] = "L2 Norm of the Divergence Error is: %g\n H(div) elements working correctly: %s\n"; 408c4762a1bSJed Brown 409c4762a1bSJed Brown /* Initialize PETSc */ 410327415f7SBarry Smith PetscFunctionBeginUser; 4119566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 4129566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 413c4762a1bSJed Brown 414c4762a1bSJed Brown /* Set up the system, we need to create a solver and a mesh and then assign 415c4762a1bSJed Brown * the correct spaces into the mesh*/ 4169566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 4179566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &mesh)); 4189566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, mesh)); 4199566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(mesh, SetupProblem, &user)); 4209566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(mesh, &user, &user, &user)); 4219566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 422c4762a1bSJed Brown 423c4762a1bSJed Brown /* Grab field IS so that we can view the solution by field */ 4249566063dSJacob Faibussowitsch PetscCall(DMCreateFieldIS(mesh, NULL, NULL, &fieldIS)); 425c4762a1bSJed Brown 426c4762a1bSJed Brown /* Create a vector to store the SNES solution, solve the system and grab the 427c4762a1bSJed Brown * solution from SNES */ 4289566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(mesh, &computed)); 4299566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)computed, "computedSolution")); 4309566063dSJacob Faibussowitsch PetscCall(VecSet(computed, 0.0)); 4319566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, computed)); 4329566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &computed)); 4339566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(computed, NULL, "-computedSolution_view")); 434c4762a1bSJed Brown 435c4762a1bSJed Brown /* Now we pull out the portion of the vector corresponding to the 3rd field 436c4762a1bSJed Brown * which is the error between \div{u} computed in a higher dimensional space 437c4762a1bSJed Brown * and p = \div{u} computed in a low dimension space. We report the L2 norm of 438c4762a1bSJed Brown * this vector which should be zero if the H(div) spaces are implemented 439c4762a1bSJed Brown * correctly. */ 4409566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(computed, fieldIS[2], &divErr)); 4419566063dSJacob Faibussowitsch PetscCall(VecNorm(divErr, NORM_2, &divErrNorm)); 4429566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(computed, fieldIS[2], &divErr)); 443c4762a1bSJed Brown exampleSuccess = (PetscBool)(divErrNorm <= errTol); 444c4762a1bSJed Brown 4459566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, stdFormat, divErrNorm, exampleSuccess ? "true" : "false")); 446c4762a1bSJed Brown 447c4762a1bSJed Brown /* Tear down */ 4489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&divErr)); 4499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&computed)); 450*48a46eb9SPierre Jolivet for (i = 0; i < 3; ++i) PetscCall(ISDestroy(&fieldIS[i])); 4519566063dSJacob Faibussowitsch PetscCall(PetscFree(fieldIS)); 4529566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 4539566063dSJacob Faibussowitsch PetscCall(DMDestroy(&mesh)); 4549566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 455b122ec5aSJacob Faibussowitsch return 0; 456c4762a1bSJed Brown } 457c4762a1bSJed Brown 458c4762a1bSJed Brown /*TEST 459c4762a1bSJed Brown testset: 460c4762a1bSJed Brown suffix: 2d_bdm 461c4762a1bSJed Brown requires: triangle 46230602db0SMatthew G. Knepley args: -velocity_petscfe_default_quadrature_order 1 \ 463c4762a1bSJed Brown -velocity_petscspace_degree 1 \ 464c4762a1bSJed Brown -velocity_petscdualspace_type bdm \ 465c4762a1bSJed Brown -divErr_petscspace_degree 1 \ 466c4762a1bSJed Brown -divErr_petscdualspace_lagrange_continuity false \ 467c4762a1bSJed Brown -snes_error_if_not_converged \ 468c4762a1bSJed Brown -ksp_rtol 1e-10 \ 469c4762a1bSJed Brown -ksp_error_if_not_converged \ 470c4762a1bSJed Brown -pc_type fieldsplit\ 471c4762a1bSJed Brown -pc_fieldsplit_detect_saddle_point\ 472c4762a1bSJed Brown -pc_fieldsplit_type schur\ 473c4762a1bSJed Brown -pc_fieldsplit_schur_precondition full 474c4762a1bSJed Brown test: 475c4762a1bSJed Brown suffix: linear 476c4762a1bSJed Brown args: -sol_form linear -mesh_transform none 477c4762a1bSJed Brown test: 478c4762a1bSJed Brown suffix: sinusoidal 479c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform none 480c4762a1bSJed Brown test: 481c4762a1bSJed Brown suffix: sinusoidal_skew 482c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew 483c4762a1bSJed Brown test: 484c4762a1bSJed Brown suffix: sinusoidal_perturb 485c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform perturb 486c4762a1bSJed Brown test: 487c4762a1bSJed Brown suffix: sinusoidal_skew_perturb 488c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew_perturb 489c4762a1bSJed Brown 490c4762a1bSJed Brown testset: 491c4762a1bSJed Brown TODO: broken 492c4762a1bSJed Brown suffix: 2d_bdmq 49330602db0SMatthew G. Knepley args: -dm_plex_simplex false \ 494c4762a1bSJed Brown -velocity_petscspace_degree 1 \ 495c4762a1bSJed Brown -velocity_petscdualspace_type bdm \ 496d092c84bSBrandon Whitchurch -velocity_petscdualspace_lagrange_tensor 1 \ 497c4762a1bSJed Brown -divErr_petscspace_degree 1 \ 498c4762a1bSJed Brown -divErr_petscdualspace_lagrange_continuity false \ 499c4762a1bSJed Brown -snes_error_if_not_converged \ 500c4762a1bSJed Brown -ksp_rtol 1e-10 \ 501c4762a1bSJed Brown -ksp_error_if_not_converged \ 502c4762a1bSJed Brown -pc_type fieldsplit\ 503c4762a1bSJed Brown -pc_fieldsplit_detect_saddle_point\ 504c4762a1bSJed Brown -pc_fieldsplit_type schur\ 505c4762a1bSJed Brown -pc_fieldsplit_schur_precondition full 506c4762a1bSJed Brown test: 507c4762a1bSJed Brown suffix: linear 508c4762a1bSJed Brown args: -sol_form linear -mesh_transform none 509c4762a1bSJed Brown test: 510c4762a1bSJed Brown suffix: sinusoidal 511c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform none 512c4762a1bSJed Brown test: 513c4762a1bSJed Brown suffix: sinusoidal_skew 514c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew 515c4762a1bSJed Brown test: 516c4762a1bSJed Brown suffix: sinusoidal_perturb 517c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform perturb 518c4762a1bSJed Brown test: 519c4762a1bSJed Brown suffix: sinusoidal_skew_perturb 520c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew_perturb 521c4762a1bSJed Brown 522c4762a1bSJed Brown testset: 523c4762a1bSJed Brown suffix: 3d_bdm 524c4762a1bSJed Brown requires: ctetgen 52530602db0SMatthew G. Knepley args: -dm_plex_dim 3 \ 526c4762a1bSJed Brown -velocity_petscspace_degree 1 \ 527c4762a1bSJed Brown -velocity_petscdualspace_type bdm \ 528c4762a1bSJed Brown -divErr_petscspace_degree 1 \ 529c4762a1bSJed Brown -divErr_petscdualspace_lagrange_continuity false \ 530c4762a1bSJed Brown -snes_error_if_not_converged \ 531c4762a1bSJed Brown -ksp_rtol 1e-10 \ 532c4762a1bSJed Brown -ksp_error_if_not_converged \ 533c4762a1bSJed Brown -pc_type fieldsplit \ 534c4762a1bSJed Brown -pc_fieldsplit_detect_saddle_point \ 535c4762a1bSJed Brown -pc_fieldsplit_type schur \ 536c4762a1bSJed Brown -pc_fieldsplit_schur_precondition full 537c4762a1bSJed Brown test: 538c4762a1bSJed Brown suffix: linear 539c4762a1bSJed Brown args: -sol_form linear -mesh_transform none 540c4762a1bSJed Brown test: 541c4762a1bSJed Brown suffix: sinusoidal 542c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform none 543c4762a1bSJed Brown test: 544c4762a1bSJed Brown suffix: sinusoidal_skew 545c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew 546c4762a1bSJed Brown test: 547c4762a1bSJed Brown suffix: sinusoidal_perturb 548c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform perturb 549c4762a1bSJed Brown test: 550c4762a1bSJed Brown suffix: sinusoidal_skew_perturb 551c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew_perturb 552c4762a1bSJed Brown 553c4762a1bSJed Brown testset: 554c4762a1bSJed Brown TODO: broken 555c4762a1bSJed Brown suffix: 3d_bdmq 556c4762a1bSJed Brown requires: ctetgen 55730602db0SMatthew G. Knepley args: -dm_plex_dim 3 \ 55830602db0SMatthew G. Knepley -dm_plex_simplex false \ 559c4762a1bSJed Brown -velocity_petscspace_degree 1 \ 560c4762a1bSJed Brown -velocity_petscdualspace_type bdm \ 561d092c84bSBrandon Whitchurch -velocity_petscdualspace_lagrange_tensor 1 \ 562c4762a1bSJed Brown -divErr_petscspace_degree 1 \ 563c4762a1bSJed Brown -divErr_petscdualspace_lagrange_continuity false \ 564c4762a1bSJed Brown -snes_error_if_not_converged \ 565c4762a1bSJed Brown -ksp_rtol 1e-10 \ 566c4762a1bSJed Brown -ksp_error_if_not_converged \ 567c4762a1bSJed Brown -pc_type fieldsplit \ 568c4762a1bSJed Brown -pc_fieldsplit_detect_saddle_point \ 569c4762a1bSJed Brown -pc_fieldsplit_type schur \ 570c4762a1bSJed Brown -pc_fieldsplit_schur_precondition full 571c4762a1bSJed Brown test: 572c4762a1bSJed Brown suffix: linear 573c4762a1bSJed Brown args: -sol_form linear -mesh_transform none 574c4762a1bSJed Brown test: 575c4762a1bSJed Brown suffix: sinusoidal 576c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform none 577c4762a1bSJed Brown test: 578c4762a1bSJed Brown suffix: sinusoidal_skew 579c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew 580c4762a1bSJed Brown test: 581c4762a1bSJed Brown suffix: sinusoidal_perturb 582c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform perturb 583c4762a1bSJed Brown test: 584c4762a1bSJed Brown suffix: sinusoidal_skew_perturb 585c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew_perturb 586d092c84bSBrandon Whitchurch 587d092c84bSBrandon Whitchurch test: 588d092c84bSBrandon Whitchurch suffix: quad_rt_0 58930602db0SMatthew G. Knepley args: -dm_plex_simplex false -mesh_transform skew \ 590d092c84bSBrandon Whitchurch -divErr_petscspace_degree 1 \ 591d092c84bSBrandon Whitchurch -divErr_petscdualspace_lagrange_continuity false \ 592d092c84bSBrandon Whitchurch -snes_error_if_not_converged \ 593d092c84bSBrandon Whitchurch -ksp_rtol 1e-10 \ 594d092c84bSBrandon Whitchurch -ksp_error_if_not_converged \ 595d092c84bSBrandon Whitchurch -pc_type fieldsplit\ 596d092c84bSBrandon Whitchurch -pc_fieldsplit_detect_saddle_point\ 597d092c84bSBrandon Whitchurch -pc_fieldsplit_type schur\ 598d092c84bSBrandon Whitchurch -pc_fieldsplit_schur_precondition full \ 599d092c84bSBrandon Whitchurch -velocity_petscfe_default_quadrature_order 1 \ 600d092c84bSBrandon Whitchurch -velocity_petscspace_type sum \ 601d092c84bSBrandon Whitchurch -velocity_petscspace_variables 2 \ 602d092c84bSBrandon Whitchurch -velocity_petscspace_components 2 \ 603d092c84bSBrandon Whitchurch -velocity_petscspace_sum_spaces 2 \ 604d092c84bSBrandon Whitchurch -velocity_petscspace_sum_concatenate true \ 605417c287bSToby Isaac -velocity_sumcomp_0_petscspace_variables 2 \ 606417c287bSToby Isaac -velocity_sumcomp_0_petscspace_type tensor \ 607417c287bSToby Isaac -velocity_sumcomp_0_petscspace_tensor_spaces 2 \ 608417c287bSToby Isaac -velocity_sumcomp_0_petscspace_tensor_uniform false \ 609417c287bSToby Isaac -velocity_sumcomp_0_tensorcomp_0_petscspace_degree 1 \ 610417c287bSToby Isaac -velocity_sumcomp_0_tensorcomp_1_petscspace_degree 0 \ 611417c287bSToby Isaac -velocity_sumcomp_1_petscspace_variables 2 \ 612417c287bSToby Isaac -velocity_sumcomp_1_petscspace_type tensor \ 613417c287bSToby Isaac -velocity_sumcomp_1_petscspace_tensor_spaces 2 \ 614417c287bSToby Isaac -velocity_sumcomp_1_petscspace_tensor_uniform false \ 615417c287bSToby Isaac -velocity_sumcomp_1_tensorcomp_0_petscspace_degree 0 \ 616417c287bSToby Isaac -velocity_sumcomp_1_tensorcomp_1_petscspace_degree 1 \ 617d092c84bSBrandon Whitchurch -velocity_petscdualspace_form_degree -1 \ 618d092c84bSBrandon Whitchurch -velocity_petscdualspace_order 1 \ 619d092c84bSBrandon Whitchurch -velocity_petscdualspace_lagrange_trimmed true 620c4762a1bSJed Brown TEST*/ 621