1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 2 #include <petscdmplex.h> 3 4 static PetscErrorCode PetscDualSpaceLagrangeGetTensor_Lagrange(PetscDualSpace sp, PetscBool *tensor) 5 { 6 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 7 8 PetscFunctionBegin; 9 *tensor = lag->tensorSpace; 10 PetscFunctionReturn(0); 11 } 12 13 static PetscErrorCode PetscDualSpaceLagrangeSetTensor_Lagrange(PetscDualSpace sp, PetscBool tensor) 14 { 15 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 16 17 PetscFunctionBegin; 18 lag->tensorSpace = tensor; 19 PetscFunctionReturn(0); 20 } 21 22 /* 23 LatticePointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to at most 'max'. 24 Ordering is lexicographic with lowest index as least significant in ordering. 25 e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,0}. 26 27 Input Parameters: 28 + len - The length of the tuple 29 . max - The maximum sum 30 - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition 31 32 Output Parameter: 33 . tup - A tuple of len integers whos sum is at most 'max' 34 */ 35 static PetscErrorCode LatticePointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[]) 36 { 37 PetscFunctionBegin; 38 while (len--) { 39 max -= tup[len]; 40 if (!max) { 41 tup[len] = 0; 42 break; 43 } 44 } 45 tup[++len]++; 46 PetscFunctionReturn(0); 47 } 48 49 /* 50 TensorPointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that are all less than or equal to 'max'. 51 Ordering is lexicographic with lowest index as least significant in ordering. 52 e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,1}, {0,2}, {1,2}, {2,2}. 53 54 Input Parameters: 55 + len - The length of the tuple 56 . max - The maximum value 57 - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition 58 59 Output Parameter: 60 . tup - A tuple of len integers whos sum is at most 'max' 61 */ 62 static PetscErrorCode TensorPointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[]) 63 { 64 PetscInt i; 65 66 PetscFunctionBegin; 67 for (i = 0; i < len; i++) { 68 if (tup[i] < max) { 69 break; 70 } else { 71 tup[i] = 0; 72 } 73 } 74 tup[i]++; 75 PetscFunctionReturn(0); 76 } 77 78 79 #define BaryIndex(perEdge,a,b,c) (((b)*(2*perEdge+1-(b)))/2)+(c) 80 81 #define CartIndex(perEdge,a,b) (perEdge*(a)+b) 82 83 static PetscErrorCode PetscDualSpaceGetSymmetries_Lagrange(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips) 84 { 85 86 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 87 PetscInt dim, order, p, Nc; 88 PetscErrorCode ierr; 89 90 PetscFunctionBegin; 91 ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr); 92 ierr = PetscDualSpaceGetNumComponents(sp,&Nc);CHKERRQ(ierr); 93 ierr = DMGetDimension(sp->dm,&dim);CHKERRQ(ierr); 94 if (!dim || !lag->continuous || order < 3) PetscFunctionReturn(0); 95 if (dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Lagrange symmetries not implemented for dim = %D > 3",dim); 96 if (!lag->symmetries) { /* store symmetries */ 97 PetscDualSpace hsp; 98 DM K; 99 PetscInt numPoints = 1, d; 100 PetscInt numFaces; 101 PetscInt ***symmetries; 102 const PetscInt ***hsymmetries; 103 104 if (lag->simplexCell) { 105 numFaces = 1 + dim; 106 for (d = 0; d < dim; d++) numPoints = numPoints * 2 + 1; 107 } 108 else { 109 numPoints = PetscPowInt(3,dim); 110 numFaces = 2 * dim; 111 } 112 ierr = PetscCalloc1(numPoints,&symmetries);CHKERRQ(ierr); 113 if (0 < dim && dim < 3) { /* compute self symmetries */ 114 PetscInt **cellSymmetries; 115 116 lag->numSelfSym = 2 * numFaces; 117 lag->selfSymOff = numFaces; 118 ierr = PetscCalloc1(2*numFaces,&cellSymmetries);CHKERRQ(ierr); 119 /* we want to be able to index symmetries directly with the orientations, which range from [-numFaces,numFaces) */ 120 symmetries[0] = &cellSymmetries[numFaces]; 121 if (dim == 1) { 122 PetscInt dofPerEdge = order - 1; 123 124 if (dofPerEdge > 1) { 125 PetscInt i, j, *reverse; 126 127 ierr = PetscMalloc1(dofPerEdge*Nc,&reverse);CHKERRQ(ierr); 128 for (i = 0; i < dofPerEdge; i++) { 129 for (j = 0; j < Nc; j++) { 130 reverse[i*Nc + j] = Nc * (dofPerEdge - 1 - i) + j; 131 } 132 } 133 symmetries[0][-2] = reverse; 134 135 /* yes, this is redundant, but it makes it easier to cleanup if I don't have to worry about what not to free */ 136 ierr = PetscMalloc1(dofPerEdge*Nc,&reverse);CHKERRQ(ierr); 137 for (i = 0; i < dofPerEdge; i++) { 138 for (j = 0; j < Nc; j++) { 139 reverse[i*Nc + j] = Nc * (dofPerEdge - 1 - i) + j; 140 } 141 } 142 symmetries[0][1] = reverse; 143 } 144 } else { 145 PetscInt dofPerEdge = lag->simplexCell ? (order - 2) : (order - 1), s; 146 PetscInt dofPerFace; 147 148 if (dofPerEdge > 1) { 149 for (s = -numFaces; s < numFaces; s++) { 150 PetscInt *sym, i, j, k, l; 151 152 if (!s) continue; 153 if (lag->simplexCell) { 154 dofPerFace = (dofPerEdge * (dofPerEdge + 1))/2; 155 ierr = PetscMalloc1(Nc*dofPerFace,&sym);CHKERRQ(ierr); 156 for (j = 0, l = 0; j < dofPerEdge; j++) { 157 for (k = 0; k < dofPerEdge - j; k++, l++) { 158 i = dofPerEdge - 1 - j - k; 159 switch (s) { 160 case -3: 161 sym[Nc*l] = BaryIndex(dofPerEdge,i,k,j); 162 break; 163 case -2: 164 sym[Nc*l] = BaryIndex(dofPerEdge,j,i,k); 165 break; 166 case -1: 167 sym[Nc*l] = BaryIndex(dofPerEdge,k,j,i); 168 break; 169 case 1: 170 sym[Nc*l] = BaryIndex(dofPerEdge,k,i,j); 171 break; 172 case 2: 173 sym[Nc*l] = BaryIndex(dofPerEdge,j,k,i); 174 break; 175 } 176 } 177 } 178 } else { 179 dofPerFace = dofPerEdge * dofPerEdge; 180 ierr = PetscMalloc1(Nc*dofPerFace,&sym);CHKERRQ(ierr); 181 for (j = 0, l = 0; j < dofPerEdge; j++) { 182 for (k = 0; k < dofPerEdge; k++, l++) { 183 switch (s) { 184 case -4: 185 sym[Nc*l] = CartIndex(dofPerEdge,k,j); 186 break; 187 case -3: 188 sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - j),k); 189 break; 190 case -2: 191 sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - k),(dofPerEdge - 1 - j)); 192 break; 193 case -1: 194 sym[Nc*l] = CartIndex(dofPerEdge,j,(dofPerEdge - 1 - k)); 195 break; 196 case 1: 197 sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - k),j); 198 break; 199 case 2: 200 sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - j),(dofPerEdge - 1 - k)); 201 break; 202 case 3: 203 sym[Nc*l] = CartIndex(dofPerEdge,k,(dofPerEdge - 1 - j)); 204 break; 205 } 206 } 207 } 208 } 209 for (i = 0; i < dofPerFace; i++) { 210 sym[Nc*i] *= Nc; 211 for (j = 1; j < Nc; j++) { 212 sym[Nc*i+j] = sym[Nc*i] + j; 213 } 214 } 215 symmetries[0][s] = sym; 216 } 217 } 218 } 219 } 220 ierr = PetscDualSpaceGetHeightSubspace(sp,1,&hsp);CHKERRQ(ierr); 221 ierr = PetscDualSpaceGetSymmetries(hsp,&hsymmetries,NULL);CHKERRQ(ierr); 222 if (hsymmetries) { 223 PetscBool *seen; 224 const PetscInt *cone; 225 PetscInt KclosureSize, *Kclosure = NULL; 226 227 ierr = PetscDualSpaceGetDM(sp,&K);CHKERRQ(ierr); 228 ierr = PetscCalloc1(numPoints,&seen);CHKERRQ(ierr); 229 ierr = DMPlexGetCone(K,0,&cone);CHKERRQ(ierr); 230 ierr = DMPlexGetTransitiveClosure(K,0,PETSC_TRUE,&KclosureSize,&Kclosure);CHKERRQ(ierr); 231 for (p = 0; p < numFaces; p++) { 232 PetscInt closureSize, *closure = NULL, q; 233 234 ierr = DMPlexGetTransitiveClosure(K,cone[p],PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 235 for (q = 0; q < closureSize; q++) { 236 PetscInt point = closure[2*q], r; 237 238 if(!seen[point]) { 239 for (r = 0; r < KclosureSize; r++) { 240 if (Kclosure[2 * r] == point) break; 241 } 242 seen[point] = PETSC_TRUE; 243 symmetries[r] = (PetscInt **) hsymmetries[q]; 244 } 245 } 246 ierr = DMPlexRestoreTransitiveClosure(K,cone[p],PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 247 } 248 ierr = DMPlexRestoreTransitiveClosure(K,0,PETSC_TRUE,&KclosureSize,&Kclosure);CHKERRQ(ierr); 249 ierr = PetscFree(seen);CHKERRQ(ierr); 250 } 251 lag->symmetries = symmetries; 252 } 253 if (perms) *perms = (const PetscInt ***) lag->symmetries; 254 PetscFunctionReturn(0); 255 } 256 257 /*@C 258 PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis 259 260 Not collective 261 262 Input Parameter: 263 . sp - the PetscDualSpace object 264 265 Output Parameters: 266 + perms - Permutations of the local degrees of freedom, parameterized by the point orientation 267 - flips - Sign reversal of the local degrees of freedom, parameterized by the point orientation 268 269 Note: The permutation and flip arrays are organized in the following way 270 $ perms[p][ornt][dof # on point] = new local dof # 271 $ flips[p][ornt][dof # on point] = reversal or not 272 273 Level: developer 274 275 .seealso: PetscDualSpaceSetSymmetries() 276 @*/ 277 PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips) 278 { 279 PetscErrorCode ierr; 280 281 PetscFunctionBegin; 282 PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1); 283 if (perms) { 284 PetscValidPointer(perms,2); 285 *perms = NULL; 286 } 287 if (flips) { 288 PetscValidPointer(flips,3); 289 *flips = NULL; 290 } 291 if (sp->ops->getsymmetries) { 292 ierr = (sp->ops->getsymmetries)(sp,perms,flips);CHKERRQ(ierr); 293 } 294 PetscFunctionReturn(0); 295 } 296 297 static PetscErrorCode PetscDualSpaceLagrangeView_Ascii(PetscDualSpace sp, PetscViewer viewer) 298 { 299 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 300 PetscErrorCode ierr; 301 302 PetscFunctionBegin; 303 ierr = PetscViewerASCIIPrintf(viewer, "%s %sLagrange dual space\n", lag->continuous ? "Continuous" : "Discontinuous", lag->tensorSpace ? "tensor " : "");CHKERRQ(ierr); 304 PetscFunctionReturn(0); 305 } 306 307 PetscErrorCode PetscDualSpaceView_Lagrange(PetscDualSpace sp, PetscViewer viewer) 308 { 309 PetscBool iascii; 310 PetscErrorCode ierr; 311 312 PetscFunctionBegin; 313 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 314 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 315 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 316 if (iascii) {ierr = PetscDualSpaceLagrangeView_Ascii(sp, viewer);CHKERRQ(ierr);} 317 PetscFunctionReturn(0); 318 } 319 320 static PetscErrorCode PetscDualSpaceGetDimension_SingleCell_Lagrange(PetscDualSpace sp, PetscInt order, PetscInt *dim) 321 { 322 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 323 PetscReal D = 1.0; 324 PetscInt n, i; 325 PetscErrorCode ierr; 326 327 PetscFunctionBegin; 328 *dim = -1; /* Ensure that the compiler knows *dim is set. */ 329 ierr = DMGetDimension(sp->dm, &n);CHKERRQ(ierr); 330 if (!lag->tensorSpace) { 331 for (i = 1; i <= n; ++i) { 332 D *= ((PetscReal) (order+i))/i; 333 } 334 *dim = (PetscInt) (D + 0.5); 335 } else { 336 *dim = 1; 337 for (i = 0; i < n; ++i) *dim *= (order+1); 338 } 339 *dim *= sp->Nc; 340 PetscFunctionReturn(0); 341 } 342 343 static PetscErrorCode PetscDualSpaceCreateHeightSubspace_Lagrange(PetscDualSpace sp, PetscInt height, PetscDualSpace *bdsp) 344 { 345 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 346 PetscBool continuous, tensor; 347 PetscInt order; 348 PetscErrorCode ierr; 349 350 PetscFunctionBegin; 351 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 352 PetscValidPointer(bdsp,2); 353 ierr = PetscDualSpaceLagrangeGetContinuity(sp,&continuous);CHKERRQ(ierr); 354 ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr); 355 if (height == 0) { 356 ierr = PetscObjectReference((PetscObject)sp);CHKERRQ(ierr); 357 *bdsp = sp; 358 } else if (continuous == PETSC_FALSE || !order) { 359 *bdsp = NULL; 360 } else { 361 DM dm, K; 362 PetscInt dim; 363 364 ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 365 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 366 if (height > dim || height < 0) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Asked for dual space at height %d for dimension %d reference element\n",height,dim);} 367 ierr = PetscDualSpaceDuplicate(sp,bdsp);CHKERRQ(ierr); 368 ierr = PetscDualSpaceCreateReferenceCell(*bdsp, dim-height, lag->simplexCell, &K);CHKERRQ(ierr); 369 ierr = PetscDualSpaceSetDM(*bdsp, K);CHKERRQ(ierr); 370 ierr = DMDestroy(&K);CHKERRQ(ierr); 371 ierr = PetscDualSpaceLagrangeGetTensor(sp,&tensor);CHKERRQ(ierr); 372 ierr = PetscDualSpaceLagrangeSetTensor(*bdsp,tensor);CHKERRQ(ierr); 373 ierr = PetscDualSpaceSetUp(*bdsp);CHKERRQ(ierr); 374 } 375 PetscFunctionReturn(0); 376 } 377 378 PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp) 379 { 380 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 381 DM dm = sp->dm; 382 PetscInt order = sp->order; 383 PetscInt Nc = sp->Nc; 384 PetscBool continuous; 385 PetscSection csection; 386 Vec coordinates; 387 PetscReal *qpoints, *qweights; 388 PetscInt depth, dim, pdimMax, pStart, pEnd, p, *pStratStart, *pStratEnd, coneSize, d, f = 0, c; 389 PetscBool simplex, tensorSpace; 390 PetscErrorCode ierr; 391 392 PetscFunctionBegin; 393 /* Classify element type */ 394 if (!order) lag->continuous = PETSC_FALSE; 395 continuous = lag->continuous; 396 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 397 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 398 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 399 ierr = PetscCalloc1(dim+1, &lag->numDof);CHKERRQ(ierr); 400 ierr = PetscMalloc2(depth+1,&pStratStart,depth+1,&pStratEnd);CHKERRQ(ierr); 401 for (d = 0; d <= depth; ++d) {ierr = DMPlexGetDepthStratum(dm, d, &pStratStart[d], &pStratEnd[d]);CHKERRQ(ierr);} 402 ierr = DMPlexGetConeSize(dm, pStratStart[depth], &coneSize);CHKERRQ(ierr); 403 ierr = DMGetCoordinateSection(dm, &csection);CHKERRQ(ierr); 404 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 405 if (depth == 1) { 406 if (coneSize == dim+1) simplex = PETSC_TRUE; 407 else if (coneSize == 1 << dim) simplex = PETSC_FALSE; 408 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support simplices and tensor product cells"); 409 } else if (depth == dim) { 410 if (coneSize == dim+1) simplex = PETSC_TRUE; 411 else if (coneSize == 2 * dim) simplex = PETSC_FALSE; 412 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support simplices and tensor product cells"); 413 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support cell-vertex meshes or interpolated meshes"); 414 lag->simplexCell = simplex; 415 if (dim > 1 && continuous && lag->simplexCell == lag->tensorSpace) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Mismatching simplex/tensor cells and spaces only allowed for discontinuous elements"); 416 tensorSpace = lag->tensorSpace; 417 lag->height = 0; 418 lag->subspaces = NULL; 419 if (continuous && sp->order > 0 && dim > 0) { 420 PetscInt i; 421 422 lag->height = dim; 423 ierr = PetscMalloc1(dim,&lag->subspaces);CHKERRQ(ierr); 424 ierr = PetscDualSpaceCreateHeightSubspace_Lagrange(sp,1,&lag->subspaces[0]);CHKERRQ(ierr); 425 ierr = PetscDualSpaceSetUp(lag->subspaces[0]);CHKERRQ(ierr); 426 for (i = 1; i < dim; i++) { 427 ierr = PetscDualSpaceGetHeightSubspace(lag->subspaces[i-1],1,&lag->subspaces[i]);CHKERRQ(ierr); 428 ierr = PetscObjectReference((PetscObject)(lag->subspaces[i]));CHKERRQ(ierr); 429 } 430 } 431 ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, sp->order, &pdimMax);CHKERRQ(ierr); 432 pdimMax *= (pStratEnd[depth] - pStratStart[depth]); 433 ierr = PetscMalloc1(pdimMax, &sp->functional);CHKERRQ(ierr); 434 if (!dim) { 435 for (c = 0; c < Nc; ++c) { 436 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 437 ierr = PetscCalloc1(Nc, &qweights);CHKERRQ(ierr); 438 ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr); 439 ierr = PetscQuadratureSetData(sp->functional[f], 0, Nc, 1, NULL, qweights);CHKERRQ(ierr); 440 qweights[c] = 1.0; 441 ++f; 442 lag->numDof[0]++; 443 } 444 } else { 445 PetscInt *tup; 446 PetscReal *v0, *hv0, *J, *invJ, detJ, hdetJ; 447 PetscSection section; 448 449 ierr = PetscSectionCreate(PETSC_COMM_SELF,§ion);CHKERRQ(ierr); 450 ierr = PetscSectionSetChart(section,pStart,pEnd);CHKERRQ(ierr); 451 ierr = PetscCalloc5(dim+1,&tup,dim,&v0,dim,&hv0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 452 for (p = pStart; p < pEnd; p++) { 453 PetscInt pointDim, d, nFunc = 0; 454 PetscDualSpace hsp; 455 456 ierr = DMPlexComputeCellGeometryFEM(dm, p, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 457 for (d = 0; d < depth; d++) {if (p >= pStratStart[d] && p < pStratEnd[d]) break;} 458 pointDim = (depth == 1 && d == 1) ? dim : d; 459 hsp = ((pointDim < dim) && lag->subspaces) ? lag->subspaces[dim - pointDim - 1] : NULL; 460 if (hsp) { 461 PetscDualSpace_Lag *hlag = (PetscDualSpace_Lag *) hsp->data; 462 DM hdm; 463 464 ierr = PetscDualSpaceGetDM(hsp,&hdm);CHKERRQ(ierr); 465 ierr = DMPlexComputeCellGeometryFEM(hdm, 0, NULL, hv0, NULL, NULL, &hdetJ);CHKERRQ(ierr); 466 nFunc = lag->numDof[pointDim] = hlag->numDof[pointDim]; 467 } 468 if (pointDim == dim) { 469 /* Cells, create for self */ 470 PetscInt orderEff = continuous ? (!tensorSpace ? order-1-dim : order-2) : order; 471 PetscReal denom = continuous ? order : (!tensorSpace ? order+1+dim : order+2); 472 PetscReal numer = (!simplex || !tensorSpace) ? 2. : (2./dim); 473 PetscReal dx = numer/denom; 474 PetscInt cdim, d, d2; 475 476 if (orderEff < 0) continue; 477 ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, orderEff, &cdim);CHKERRQ(ierr); 478 ierr = PetscMemzero(tup,(dim+1)*sizeof(PetscInt));CHKERRQ(ierr); 479 if (!tensorSpace) { 480 while (!tup[dim]) { 481 for (c = 0; c < Nc; ++c) { 482 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 483 ierr = PetscMalloc1(dim, &qpoints);CHKERRQ(ierr); 484 ierr = PetscCalloc1(Nc, &qweights);CHKERRQ(ierr); 485 ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr); 486 ierr = PetscQuadratureSetData(sp->functional[f], dim, Nc, 1, qpoints, qweights);CHKERRQ(ierr); 487 for (d = 0; d < dim; ++d) { 488 qpoints[d] = v0[d]; 489 for (d2 = 0; d2 < dim; ++d2) qpoints[d] += J[d*dim+d2]*((tup[d2]+1)*dx); 490 } 491 qweights[c] = 1.0; 492 ++f; 493 } 494 ierr = LatticePointLexicographic_Internal(dim, orderEff, tup);CHKERRQ(ierr); 495 } 496 } else { 497 while (!tup[dim]) { 498 for (c = 0; c < Nc; ++c) { 499 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 500 ierr = PetscMalloc1(dim, &qpoints);CHKERRQ(ierr); 501 ierr = PetscCalloc1(Nc, &qweights);CHKERRQ(ierr); 502 ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr); 503 ierr = PetscQuadratureSetData(sp->functional[f], dim, Nc, 1, qpoints, qweights);CHKERRQ(ierr); 504 for (d = 0; d < dim; ++d) { 505 qpoints[d] = v0[d]; 506 for (d2 = 0; d2 < dim; ++d2) qpoints[d] += J[d*dim+d2]*((tup[d2]+1)*dx); 507 } 508 qweights[c] = 1.0; 509 ++f; 510 } 511 ierr = TensorPointLexicographic_Internal(dim, orderEff, tup);CHKERRQ(ierr); 512 } 513 } 514 lag->numDof[dim] = cdim; 515 } else { /* transform functionals from subspaces */ 516 PetscInt q; 517 518 for (q = 0; q < nFunc; q++, f++) { 519 PetscQuadrature fn; 520 PetscInt fdim, Nc, c, nPoints, i; 521 const PetscReal *points; 522 const PetscReal *weights; 523 PetscReal *qpoints; 524 PetscReal *qweights; 525 526 ierr = PetscDualSpaceGetFunctional(hsp, q, &fn);CHKERRQ(ierr); 527 ierr = PetscQuadratureGetData(fn,&fdim,&Nc,&nPoints,&points,&weights);CHKERRQ(ierr); 528 if (fdim != pointDim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected height dual space dim %D, got %D",pointDim,fdim); 529 ierr = PetscMalloc1(nPoints * dim, &qpoints);CHKERRQ(ierr); 530 ierr = PetscCalloc1(nPoints * Nc, &qweights);CHKERRQ(ierr); 531 for (i = 0; i < nPoints; i++) { 532 PetscInt j, k; 533 PetscReal *qp = &qpoints[i * dim]; 534 535 for (c = 0; c < Nc; ++c) qweights[i*Nc+c] = weights[i*Nc+c]; 536 for (j = 0; j < dim; ++j) qp[j] = v0[j]; 537 for (j = 0; j < dim; ++j) { 538 for (k = 0; k < pointDim; k++) qp[j] += J[dim * j + k] * (points[pointDim * i + k] - hv0[k]); 539 } 540 } 541 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 542 ierr = PetscQuadratureSetOrder(sp->functional[f],0);CHKERRQ(ierr); 543 ierr = PetscQuadratureSetData(sp->functional[f],dim,Nc,nPoints,qpoints,qweights);CHKERRQ(ierr); 544 } 545 } 546 ierr = PetscSectionSetDof(section,p,lag->numDof[pointDim]);CHKERRQ(ierr); 547 } 548 ierr = PetscFree5(tup,v0,hv0,J,invJ);CHKERRQ(ierr); 549 ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 550 { /* reorder to closure order */ 551 PetscInt *key, count; 552 PetscQuadrature *reorder = NULL; 553 554 ierr = PetscCalloc1(f,&key);CHKERRQ(ierr); 555 ierr = PetscMalloc1(f*sp->Nc,&reorder);CHKERRQ(ierr); 556 557 for (p = pStratStart[depth], count = 0; p < pStratEnd[depth]; p++) { 558 PetscInt *closure = NULL, closureSize, c; 559 560 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 561 for (c = 0; c < closureSize; c++) { 562 PetscInt point = closure[2 * c], dof, off, i; 563 564 ierr = PetscSectionGetDof(section,point,&dof);CHKERRQ(ierr); 565 ierr = PetscSectionGetOffset(section,point,&off);CHKERRQ(ierr); 566 for (i = 0; i < dof; i++) { 567 PetscInt fi = i + off; 568 if (!key[fi]) { 569 key[fi] = 1; 570 reorder[count++] = sp->functional[fi]; 571 } 572 } 573 } 574 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 575 } 576 ierr = PetscFree(sp->functional);CHKERRQ(ierr); 577 sp->functional = reorder; 578 ierr = PetscFree(key);CHKERRQ(ierr); 579 } 580 ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 581 } 582 if (pStratEnd[depth] == 1 && f != pdimMax) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of dual basis vectors %d not equal to dimension %d", f, pdimMax); 583 ierr = PetscFree2(pStratStart, pStratEnd);CHKERRQ(ierr); 584 if (f > pdimMax) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of dual basis vectors %d is greater than dimension %d", f, pdimMax); 585 PetscFunctionReturn(0); 586 } 587 588 PetscErrorCode PetscDualSpaceDestroy_Lagrange(PetscDualSpace sp) 589 { 590 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 591 PetscInt i; 592 PetscErrorCode ierr; 593 594 PetscFunctionBegin; 595 if (lag->symmetries) { 596 PetscInt **selfSyms = lag->symmetries[0]; 597 598 if (selfSyms) { 599 PetscInt i, **allocated = &selfSyms[-lag->selfSymOff]; 600 601 for (i = 0; i < lag->numSelfSym; i++) { 602 ierr = PetscFree(allocated[i]);CHKERRQ(ierr); 603 } 604 ierr = PetscFree(allocated);CHKERRQ(ierr); 605 } 606 ierr = PetscFree(lag->symmetries);CHKERRQ(ierr); 607 } 608 for (i = 0; i < lag->height; i++) { 609 ierr = PetscDualSpaceDestroy(&lag->subspaces[i]);CHKERRQ(ierr); 610 } 611 ierr = PetscFree(lag->subspaces);CHKERRQ(ierr); 612 ierr = PetscFree(lag->numDof);CHKERRQ(ierr); 613 ierr = PetscFree(lag);CHKERRQ(ierr); 614 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", NULL);CHKERRQ(ierr); 615 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", NULL);CHKERRQ(ierr); 616 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", NULL);CHKERRQ(ierr); 617 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", NULL);CHKERRQ(ierr); 618 PetscFunctionReturn(0); 619 } 620 621 PetscErrorCode PetscDualSpaceDuplicate_Lagrange(PetscDualSpace sp, PetscDualSpace *spNew) 622 { 623 PetscInt order, Nc; 624 PetscBool cont, tensor; 625 PetscErrorCode ierr; 626 627 PetscFunctionBegin; 628 ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) sp), spNew);CHKERRQ(ierr); 629 ierr = PetscDualSpaceSetType(*spNew, PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 630 ierr = PetscDualSpaceGetOrder(sp, &order);CHKERRQ(ierr); 631 ierr = PetscDualSpaceSetOrder(*spNew, order);CHKERRQ(ierr); 632 ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr); 633 ierr = PetscDualSpaceSetNumComponents(*spNew, Nc);CHKERRQ(ierr); 634 ierr = PetscDualSpaceLagrangeGetContinuity(sp, &cont);CHKERRQ(ierr); 635 ierr = PetscDualSpaceLagrangeSetContinuity(*spNew, cont);CHKERRQ(ierr); 636 ierr = PetscDualSpaceLagrangeGetTensor(sp, &tensor);CHKERRQ(ierr); 637 ierr = PetscDualSpaceLagrangeSetTensor(*spNew, tensor);CHKERRQ(ierr); 638 PetscFunctionReturn(0); 639 } 640 641 PetscErrorCode PetscDualSpaceSetFromOptions_Lagrange(PetscOptionItems *PetscOptionsObject,PetscDualSpace sp) 642 { 643 PetscBool continuous, tensor, flg; 644 PetscErrorCode ierr; 645 646 PetscFunctionBegin; 647 ierr = PetscDualSpaceLagrangeGetContinuity(sp, &continuous);CHKERRQ(ierr); 648 ierr = PetscDualSpaceLagrangeGetTensor(sp, &tensor);CHKERRQ(ierr); 649 ierr = PetscOptionsHead(PetscOptionsObject,"PetscDualSpace Lagrange Options");CHKERRQ(ierr); 650 ierr = PetscOptionsBool("-petscdualspace_lagrange_continuity", "Flag for continuous element", "PetscDualSpaceLagrangeSetContinuity", continuous, &continuous, &flg);CHKERRQ(ierr); 651 if (flg) {ierr = PetscDualSpaceLagrangeSetContinuity(sp, continuous);CHKERRQ(ierr);} 652 ierr = PetscOptionsBool("-petscdualspace_lagrange_tensor", "Flag for tensor dual space", "PetscDualSpaceLagrangeSetContinuity", tensor, &tensor, &flg);CHKERRQ(ierr); 653 if (flg) {ierr = PetscDualSpaceLagrangeSetTensor(sp, tensor);CHKERRQ(ierr);} 654 ierr = PetscOptionsTail();CHKERRQ(ierr); 655 PetscFunctionReturn(0); 656 } 657 658 PetscErrorCode PetscDualSpaceGetDimension_Lagrange(PetscDualSpace sp, PetscInt *dim) 659 { 660 DM K; 661 const PetscInt *numDof; 662 PetscInt spatialDim, Nc, size = 0, d; 663 PetscErrorCode ierr; 664 665 PetscFunctionBegin; 666 ierr = PetscDualSpaceGetDM(sp, &K);CHKERRQ(ierr); 667 ierr = PetscDualSpaceGetNumDof(sp, &numDof);CHKERRQ(ierr); 668 ierr = DMGetDimension(K, &spatialDim);CHKERRQ(ierr); 669 ierr = DMPlexGetHeightStratum(K, 0, NULL, &Nc);CHKERRQ(ierr); 670 if (Nc == 1) {ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, sp->order, dim);CHKERRQ(ierr); PetscFunctionReturn(0);} 671 for (d = 0; d <= spatialDim; ++d) { 672 PetscInt pStart, pEnd; 673 674 ierr = DMPlexGetDepthStratum(K, d, &pStart, &pEnd);CHKERRQ(ierr); 675 size += (pEnd-pStart)*numDof[d]; 676 } 677 *dim = size; 678 PetscFunctionReturn(0); 679 } 680 681 PetscErrorCode PetscDualSpaceGetNumDof_Lagrange(PetscDualSpace sp, const PetscInt **numDof) 682 { 683 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 684 685 PetscFunctionBegin; 686 *numDof = lag->numDof; 687 PetscFunctionReturn(0); 688 } 689 690 static PetscErrorCode PetscDualSpaceLagrangeGetContinuity_Lagrange(PetscDualSpace sp, PetscBool *continuous) 691 { 692 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 693 694 PetscFunctionBegin; 695 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 696 PetscValidPointer(continuous, 2); 697 *continuous = lag->continuous; 698 PetscFunctionReturn(0); 699 } 700 701 static PetscErrorCode PetscDualSpaceLagrangeSetContinuity_Lagrange(PetscDualSpace sp, PetscBool continuous) 702 { 703 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 704 705 PetscFunctionBegin; 706 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 707 lag->continuous = continuous; 708 PetscFunctionReturn(0); 709 } 710 711 /*@ 712 PetscDualSpaceLagrangeGetContinuity - Retrieves the flag for element continuity 713 714 Not Collective 715 716 Input Parameter: 717 . sp - the PetscDualSpace 718 719 Output Parameter: 720 . continuous - flag for element continuity 721 722 Level: intermediate 723 724 .keywords: PetscDualSpace, Lagrange, continuous, discontinuous 725 .seealso: PetscDualSpaceLagrangeSetContinuity() 726 @*/ 727 PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace sp, PetscBool *continuous) 728 { 729 PetscErrorCode ierr; 730 731 PetscFunctionBegin; 732 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 733 PetscValidPointer(continuous, 2); 734 ierr = PetscTryMethod(sp, "PetscDualSpaceLagrangeGetContinuity_C", (PetscDualSpace,PetscBool*),(sp,continuous));CHKERRQ(ierr); 735 PetscFunctionReturn(0); 736 } 737 738 /*@ 739 PetscDualSpaceLagrangeSetContinuity - Indicate whether the element is continuous 740 741 Logically Collective on PetscDualSpace 742 743 Input Parameters: 744 + sp - the PetscDualSpace 745 - continuous - flag for element continuity 746 747 Options Database: 748 . -petscdualspace_lagrange_continuity <bool> 749 750 Level: intermediate 751 752 .keywords: PetscDualSpace, Lagrange, continuous, discontinuous 753 .seealso: PetscDualSpaceLagrangeGetContinuity() 754 @*/ 755 PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace sp, PetscBool continuous) 756 { 757 PetscErrorCode ierr; 758 759 PetscFunctionBegin; 760 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 761 PetscValidLogicalCollectiveBool(sp, continuous, 2); 762 ierr = PetscTryMethod(sp, "PetscDualSpaceLagrangeSetContinuity_C", (PetscDualSpace,PetscBool),(sp,continuous));CHKERRQ(ierr); 763 PetscFunctionReturn(0); 764 } 765 766 PetscErrorCode PetscDualSpaceGetHeightSubspace_Lagrange(PetscDualSpace sp, PetscInt height, PetscDualSpace *bdsp) 767 { 768 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 769 PetscErrorCode ierr; 770 771 PetscFunctionBegin; 772 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 773 PetscValidPointer(bdsp,2); 774 if (height == 0) { 775 *bdsp = sp; 776 } 777 else { 778 DM dm; 779 PetscInt dim; 780 781 ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 782 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 783 if (height > dim || height < 0) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Asked for dual space at height %d for dimension %d reference element\n",height,dim);} 784 if (height <= lag->height) { 785 *bdsp = lag->subspaces[height-1]; 786 } 787 else { 788 *bdsp = NULL; 789 } 790 } 791 PetscFunctionReturn(0); 792 } 793 794 PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp) 795 { 796 PetscFunctionBegin; 797 sp->ops->setfromoptions = PetscDualSpaceSetFromOptions_Lagrange; 798 sp->ops->setup = PetscDualSpaceSetUp_Lagrange; 799 sp->ops->view = PetscDualSpaceView_Lagrange; 800 sp->ops->destroy = PetscDualSpaceDestroy_Lagrange; 801 sp->ops->duplicate = PetscDualSpaceDuplicate_Lagrange; 802 sp->ops->getdimension = PetscDualSpaceGetDimension_Lagrange; 803 sp->ops->getnumdof = PetscDualSpaceGetNumDof_Lagrange; 804 sp->ops->getheightsubspace = PetscDualSpaceGetHeightSubspace_Lagrange; 805 sp->ops->getsymmetries = PetscDualSpaceGetSymmetries_Lagrange; 806 sp->ops->apply = PetscDualSpaceApplyDefault; 807 sp->ops->applyall = PetscDualSpaceApplyAllDefault; 808 sp->ops->createallpoints = PetscDualSpaceCreateAllPointsDefault; 809 PetscFunctionReturn(0); 810 } 811 812 /*MC 813 PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals 814 815 Level: intermediate 816 817 .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType() 818 M*/ 819 820 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp) 821 { 822 PetscDualSpace_Lag *lag; 823 PetscErrorCode ierr; 824 825 PetscFunctionBegin; 826 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 827 ierr = PetscNewLog(sp,&lag);CHKERRQ(ierr); 828 sp->data = lag; 829 830 lag->numDof = NULL; 831 lag->simplexCell = PETSC_TRUE; 832 lag->tensorSpace = PETSC_FALSE; 833 lag->continuous = PETSC_TRUE; 834 835 ierr = PetscDualSpaceInitialize_Lagrange(sp);CHKERRQ(ierr); 836 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", PetscDualSpaceLagrangeGetContinuity_Lagrange);CHKERRQ(ierr); 837 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", PetscDualSpaceLagrangeSetContinuity_Lagrange);CHKERRQ(ierr); 838 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", PetscDualSpaceLagrangeGetTensor_Lagrange);CHKERRQ(ierr); 839 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", PetscDualSpaceLagrangeSetTensor_Lagrange);CHKERRQ(ierr); 840 PetscFunctionReturn(0); 841 } 842 843