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