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