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 /* 11*0e3d61c9SBarry Smith We are using the system 12*0e3d61c9SBarry Smith 13*0e3d61c9SBarry Smith \vec{u} = \vec{\hat{u}} 14*0e3d61c9SBarry Smith p = \div{\vec{u}} in low degree approximation space 15*0e3d61c9SBarry Smith d = \div{\vec{u}} - p == 0 in higher degree approximation space 16*0e3d61c9SBarry Smith 17*0e3d61c9SBarry Smith That is, we are using the field d to compute the error between \div{\vec{u}} 18*0e3d61c9SBarry Smith computed in a space 1 degree higher than p and the value of p which is 19*0e3d61c9SBarry Smith \div{u} computed in the low degree space. If H-div 20*0e3d61c9SBarry Smith elements are implemented correctly then this should be identically zero since 21*0e3d61c9SBarry Smith the divergence of a function in H(div) should be exactly representable in L_2 22*0e3d61c9SBarry Smith by definition. 23c4762a1bSJed Brown */ 24c4762a1bSJed Brown static PetscErrorCode zero_func(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar *u,void *ctx) 25c4762a1bSJed Brown { 26c4762a1bSJed Brown PetscInt c; 27c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = 0; 28c4762a1bSJed Brown return 0; 29c4762a1bSJed Brown } 30c4762a1bSJed Brown /* Linear Exact Functions 31c4762a1bSJed Brown \vec{u} = \vec{x}; 32c4762a1bSJed Brown p = dim; 33c4762a1bSJed Brown */ 34c4762a1bSJed Brown static PetscErrorCode linear_u(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar *u,void *ctx) 35c4762a1bSJed Brown { 36c4762a1bSJed Brown PetscInt c; 37c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = x[c]; 38c4762a1bSJed Brown return 0; 39c4762a1bSJed Brown } 40c4762a1bSJed Brown static PetscErrorCode linear_p(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar *u,void *ctx) 41c4762a1bSJed Brown { 42c4762a1bSJed Brown u[0] = dim; 43c4762a1bSJed Brown return 0; 44c4762a1bSJed Brown } 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* Sinusoidal Exact Functions 47c4762a1bSJed Brown * u_i = \sin{2*\pi*x_i} 48c4762a1bSJed Brown * p = \Sum_{i=1}^{dim} 2*\pi*cos{2*\pi*x_i} 49c4762a1bSJed Brown * */ 50c4762a1bSJed Brown 51c4762a1bSJed Brown static PetscErrorCode sinusoid_u(PetscInt dim,PetscReal time,const PetscReal 52c4762a1bSJed Brown x[],PetscInt Nc,PetscScalar *u,void *ctx) 53c4762a1bSJed Brown { 54c4762a1bSJed Brown PetscInt c; 55c4762a1bSJed Brown for (c = 0; c< Nc; ++c) u[c] = PetscSinReal(2*PETSC_PI*x[c]); 56c4762a1bSJed Brown return 0; 57c4762a1bSJed Brown } 58c4762a1bSJed Brown static PetscErrorCode sinusoid_p(PetscInt dim,PetscReal time,const PetscReal 59c4762a1bSJed Brown x[],PetscInt Nc,PetscScalar *u,void *ctx) 60c4762a1bSJed Brown { 61c4762a1bSJed Brown PetscInt d; 62c4762a1bSJed Brown u[0]=0; 63c4762a1bSJed Brown for (d=0; d<dim; ++d) u[0] += 2*PETSC_PI*PetscCosReal(2*PETSC_PI*x[d]); 64c4762a1bSJed Brown return 0; 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* Pointwise residual for u = u*. Need one of these for each possible u* */ 68c4762a1bSJed Brown 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[]) 69c4762a1bSJed Brown { 70c4762a1bSJed Brown PetscInt i; 71c4762a1bSJed Brown PetscScalar *u_rhs; 72c4762a1bSJed Brown 73c4762a1bSJed Brown PetscCalloc1(dim,&u_rhs); 74c4762a1bSJed Brown (void) linear_u(dim,t,x,dim,u_rhs,NULL); 75c4762a1bSJed Brown for (i = 0; i < dim; ++i) f0[i] = u[uOff[0]+i]-u_rhs[i]; 76c4762a1bSJed Brown PetscFree(u_rhs); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79c4762a1bSJed Brown 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[]) 80c4762a1bSJed Brown { 81c4762a1bSJed Brown PetscInt i; 82c4762a1bSJed Brown PetscScalar *u_rhs; 83c4762a1bSJed Brown 84c4762a1bSJed Brown PetscCalloc1(dim,&u_rhs); 85c4762a1bSJed Brown (void) sinusoid_u(dim,t,x,dim,u_rhs,NULL); 86c4762a1bSJed Brown for (i = 0; i < dim; ++i) f0[i] = u[uOff[0]+i]-u_rhs[i]; 87c4762a1bSJed Brown PetscFree(u_rhs); 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Residual function for enforcing p = \div{u}. */ 91c4762a1bSJed Brown 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[]) 92c4762a1bSJed Brown { 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[1]] - divu; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* Residual function for p_err = \div{u} - p. */ 102c4762a1bSJed Brown 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[]) 103c4762a1bSJed Brown { 104c4762a1bSJed Brown PetscInt i; 105c4762a1bSJed Brown PetscScalar divu; 106c4762a1bSJed Brown 107c4762a1bSJed Brown divu = 0.; 108c4762a1bSJed Brown for (i = 0; i < dim; ++i) divu += u_x[uOff_x[0] + i*dim +i]; 109c4762a1bSJed Brown f0[0] = u[uOff[2]] - u[uOff[1]] + divu; 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* Boundary residual for the embedding system. Need one for each form of 113c4762a1bSJed Brown * solution. These enforce u = \hat{u} at the boundary. */ 114c4762a1bSJed Brown 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[]) 115c4762a1bSJed Brown { 116c4762a1bSJed Brown PetscInt d; 117c4762a1bSJed Brown PetscScalar *u_rhs; 118c4762a1bSJed Brown PetscCalloc1(dim,&u_rhs); 119c4762a1bSJed Brown (void) sinusoid_u(dim,t,x,dim,u_rhs,NULL); 120c4762a1bSJed Brown 121c4762a1bSJed Brown for (d=0; d<dim; ++d) f0[d] = u_rhs[d]; 122c4762a1bSJed Brown PetscFree(u_rhs); 123c4762a1bSJed Brown 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown 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[]) 127c4762a1bSJed Brown { 128c4762a1bSJed Brown PetscInt d; 129c4762a1bSJed Brown PetscScalar *u_rhs; 130c4762a1bSJed Brown PetscCalloc1(dim,&u_rhs); 131c4762a1bSJed Brown (void) linear_u(dim,t,x,dim,u_rhs,NULL); 132c4762a1bSJed Brown 133c4762a1bSJed Brown for (d=0; d<dim; ++d) f0[d] = u_rhs[d]; 134c4762a1bSJed Brown PetscFree(u_rhs); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown /* Jacobian functions. For the following, v is the test function associated with 137c4762a1bSJed Brown * u, q the test function associated with p, and w the test function associated 138c4762a1bSJed Brown * with d. */ 139c4762a1bSJed Brown /* <v, u> */ 140c4762a1bSJed Brown 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[]) 141c4762a1bSJed Brown { 142c4762a1bSJed Brown PetscInt c; 143c4762a1bSJed Brown for (c = 0; c < dim; ++c) g0[c * dim + c] = 1.0; 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* <q, p> */ 147c4762a1bSJed Brown 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[]) 148c4762a1bSJed Brown { 149c4762a1bSJed Brown PetscInt d; 150c4762a1bSJed Brown for (d=0; d< dim; ++d) g0[d * dim + d] = 1.0; 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* -<q, \div{u}> For the embedded system. This is different from the method of 154c4762a1bSJed Brown * manufactured solution because instead of computing <q,\div{u}> - <q,f> we 155c4762a1bSJed Brown * need <q,p> - <q,\div{u}.*/ 156c4762a1bSJed Brown 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[]) 157c4762a1bSJed Brown { 158c4762a1bSJed Brown PetscInt d; 159c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d * dim + d] = -1.0; 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* <w, p> This is only used by the embedded system. Where we need to compute 163c4762a1bSJed Brown * <w,d> - <w,p> + <w, \div{u}>*/ 164c4762a1bSJed Brown 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[]) 165c4762a1bSJed Brown { 166c4762a1bSJed Brown PetscInt d; 167c4762a1bSJed Brown for (d=0; d< dim; ++d) g0[d * dim + d] = -1.0; 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* <w, d> */ 171c4762a1bSJed Brown 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[]) 172c4762a1bSJed Brown { 173c4762a1bSJed Brown PetscInt c; 174c4762a1bSJed Brown for (c = 0; c < dim; ++c) g0[c*dim+c] = 1.0; 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* <w, \div{u}> for the embedded system. */ 178c4762a1bSJed Brown 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[]) 179c4762a1bSJed Brown { 180c4762a1bSJed Brown PetscInt d; 181c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* Enum and string array for selecting mesh perturbation options */ 185c4762a1bSJed Brown typedef enum { NONE = 0,PERTURB = 1,SKEW = 2,SKEW_PERTURB = 3 } Transform; 186c4762a1bSJed Brown const char* const TransformTypes[] = {"none","perturb","skew","skew_perturb","Perturbation","",NULL}; 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* Enum and string array for selecting the form of the exact solution*/ 189c4762a1bSJed Brown typedef enum 190c4762a1bSJed Brown { LINEAR = 0,SINUSOIDAL = 1 } Solution; 191c4762a1bSJed Brown const char* const SolutionTypes[] = {"linear","sinusoidal","Solution","",NULL}; 192c4762a1bSJed Brown 193c4762a1bSJed Brown typedef struct 194c4762a1bSJed Brown { 195c4762a1bSJed Brown Transform mesh_transform; 196c4762a1bSJed Brown Solution sol_form; 197c4762a1bSJed Brown } UserCtx; 198c4762a1bSJed Brown 199c4762a1bSJed Brown /* Process command line options and initialize the UserCtx struct */ 200c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm,UserCtx * user) 201c4762a1bSJed Brown { 202c4762a1bSJed Brown PetscErrorCode ierr; 203c4762a1bSJed Brown 204c4762a1bSJed Brown PetscFunctionBegin; 205c4762a1bSJed Brown /* Default to 2D, unperturbed triangle mesh and Linear solution.*/ 206c4762a1bSJed Brown user->mesh_transform = NONE; 207c4762a1bSJed Brown user->sol_form = LINEAR; 208c4762a1bSJed Brown 209c4762a1bSJed Brown ierr = PetscOptionsBegin(comm,"","H-div Test Options","DMPLEX");CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 212c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 213c4762a1bSJed Brown PetscFunctionReturn(0); 214c4762a1bSJed Brown } 215c4762a1bSJed Brown 216c4762a1bSJed Brown /* Perturb the position of each mesh vertex by a small amount.*/ 217c4762a1bSJed Brown static PetscErrorCode PerturbMesh(DM *mesh,PetscScalar *coordVals,PetscInt npoints,PetscInt dim) 218c4762a1bSJed Brown { 219c4762a1bSJed Brown PetscInt i,j,k; 220c4762a1bSJed Brown PetscErrorCode ierr; 221d092c84bSBrandon Whitchurch PetscReal minCoords[3],maxCoords[3],maxPert[3],randVal,amp; 222c4762a1bSJed Brown PetscRandom ran; 223c4762a1bSJed Brown 224c4762a1bSJed Brown PetscFunctionBegin; 225c4762a1bSJed Brown ierr = DMGetCoordinateDim(*mesh,&dim);CHKERRQ(ierr); 226c4762a1bSJed Brown ierr = DMGetLocalBoundingBox(*mesh,minCoords,maxCoords);CHKERRQ(ierr); 227c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&ran);CHKERRQ(ierr); 228c4762a1bSJed Brown 229c4762a1bSJed Brown /* Compute something approximately equal to half an edge length. This is the 230c4762a1bSJed Brown * most we can perturb points and gaurantee that there won't be any topology 231c4762a1bSJed Brown * issues. */ 232d092c84bSBrandon Whitchurch for (k = 0; k < dim; ++k) maxPert[k] = 0.025 * (maxCoords[k] - minCoords[k]) / (PetscPowReal(npoints,1. / dim) - 1); 233c4762a1bSJed Brown /* For each mesh vertex */ 234c4762a1bSJed Brown for (i = 0; i < npoints; ++i) { 235c4762a1bSJed Brown /* For each coordinate of the vertex */ 236c4762a1bSJed Brown for (j = 0; j < dim; ++j) { 237c4762a1bSJed Brown /* Generate a random amplitude in [-0.5*maxPert, 0.5*maxPert] */ 238c4762a1bSJed Brown ierr = PetscRandomGetValueReal(ran,&randVal);CHKERRQ(ierr); 239c4762a1bSJed Brown amp = maxPert[j] * (randVal - 0.5); 240c4762a1bSJed Brown /* Add the perturbation to the vertex*/ 241d092c84bSBrandon Whitchurch coordVals[dim * i + j] += amp; 242c4762a1bSJed Brown } 243c4762a1bSJed Brown } 244c4762a1bSJed Brown 245c4762a1bSJed Brown PetscRandomDestroy(&ran); 246c4762a1bSJed Brown PetscFunctionReturn(0); 247c4762a1bSJed Brown } 248c4762a1bSJed Brown 249c4762a1bSJed Brown /* Apply a global skew transformation to the mesh. */ 250c4762a1bSJed Brown static PetscErrorCode SkewMesh(DM * mesh,PetscScalar * coordVals,PetscInt npoints,PetscInt dim) 251c4762a1bSJed Brown { 252c4762a1bSJed Brown PetscInt i,j,k,l; 253c4762a1bSJed Brown PetscErrorCode ierr; 254c4762a1bSJed Brown PetscScalar * transMat; 255c4762a1bSJed Brown PetscScalar tmpcoord[3]; 256c4762a1bSJed Brown PetscRandom ran; 257c4762a1bSJed Brown PetscReal randVal; 258c4762a1bSJed Brown 259c4762a1bSJed Brown PetscFunctionBegin; 260c4762a1bSJed Brown ierr = PetscCalloc1(dim * dim,&transMat);CHKERRQ(ierr); 261c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&ran);CHKERRQ(ierr); 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* Make a matrix representing a skew transformation */ 264c4762a1bSJed Brown for (i = 0; i < dim; ++i) { 265c4762a1bSJed Brown for (j = 0; j < dim; ++j) { 266c4762a1bSJed Brown ierr = PetscRandomGetValueReal(ran,&randVal);CHKERRQ(ierr); 267d092c84bSBrandon Whitchurch if (i == j) transMat[i * dim + j] = 1.; 268c4762a1bSJed Brown else if (j < i) transMat[i * dim + j] = 2 * (j + i)*randVal; 269c4762a1bSJed Brown else transMat[i * dim + j] = 0; 270c4762a1bSJed Brown } 271c4762a1bSJed Brown } 272c4762a1bSJed Brown 273c4762a1bSJed Brown /* Multiply each coordinate vector by our tranformation.*/ 274c4762a1bSJed Brown for (i = 0; i < npoints; ++i) { 275c4762a1bSJed Brown for (j = 0; j < dim; ++j) { 276c4762a1bSJed Brown tmpcoord[j] = 0; 277c4762a1bSJed Brown for (k = 0; k < dim; ++k) tmpcoord[j] += coordVals[dim * i + k] * transMat[dim * k + j]; 278c4762a1bSJed Brown } 279c4762a1bSJed Brown for (l = 0; l < dim; ++l) coordVals[dim * i + l] = tmpcoord[l]; 280c4762a1bSJed Brown } 281c4762a1bSJed Brown ierr = PetscFree(transMat);CHKERRQ(ierr); 282c4762a1bSJed Brown ierr = PetscRandomDestroy(&ran);CHKERRQ(ierr); 283c4762a1bSJed Brown PetscFunctionReturn(0); 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* Accesses the mesh coordinate array and performs the transformation operations 287c4762a1bSJed Brown * specified by the user options */ 288c4762a1bSJed Brown static PetscErrorCode TransformMesh(UserCtx * user,DM * mesh) 289c4762a1bSJed Brown { 290c4762a1bSJed Brown PetscErrorCode ierr; 291c4762a1bSJed Brown PetscInt dim,npoints; 292c4762a1bSJed Brown PetscScalar * coordVals; 293c4762a1bSJed Brown Vec coords; 294c4762a1bSJed Brown 295c4762a1bSJed Brown PetscFunctionBegin; 296c4762a1bSJed Brown ierr = DMGetCoordinates(*mesh,&coords);CHKERRQ(ierr); 297c4762a1bSJed Brown ierr = VecGetArray(coords,&coordVals);CHKERRQ(ierr); 298c4762a1bSJed Brown ierr = VecGetLocalSize(coords,&npoints);CHKERRQ(ierr); 299c4762a1bSJed Brown ierr = DMGetCoordinateDim(*mesh,&dim);CHKERRQ(ierr); 300c4762a1bSJed Brown npoints = npoints / dim; 301c4762a1bSJed Brown 302c4762a1bSJed Brown switch (user->mesh_transform) { 303c4762a1bSJed Brown case PERTURB: 304c4762a1bSJed Brown ierr = PerturbMesh(mesh,coordVals,npoints,dim);CHKERRQ(ierr); 305c4762a1bSJed Brown break; 306c4762a1bSJed Brown case SKEW: 307c4762a1bSJed Brown ierr = SkewMesh(mesh,coordVals,npoints,dim);CHKERRQ(ierr); 308c4762a1bSJed Brown break; 309c4762a1bSJed Brown case SKEW_PERTURB: 310c4762a1bSJed Brown ierr = SkewMesh(mesh,coordVals,npoints,dim);CHKERRQ(ierr); 311c4762a1bSJed Brown ierr = PerturbMesh(mesh,coordVals,npoints,dim);CHKERRQ(ierr); 312c4762a1bSJed Brown break; 313c4762a1bSJed Brown default: 314c4762a1bSJed Brown PetscFunctionReturn(-1); 315c4762a1bSJed Brown } 316c4762a1bSJed Brown ierr = VecRestoreArray(coords,&coordVals);CHKERRQ(ierr); 317c4762a1bSJed Brown ierr = DMSetCoordinates(*mesh,coords);CHKERRQ(ierr); 318c4762a1bSJed Brown PetscFunctionReturn(0); 319c4762a1bSJed Brown } 320c4762a1bSJed Brown 321c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm,UserCtx * user,DM * mesh) 322c4762a1bSJed Brown { 323c4762a1bSJed Brown PetscErrorCode ierr; 324c4762a1bSJed Brown 325c4762a1bSJed Brown PetscFunctionBegin; 32630602db0SMatthew G. Knepley ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 32730602db0SMatthew G. Knepley ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 32830602db0SMatthew G. Knepley ierr = DMSetFromOptions(*mesh);CHKERRQ(ierr); 329c4762a1bSJed Brown 330c4762a1bSJed Brown /* Perform any mesh transformations if specified by user */ 331c4762a1bSJed Brown if (user->mesh_transform != NONE) { 332c4762a1bSJed Brown ierr = TransformMesh(user,mesh);CHKERRQ(ierr); 333c4762a1bSJed Brown } 334c4762a1bSJed Brown 335c4762a1bSJed Brown /* Get any other mesh options from the command line */ 336c4762a1bSJed Brown ierr = DMSetApplicationContext(*mesh,user);CHKERRQ(ierr); 337c4762a1bSJed Brown ierr = DMViewFromOptions(*mesh,NULL,"-dm_view");CHKERRQ(ierr); 338c4762a1bSJed Brown PetscFunctionReturn(0); 339c4762a1bSJed Brown } 340c4762a1bSJed Brown 341c4762a1bSJed Brown /* Setup the system of equations that we wish to solve */ 342c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm,UserCtx * user) 343c4762a1bSJed Brown { 344c4762a1bSJed Brown PetscDS prob; 34545480ffeSMatthew G. Knepley DMLabel label; 346c4762a1bSJed Brown PetscErrorCode ierr; 347c4762a1bSJed Brown const PetscInt id=1; 348c4762a1bSJed Brown 349c4762a1bSJed Brown PetscFunctionBegin; 350c4762a1bSJed Brown ierr = DMGetDS(dm,&prob);CHKERRQ(ierr); 351c4762a1bSJed Brown /* All of these are independent of the user's choice of solution */ 352c4762a1bSJed Brown ierr = PetscDSSetResidual(prob,1,f0_q,NULL);CHKERRQ(ierr); 353c4762a1bSJed Brown ierr = PetscDSSetResidual(prob,2,f0_w,NULL);CHKERRQ(ierr); 354c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob,0,0,g0_vu,NULL,NULL,NULL);CHKERRQ(ierr); 355c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob,1,0,NULL,g1_qu,NULL,NULL);CHKERRQ(ierr); 356c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob,1,1,g0_qp,NULL,NULL,NULL);CHKERRQ(ierr); 357c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob,2,0,NULL,g1_wu,NULL,NULL);CHKERRQ(ierr); 358c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob,2,1,g0_wp,NULL,NULL,NULL);CHKERRQ(ierr); 359c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob,2,2,g0_wd,NULL,NULL,NULL);CHKERRQ(ierr); 360c4762a1bSJed Brown 361c4762a1bSJed Brown /* Field 2 is the error between \div{u} and pressure in a higher dimensional 362c4762a1bSJed Brown * space. If all is right this should be machine zero. */ 363c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob,2,zero_func,NULL);CHKERRQ(ierr); 364c4762a1bSJed Brown 365c4762a1bSJed Brown switch (user->sol_form) { 366c4762a1bSJed Brown case LINEAR: 367c4762a1bSJed Brown ierr = PetscDSSetResidual(prob,0,f0_v_linear,NULL);CHKERRQ(ierr); 368c4762a1bSJed Brown ierr = PetscDSSetBdResidual(prob,0,f0_bd_u_linear,NULL);CHKERRQ(ierr); 369c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob,0,linear_u,NULL);CHKERRQ(ierr); 370c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob,1,linear_p,NULL);CHKERRQ(ierr); 371c4762a1bSJed Brown break; 372c4762a1bSJed Brown case SINUSOIDAL: 373c4762a1bSJed Brown ierr = PetscDSSetResidual(prob,0,f0_v_sinusoid,NULL);CHKERRQ(ierr); 374c4762a1bSJed Brown ierr = PetscDSSetBdResidual(prob,0,f0_bd_u_sinusoid,NULL);CHKERRQ(ierr); 375c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob,0,sinusoid_u,NULL);CHKERRQ(ierr); 376c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob,1,sinusoid_p,NULL);CHKERRQ(ierr); 377c4762a1bSJed Brown break; 378c4762a1bSJed Brown default: 379c4762a1bSJed Brown PetscFunctionReturn(-1); 380c4762a1bSJed Brown } 381c4762a1bSJed Brown 38245480ffeSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 38345480ffeSMatthew G. Knepley ierr = PetscDSAddBoundary(prob,DM_BC_NATURAL,"Boundary Integral",label,1,&id,0,0,NULL,(void (*)(void))NULL,NULL,user,NULL);CHKERRQ(ierr); 384c4762a1bSJed Brown PetscFunctionReturn(0); 385c4762a1bSJed Brown } 386c4762a1bSJed Brown 387c4762a1bSJed Brown /* Create the finite element spaces we will use for this system */ 388c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM mesh,PetscErrorCode (*setup)(DM,UserCtx*),UserCtx *user) 389c4762a1bSJed Brown { 390c4762a1bSJed Brown DM cdm = mesh; 391c4762a1bSJed Brown PetscFE fevel,fepres,fedivErr; 39230602db0SMatthew G. Knepley PetscInt dim; 39330602db0SMatthew G. Knepley PetscBool simplex; 394c4762a1bSJed Brown PetscErrorCode ierr; 395c4762a1bSJed Brown 396c4762a1bSJed Brown PetscFunctionBegin; 39730602db0SMatthew G. Knepley ierr = DMGetDimension(mesh, &dim);CHKERRQ(ierr); 39830602db0SMatthew G. Knepley ierr = DMPlexIsSimplex(mesh, &simplex);CHKERRQ(ierr); 399c4762a1bSJed Brown /* Create FE objects and give them names so that options can be set from 400c4762a1bSJed Brown * command line */ 40130602db0SMatthew G. Knepley ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) mesh),dim,dim,simplex,"velocity_",-1,&fevel);CHKERRQ(ierr); 402c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fevel,"velocity");CHKERRQ(ierr); 403c4762a1bSJed Brown 40430602db0SMatthew G. Knepley ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) mesh),dim,1,simplex,"pressure_",-1,&fepres);CHKERRQ(ierr); 405c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fepres,"pressure");CHKERRQ(ierr); 406c4762a1bSJed Brown 407c4762a1bSJed Brown ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) 40830602db0SMatthew G. Knepley mesh),dim,1,simplex,"divErr_",-1,&fedivErr);CHKERRQ(ierr); 409c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fedivErr,"divErr");CHKERRQ(ierr); 410c4762a1bSJed Brown 411c4762a1bSJed Brown ierr = PetscFECopyQuadrature(fevel,fepres);CHKERRQ(ierr); 412c4762a1bSJed Brown ierr = PetscFECopyQuadrature(fevel,fedivErr);CHKERRQ(ierr); 413c4762a1bSJed Brown 414c4762a1bSJed Brown /* Associate the FE objects with the mesh and setup the system */ 415c4762a1bSJed Brown ierr = DMSetField(mesh,0,NULL,(PetscObject) fevel);CHKERRQ(ierr); 416c4762a1bSJed Brown ierr = DMSetField(mesh,1,NULL,(PetscObject) fepres);CHKERRQ(ierr); 417c4762a1bSJed Brown ierr = DMSetField(mesh,2,NULL,(PetscObject) fedivErr);CHKERRQ(ierr); 418c4762a1bSJed Brown ierr = DMCreateDS(mesh);CHKERRQ(ierr); 419c4762a1bSJed Brown ierr = (*setup)(mesh,user);CHKERRQ(ierr); 420c4762a1bSJed Brown 421c4762a1bSJed Brown while (cdm) { 422c4762a1bSJed Brown ierr = DMCopyDisc(mesh,cdm);CHKERRQ(ierr); 423c4762a1bSJed Brown ierr = DMGetCoarseDM(cdm,&cdm);CHKERRQ(ierr); 424c4762a1bSJed Brown } 425c4762a1bSJed Brown 426c4762a1bSJed Brown /* The Mesh now owns the fields, so we can destroy the FEs created in this 427c4762a1bSJed Brown * function */ 428c4762a1bSJed Brown ierr = PetscFEDestroy(&fevel);CHKERRQ(ierr); 429c4762a1bSJed Brown ierr = PetscFEDestroy(&fepres);CHKERRQ(ierr); 430c4762a1bSJed Brown ierr = PetscFEDestroy(&fedivErr);CHKERRQ(ierr); 431c4762a1bSJed Brown ierr = DMDestroy(&cdm);CHKERRQ(ierr); 432c4762a1bSJed Brown PetscFunctionReturn(0); 433c4762a1bSJed Brown } 434c4762a1bSJed Brown 435c4762a1bSJed Brown int main(int argc,char **argv) 436c4762a1bSJed Brown { 437c4762a1bSJed Brown PetscInt i; 438c4762a1bSJed Brown UserCtx user; 439c4762a1bSJed Brown DM mesh; 440c4762a1bSJed Brown SNES snes; 441c4762a1bSJed Brown Vec computed,divErr; 442c4762a1bSJed Brown PetscReal divErrNorm; 443c4762a1bSJed Brown PetscErrorCode ierr; 444c4762a1bSJed Brown IS * fieldIS; 445c4762a1bSJed Brown PetscBool exampleSuccess = PETSC_FALSE; 446d092c84bSBrandon Whitchurch const PetscReal errTol = 10. * PETSC_SMALL; 447c4762a1bSJed Brown 448c4762a1bSJed Brown char stdFormat[] = "L2 Norm of the Divergence Error is: %g\n H(div) elements working correctly: %s\n"; 449c4762a1bSJed Brown 450c4762a1bSJed Brown /* Initialize PETSc */ 451c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 452c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD,&user);CHKERRQ(ierr); 453c4762a1bSJed Brown 454c4762a1bSJed Brown /* Set up the system, we need to create a solver and a mesh and then assign 455c4762a1bSJed Brown * the correct spaces into the mesh*/ 456c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 457c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD,&user,&mesh);CHKERRQ(ierr); 458c4762a1bSJed Brown ierr = SNESSetDM(snes,mesh);CHKERRQ(ierr); 459c4762a1bSJed Brown ierr = SetupDiscretization(mesh,SetupProblem,&user);CHKERRQ(ierr); 460c4762a1bSJed Brown ierr = DMPlexSetSNESLocalFEM(mesh,&user,&user,&user);CHKERRQ(ierr); 461c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 462c4762a1bSJed Brown 463c4762a1bSJed Brown /* Grab field IS so that we can view the solution by field */ 464c4762a1bSJed Brown ierr = DMCreateFieldIS(mesh,NULL,NULL,&fieldIS);CHKERRQ(ierr); 465c4762a1bSJed Brown 466c4762a1bSJed Brown /* Create a vector to store the SNES solution, solve the system and grab the 467c4762a1bSJed Brown * solution from SNES */ 468c4762a1bSJed Brown ierr = DMCreateGlobalVector(mesh,&computed);CHKERRQ(ierr); 469c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) computed,"computedSolution");CHKERRQ(ierr); 470c4762a1bSJed Brown ierr = VecSet(computed,0.0);CHKERRQ(ierr); 471c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,computed);CHKERRQ(ierr); 472c4762a1bSJed Brown ierr = SNESGetSolution(snes,&computed);CHKERRQ(ierr); 473c4762a1bSJed Brown ierr = VecViewFromOptions(computed,NULL,"-computedSolution_view");CHKERRQ(ierr); 474c4762a1bSJed Brown 475c4762a1bSJed Brown /* Now we pull out the portion of the vector corresponding to the 3rd field 476c4762a1bSJed Brown * which is the error between \div{u} computed in a higher dimensional space 477c4762a1bSJed Brown * and p = \div{u} computed in a low dimension space. We report the L2 norm of 478c4762a1bSJed Brown * this vector which should be zero if the H(div) spaces are implemented 479c4762a1bSJed Brown * correctly. */ 480c4762a1bSJed Brown ierr = VecGetSubVector(computed,fieldIS[2],&divErr);CHKERRQ(ierr); 481c4762a1bSJed Brown ierr = VecNorm(divErr,NORM_2,&divErrNorm);CHKERRQ(ierr); 482c4762a1bSJed Brown ierr = VecRestoreSubVector(computed,fieldIS[2],&divErr);CHKERRQ(ierr); 483c4762a1bSJed Brown exampleSuccess = (PetscBool)(divErrNorm <= errTol); 484c4762a1bSJed Brown 485c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,stdFormat,divErrNorm,exampleSuccess ? "true" : "false");CHKERRQ(ierr); 486c4762a1bSJed Brown 487c4762a1bSJed Brown /* Tear down */ 488c4762a1bSJed Brown ierr = VecDestroy(&divErr);CHKERRQ(ierr); 489c4762a1bSJed Brown ierr = VecDestroy(&computed);CHKERRQ(ierr); 490c4762a1bSJed Brown for (i = 0; i < 3; ++i) { 491c4762a1bSJed Brown ierr = ISDestroy(&fieldIS[i]);CHKERRQ(ierr); 492c4762a1bSJed Brown } 493c4762a1bSJed Brown ierr = PetscFree(fieldIS);CHKERRQ(ierr); 494c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 495c4762a1bSJed Brown ierr = DMDestroy(&mesh);CHKERRQ(ierr); 496c4762a1bSJed Brown ierr = PetscFinalize(); 497c4762a1bSJed Brown return ierr; 498c4762a1bSJed Brown } 499c4762a1bSJed Brown 500c4762a1bSJed Brown /*TEST 501c4762a1bSJed Brown testset: 502c4762a1bSJed Brown suffix: 2d_bdm 503c4762a1bSJed Brown requires: triangle 50430602db0SMatthew G. Knepley args: -velocity_petscfe_default_quadrature_order 1 \ 505c4762a1bSJed Brown -velocity_petscspace_degree 1 \ 506c4762a1bSJed Brown -velocity_petscdualspace_type bdm \ 507c4762a1bSJed Brown -divErr_petscspace_degree 1 \ 508c4762a1bSJed Brown -divErr_petscdualspace_lagrange_continuity false \ 509c4762a1bSJed Brown -snes_error_if_not_converged \ 510c4762a1bSJed Brown -ksp_rtol 1e-10 \ 511c4762a1bSJed Brown -ksp_error_if_not_converged \ 512c4762a1bSJed Brown -pc_type fieldsplit\ 513c4762a1bSJed Brown -pc_fieldsplit_detect_saddle_point\ 514c4762a1bSJed Brown -pc_fieldsplit_type schur\ 515c4762a1bSJed Brown -pc_fieldsplit_schur_precondition full 516c4762a1bSJed Brown test: 517c4762a1bSJed Brown suffix: linear 518c4762a1bSJed Brown args: -sol_form linear -mesh_transform none 519c4762a1bSJed Brown test: 520c4762a1bSJed Brown suffix: sinusoidal 521c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform none 522c4762a1bSJed Brown test: 523c4762a1bSJed Brown suffix: sinusoidal_skew 524c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew 525c4762a1bSJed Brown test: 526c4762a1bSJed Brown suffix: sinusoidal_perturb 527c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform perturb 528c4762a1bSJed Brown test: 529c4762a1bSJed Brown suffix: sinusoidal_skew_perturb 530c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew_perturb 531c4762a1bSJed Brown 532c4762a1bSJed Brown testset: 533c4762a1bSJed Brown TODO: broken 534c4762a1bSJed Brown suffix: 2d_bdmq 53530602db0SMatthew G. Knepley args: -dm_plex_simplex false \ 536c4762a1bSJed Brown -velocity_petscspace_degree 1 \ 537c4762a1bSJed Brown -velocity_petscdualspace_type bdm \ 538d092c84bSBrandon Whitchurch -velocity_petscdualspace_lagrange_tensor 1 \ 539c4762a1bSJed Brown -divErr_petscspace_degree 1 \ 540c4762a1bSJed Brown -divErr_petscdualspace_lagrange_continuity false \ 541c4762a1bSJed Brown -snes_error_if_not_converged \ 542c4762a1bSJed Brown -ksp_rtol 1e-10 \ 543c4762a1bSJed Brown -ksp_error_if_not_converged \ 544c4762a1bSJed Brown -pc_type fieldsplit\ 545c4762a1bSJed Brown -pc_fieldsplit_detect_saddle_point\ 546c4762a1bSJed Brown -pc_fieldsplit_type schur\ 547c4762a1bSJed Brown -pc_fieldsplit_schur_precondition full 548c4762a1bSJed Brown test: 549c4762a1bSJed Brown suffix: linear 550c4762a1bSJed Brown args: -sol_form linear -mesh_transform none 551c4762a1bSJed Brown test: 552c4762a1bSJed Brown suffix: sinusoidal 553c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform none 554c4762a1bSJed Brown test: 555c4762a1bSJed Brown suffix: sinusoidal_skew 556c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew 557c4762a1bSJed Brown test: 558c4762a1bSJed Brown suffix: sinusoidal_perturb 559c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform perturb 560c4762a1bSJed Brown test: 561c4762a1bSJed Brown suffix: sinusoidal_skew_perturb 562c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew_perturb 563c4762a1bSJed Brown 564c4762a1bSJed Brown testset: 565c4762a1bSJed Brown suffix: 3d_bdm 566c4762a1bSJed Brown requires: ctetgen 56730602db0SMatthew G. Knepley args: -dm_plex_dim 3 \ 568c4762a1bSJed Brown -velocity_petscspace_degree 1 \ 569c4762a1bSJed Brown -velocity_petscdualspace_type bdm \ 570c4762a1bSJed Brown -divErr_petscspace_degree 1 \ 571c4762a1bSJed Brown -divErr_petscdualspace_lagrange_continuity false \ 572c4762a1bSJed Brown -snes_error_if_not_converged \ 573c4762a1bSJed Brown -ksp_rtol 1e-10 \ 574c4762a1bSJed Brown -ksp_error_if_not_converged \ 575c4762a1bSJed Brown -pc_type fieldsplit \ 576c4762a1bSJed Brown -pc_fieldsplit_detect_saddle_point \ 577c4762a1bSJed Brown -pc_fieldsplit_type schur \ 578c4762a1bSJed Brown -pc_fieldsplit_schur_precondition full 579c4762a1bSJed Brown test: 580c4762a1bSJed Brown suffix: linear 581c4762a1bSJed Brown args: -sol_form linear -mesh_transform none 582c4762a1bSJed Brown test: 583c4762a1bSJed Brown suffix: sinusoidal 584c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform none 585c4762a1bSJed Brown test: 586c4762a1bSJed Brown suffix: sinusoidal_skew 587c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew 588c4762a1bSJed Brown test: 589c4762a1bSJed Brown suffix: sinusoidal_perturb 590c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform perturb 591c4762a1bSJed Brown test: 592c4762a1bSJed Brown suffix: sinusoidal_skew_perturb 593c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew_perturb 594c4762a1bSJed Brown 595c4762a1bSJed Brown testset: 596c4762a1bSJed Brown TODO: broken 597c4762a1bSJed Brown suffix: 3d_bdmq 598c4762a1bSJed Brown requires: ctetgen 59930602db0SMatthew G. Knepley args: -dm_plex_dim 3 \ 60030602db0SMatthew G. Knepley -dm_plex_simplex false \ 601c4762a1bSJed Brown -velocity_petscspace_degree 1 \ 602c4762a1bSJed Brown -velocity_petscdualspace_type bdm \ 603d092c84bSBrandon Whitchurch -velocity_petscdualspace_lagrange_tensor 1 \ 604c4762a1bSJed Brown -divErr_petscspace_degree 1 \ 605c4762a1bSJed Brown -divErr_petscdualspace_lagrange_continuity false \ 606c4762a1bSJed Brown -snes_error_if_not_converged \ 607c4762a1bSJed Brown -ksp_rtol 1e-10 \ 608c4762a1bSJed Brown -ksp_error_if_not_converged \ 609c4762a1bSJed Brown -pc_type fieldsplit \ 610c4762a1bSJed Brown -pc_fieldsplit_detect_saddle_point \ 611c4762a1bSJed Brown -pc_fieldsplit_type schur \ 612c4762a1bSJed Brown -pc_fieldsplit_schur_precondition full 613c4762a1bSJed Brown test: 614c4762a1bSJed Brown suffix: linear 615c4762a1bSJed Brown args: -sol_form linear -mesh_transform none 616c4762a1bSJed Brown test: 617c4762a1bSJed Brown suffix: sinusoidal 618c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform none 619c4762a1bSJed Brown test: 620c4762a1bSJed Brown suffix: sinusoidal_skew 621c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew 622c4762a1bSJed Brown test: 623c4762a1bSJed Brown suffix: sinusoidal_perturb 624c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform perturb 625c4762a1bSJed Brown test: 626c4762a1bSJed Brown suffix: sinusoidal_skew_perturb 627c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew_perturb 628d092c84bSBrandon Whitchurch 629d092c84bSBrandon Whitchurch test: 630d092c84bSBrandon Whitchurch suffix: quad_rt_0 63130602db0SMatthew G. Knepley args: -dm_plex_simplex false -mesh_transform skew \ 632d092c84bSBrandon Whitchurch -divErr_petscspace_degree 1 \ 633d092c84bSBrandon Whitchurch -divErr_petscdualspace_lagrange_continuity false \ 634d092c84bSBrandon Whitchurch -snes_error_if_not_converged \ 635d092c84bSBrandon Whitchurch -ksp_rtol 1e-10 \ 636d092c84bSBrandon Whitchurch -ksp_error_if_not_converged \ 637d092c84bSBrandon Whitchurch -pc_type fieldsplit\ 638d092c84bSBrandon Whitchurch -pc_fieldsplit_detect_saddle_point\ 639d092c84bSBrandon Whitchurch -pc_fieldsplit_type schur\ 640d092c84bSBrandon Whitchurch -pc_fieldsplit_schur_precondition full \ 641d092c84bSBrandon Whitchurch -velocity_petscfe_default_quadrature_order 1 \ 642d092c84bSBrandon Whitchurch -velocity_petscspace_type sum \ 643d092c84bSBrandon Whitchurch -velocity_petscspace_variables 2 \ 644d092c84bSBrandon Whitchurch -velocity_petscspace_components 2 \ 645d092c84bSBrandon Whitchurch -velocity_petscspace_sum_spaces 2 \ 646d092c84bSBrandon Whitchurch -velocity_petscspace_sum_concatenate true \ 647d092c84bSBrandon Whitchurch -velocity_subspace0_petscspace_variables 2 \ 648d092c84bSBrandon Whitchurch -velocity_subspace0_petscspace_type tensor \ 649d092c84bSBrandon Whitchurch -velocity_subspace0_petscspace_tensor_spaces 2 \ 650d092c84bSBrandon Whitchurch -velocity_subspace0_petscspace_tensor_uniform false \ 651d092c84bSBrandon Whitchurch -velocity_subspace0_subspace_0_petscspace_degree 1 \ 652d092c84bSBrandon Whitchurch -velocity_subspace0_subspace_1_petscspace_degree 0 \ 653d092c84bSBrandon Whitchurch -velocity_subspace1_petscspace_variables 2 \ 654d092c84bSBrandon Whitchurch -velocity_subspace1_petscspace_type tensor \ 655d092c84bSBrandon Whitchurch -velocity_subspace1_petscspace_tensor_spaces 2 \ 656d092c84bSBrandon Whitchurch -velocity_subspace1_petscspace_tensor_uniform false \ 657d092c84bSBrandon Whitchurch -velocity_subspace1_subspace_0_petscspace_degree 0 \ 658d092c84bSBrandon Whitchurch -velocity_subspace1_subspace_1_petscspace_degree 1 \ 659d092c84bSBrandon Whitchurch -velocity_petscdualspace_form_degree -1 \ 660d092c84bSBrandon Whitchurch -velocity_petscdualspace_order 1 \ 661d092c84bSBrandon Whitchurch -velocity_petscdualspace_lagrange_trimmed true 662c4762a1bSJed Brown TEST*/ 663