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 of order %D", lag->continuous ? "Continuous" : "Discontinuous", lag->tensorSpace ? "Tensor " : "", sp->order, sp->Nc);CHKERRQ(ierr); 304 if (sp->Nc > 1) {ierr = PetscViewerASCIIPrintf(viewer, " with %D components", sp->Nc);CHKERRQ(ierr);} 305 ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 306 PetscFunctionReturn(0); 307 } 308 309 PetscErrorCode PetscDualSpaceView_Lagrange(PetscDualSpace sp, PetscViewer viewer) 310 { 311 PetscBool iascii; 312 PetscErrorCode ierr; 313 314 PetscFunctionBegin; 315 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 316 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 317 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 318 if (iascii) {ierr = PetscDualSpaceLagrangeView_Ascii(sp, viewer);CHKERRQ(ierr);} 319 PetscFunctionReturn(0); 320 } 321 322 static PetscErrorCode PetscDualSpaceGetDimension_SingleCell_Lagrange(PetscDualSpace sp, PetscInt order, PetscInt *dim) 323 { 324 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 325 PetscReal D = 1.0; 326 PetscInt n, i; 327 PetscErrorCode ierr; 328 329 PetscFunctionBegin; 330 *dim = -1; /* Ensure that the compiler knows *dim is set. */ 331 ierr = DMGetDimension(sp->dm, &n);CHKERRQ(ierr); 332 if (!lag->tensorSpace) { 333 for (i = 1; i <= n; ++i) { 334 D *= ((PetscReal) (order+i))/i; 335 } 336 *dim = (PetscInt) (D + 0.5); 337 } else { 338 *dim = 1; 339 for (i = 0; i < n; ++i) *dim *= (order+1); 340 } 341 *dim *= sp->Nc; 342 PetscFunctionReturn(0); 343 } 344 345 static PetscErrorCode PetscDualSpaceCreateHeightSubspace_Lagrange(PetscDualSpace sp, PetscInt height, PetscDualSpace *bdsp) 346 { 347 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 348 PetscBool continuous, tensor; 349 PetscInt order; 350 PetscErrorCode ierr; 351 352 PetscFunctionBegin; 353 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 354 PetscValidPointer(bdsp,2); 355 ierr = PetscDualSpaceLagrangeGetContinuity(sp,&continuous);CHKERRQ(ierr); 356 ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr); 357 if (height == 0) { 358 ierr = PetscObjectReference((PetscObject)sp);CHKERRQ(ierr); 359 *bdsp = sp; 360 } else if (continuous == PETSC_FALSE || !order) { 361 *bdsp = NULL; 362 } else { 363 DM dm, K; 364 PetscInt dim; 365 366 ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 367 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 368 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);} 369 ierr = PetscDualSpaceDuplicate(sp,bdsp);CHKERRQ(ierr); 370 ierr = PetscDualSpaceCreateReferenceCell(*bdsp, dim-height, lag->simplexCell, &K);CHKERRQ(ierr); 371 ierr = PetscDualSpaceSetDM(*bdsp, K);CHKERRQ(ierr); 372 ierr = DMDestroy(&K);CHKERRQ(ierr); 373 ierr = PetscDualSpaceLagrangeGetTensor(sp,&tensor);CHKERRQ(ierr); 374 ierr = PetscDualSpaceLagrangeSetTensor(*bdsp,tensor);CHKERRQ(ierr); 375 ierr = PetscDualSpaceSetUp(*bdsp);CHKERRQ(ierr); 376 } 377 PetscFunctionReturn(0); 378 } 379 380 PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp) 381 { 382 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 383 DM dm = sp->dm; 384 PetscInt order = sp->order; 385 PetscInt Nc = sp->Nc; 386 PetscBool continuous; 387 PetscSection csection; 388 Vec coordinates; 389 PetscReal *qpoints, *qweights; 390 PetscInt depth, dim, pdimMax, pStart, pEnd, p, *pStratStart, *pStratEnd, coneSize, d, f = 0, c; 391 PetscBool simplex, tensorSpace; 392 PetscErrorCode ierr; 393 394 PetscFunctionBegin; 395 /* Classify element type */ 396 if (!order) lag->continuous = PETSC_FALSE; 397 continuous = lag->continuous; 398 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 399 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 400 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 401 ierr = PetscCalloc1(dim+1, &lag->numDof);CHKERRQ(ierr); 402 ierr = PetscMalloc2(depth+1,&pStratStart,depth+1,&pStratEnd);CHKERRQ(ierr); 403 for (d = 0; d <= depth; ++d) {ierr = DMPlexGetDepthStratum(dm, d, &pStratStart[d], &pStratEnd[d]);CHKERRQ(ierr);} 404 ierr = DMPlexGetConeSize(dm, pStratStart[depth], &coneSize);CHKERRQ(ierr); 405 ierr = DMGetCoordinateSection(dm, &csection);CHKERRQ(ierr); 406 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 407 if (depth == 1) { 408 if (coneSize == dim+1) simplex = PETSC_TRUE; 409 else if (coneSize == 1 << dim) simplex = PETSC_FALSE; 410 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support simplices and tensor product cells"); 411 } else if (depth == dim) { 412 if (coneSize == dim+1) simplex = PETSC_TRUE; 413 else if (coneSize == 2 * dim) simplex = PETSC_FALSE; 414 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support simplices and tensor product cells"); 415 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support cell-vertex meshes or interpolated meshes"); 416 lag->simplexCell = simplex; 417 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"); 418 tensorSpace = lag->tensorSpace; 419 lag->height = 0; 420 lag->subspaces = NULL; 421 if (continuous && sp->order > 0 && dim > 0) { 422 PetscInt i; 423 424 lag->height = dim; 425 ierr = PetscMalloc1(dim,&lag->subspaces);CHKERRQ(ierr); 426 ierr = PetscDualSpaceCreateHeightSubspace_Lagrange(sp,1,&lag->subspaces[0]);CHKERRQ(ierr); 427 ierr = PetscDualSpaceSetUp(lag->subspaces[0]);CHKERRQ(ierr); 428 for (i = 1; i < dim; i++) { 429 ierr = PetscDualSpaceGetHeightSubspace(lag->subspaces[i-1],1,&lag->subspaces[i]);CHKERRQ(ierr); 430 ierr = PetscObjectReference((PetscObject)(lag->subspaces[i]));CHKERRQ(ierr); 431 } 432 } 433 ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, sp->order, &pdimMax);CHKERRQ(ierr); 434 pdimMax *= (pStratEnd[depth] - pStratStart[depth]); 435 ierr = PetscMalloc1(pdimMax, &sp->functional);CHKERRQ(ierr); 436 if (!dim) { 437 for (c = 0; c < Nc; ++c) { 438 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 439 ierr = PetscCalloc1(Nc, &qweights);CHKERRQ(ierr); 440 ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr); 441 ierr = PetscQuadratureSetData(sp->functional[f], 0, Nc, 1, NULL, qweights);CHKERRQ(ierr); 442 qweights[c] = 1.0; 443 ++f; 444 lag->numDof[0]++; 445 } 446 } else { 447 PetscInt *tup; 448 PetscReal *v0, *hv0, *J, *invJ, detJ, hdetJ; 449 PetscSection section; 450 451 ierr = PetscSectionCreate(PETSC_COMM_SELF,§ion);CHKERRQ(ierr); 452 ierr = PetscSectionSetChart(section,pStart,pEnd);CHKERRQ(ierr); 453 ierr = PetscCalloc5(dim+1,&tup,dim,&v0,dim,&hv0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 454 for (p = pStart; p < pEnd; p++) { 455 PetscInt pointDim, d, nFunc = 0; 456 PetscDualSpace hsp; 457 458 ierr = DMPlexComputeCellGeometryFEM(dm, p, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 459 for (d = 0; d < depth; d++) {if (p >= pStratStart[d] && p < pStratEnd[d]) break;} 460 pointDim = (depth == 1 && d == 1) ? dim : d; 461 hsp = ((pointDim < dim) && lag->subspaces) ? lag->subspaces[dim - pointDim - 1] : NULL; 462 if (hsp) { 463 PetscDualSpace_Lag *hlag = (PetscDualSpace_Lag *) hsp->data; 464 DM hdm; 465 466 ierr = PetscDualSpaceGetDM(hsp,&hdm);CHKERRQ(ierr); 467 ierr = DMPlexComputeCellGeometryFEM(hdm, 0, NULL, hv0, NULL, NULL, &hdetJ);CHKERRQ(ierr); 468 nFunc = lag->numDof[pointDim] = hlag->numDof[pointDim]; 469 } 470 if (pointDim == dim) { 471 /* Cells, create for self */ 472 PetscInt orderEff = continuous ? (!tensorSpace ? order-1-dim : order-2) : order; 473 PetscReal denom = continuous ? order : (!tensorSpace ? order+1+dim : order+2); 474 PetscReal numer = (!simplex || !tensorSpace) ? 2. : (2./dim); 475 PetscReal dx = numer/denom; 476 PetscInt cdim, d, d2; 477 478 if (orderEff < 0) continue; 479 ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, orderEff, &cdim);CHKERRQ(ierr); 480 ierr = PetscMemzero(tup,(dim+1)*sizeof(PetscInt));CHKERRQ(ierr); 481 if (!tensorSpace) { 482 while (!tup[dim]) { 483 for (c = 0; c < Nc; ++c) { 484 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 485 ierr = PetscMalloc1(dim, &qpoints);CHKERRQ(ierr); 486 ierr = PetscCalloc1(Nc, &qweights);CHKERRQ(ierr); 487 ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr); 488 ierr = PetscQuadratureSetData(sp->functional[f], dim, Nc, 1, qpoints, qweights);CHKERRQ(ierr); 489 for (d = 0; d < dim; ++d) { 490 qpoints[d] = v0[d]; 491 for (d2 = 0; d2 < dim; ++d2) qpoints[d] += J[d*dim+d2]*((tup[d2]+1)*dx); 492 } 493 qweights[c] = 1.0; 494 ++f; 495 } 496 ierr = LatticePointLexicographic_Internal(dim, orderEff, tup);CHKERRQ(ierr); 497 } 498 } else { 499 while (!tup[dim]) { 500 for (c = 0; c < Nc; ++c) { 501 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 502 ierr = PetscMalloc1(dim, &qpoints);CHKERRQ(ierr); 503 ierr = PetscCalloc1(Nc, &qweights);CHKERRQ(ierr); 504 ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr); 505 ierr = PetscQuadratureSetData(sp->functional[f], dim, Nc, 1, qpoints, qweights);CHKERRQ(ierr); 506 for (d = 0; d < dim; ++d) { 507 qpoints[d] = v0[d]; 508 for (d2 = 0; d2 < dim; ++d2) qpoints[d] += J[d*dim+d2]*((tup[d2]+1)*dx); 509 } 510 qweights[c] = 1.0; 511 ++f; 512 } 513 ierr = TensorPointLexicographic_Internal(dim, orderEff, tup);CHKERRQ(ierr); 514 } 515 } 516 lag->numDof[dim] = cdim; 517 } else { /* transform functionals from subspaces */ 518 PetscInt q; 519 520 for (q = 0; q < nFunc; q++, f++) { 521 PetscQuadrature fn; 522 PetscInt fdim, Nc, c, nPoints, i; 523 const PetscReal *points; 524 const PetscReal *weights; 525 PetscReal *qpoints; 526 PetscReal *qweights; 527 528 ierr = PetscDualSpaceGetFunctional(hsp, q, &fn);CHKERRQ(ierr); 529 ierr = PetscQuadratureGetData(fn,&fdim,&Nc,&nPoints,&points,&weights);CHKERRQ(ierr); 530 if (fdim != pointDim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected height dual space dim %D, got %D",pointDim,fdim); 531 ierr = PetscMalloc1(nPoints * dim, &qpoints);CHKERRQ(ierr); 532 ierr = PetscCalloc1(nPoints * Nc, &qweights);CHKERRQ(ierr); 533 for (i = 0; i < nPoints; i++) { 534 PetscInt j, k; 535 PetscReal *qp = &qpoints[i * dim]; 536 537 for (c = 0; c < Nc; ++c) qweights[i*Nc+c] = weights[i*Nc+c]; 538 for (j = 0; j < dim; ++j) qp[j] = v0[j]; 539 for (j = 0; j < dim; ++j) { 540 for (k = 0; k < pointDim; k++) qp[j] += J[dim * j + k] * (points[pointDim * i + k] - hv0[k]); 541 } 542 } 543 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 544 ierr = PetscQuadratureSetOrder(sp->functional[f],0);CHKERRQ(ierr); 545 ierr = PetscQuadratureSetData(sp->functional[f],dim,Nc,nPoints,qpoints,qweights);CHKERRQ(ierr); 546 } 547 } 548 ierr = PetscSectionSetDof(section,p,lag->numDof[pointDim]);CHKERRQ(ierr); 549 } 550 ierr = PetscFree5(tup,v0,hv0,J,invJ);CHKERRQ(ierr); 551 ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 552 { /* reorder to closure order */ 553 PetscInt *key, count; 554 PetscQuadrature *reorder = NULL; 555 556 ierr = PetscCalloc1(f,&key);CHKERRQ(ierr); 557 ierr = PetscMalloc1(f*sp->Nc,&reorder);CHKERRQ(ierr); 558 559 for (p = pStratStart[depth], count = 0; p < pStratEnd[depth]; p++) { 560 PetscInt *closure = NULL, closureSize, c; 561 562 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 563 for (c = 0; c < closureSize; c++) { 564 PetscInt point = closure[2 * c], dof, off, i; 565 566 ierr = PetscSectionGetDof(section,point,&dof);CHKERRQ(ierr); 567 ierr = PetscSectionGetOffset(section,point,&off);CHKERRQ(ierr); 568 for (i = 0; i < dof; i++) { 569 PetscInt fi = i + off; 570 if (!key[fi]) { 571 key[fi] = 1; 572 reorder[count++] = sp->functional[fi]; 573 } 574 } 575 } 576 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 577 } 578 ierr = PetscFree(sp->functional);CHKERRQ(ierr); 579 sp->functional = reorder; 580 ierr = PetscFree(key);CHKERRQ(ierr); 581 } 582 ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 583 } 584 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); 585 ierr = PetscFree2(pStratStart, pStratEnd);CHKERRQ(ierr); 586 if (f > pdimMax) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of dual basis vectors %d is greater than dimension %d", f, pdimMax); 587 PetscFunctionReturn(0); 588 } 589 590 PetscErrorCode PetscDualSpaceDestroy_Lagrange(PetscDualSpace sp) 591 { 592 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 593 PetscInt i; 594 PetscErrorCode ierr; 595 596 PetscFunctionBegin; 597 if (lag->symmetries) { 598 PetscInt **selfSyms = lag->symmetries[0]; 599 600 if (selfSyms) { 601 PetscInt i, **allocated = &selfSyms[-lag->selfSymOff]; 602 603 for (i = 0; i < lag->numSelfSym; i++) { 604 ierr = PetscFree(allocated[i]);CHKERRQ(ierr); 605 } 606 ierr = PetscFree(allocated);CHKERRQ(ierr); 607 } 608 ierr = PetscFree(lag->symmetries);CHKERRQ(ierr); 609 } 610 for (i = 0; i < lag->height; i++) { 611 ierr = PetscDualSpaceDestroy(&lag->subspaces[i]);CHKERRQ(ierr); 612 } 613 ierr = PetscFree(lag->subspaces);CHKERRQ(ierr); 614 ierr = PetscFree(lag->numDof);CHKERRQ(ierr); 615 ierr = PetscFree(lag);CHKERRQ(ierr); 616 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", NULL);CHKERRQ(ierr); 617 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", NULL);CHKERRQ(ierr); 618 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", NULL);CHKERRQ(ierr); 619 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", NULL);CHKERRQ(ierr); 620 PetscFunctionReturn(0); 621 } 622 623 PetscErrorCode PetscDualSpaceDuplicate_Lagrange(PetscDualSpace sp, PetscDualSpace *spNew) 624 { 625 PetscInt order, Nc; 626 PetscBool cont, tensor; 627 PetscErrorCode ierr; 628 629 PetscFunctionBegin; 630 ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) sp), spNew);CHKERRQ(ierr); 631 ierr = PetscDualSpaceSetType(*spNew, PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 632 ierr = PetscDualSpaceGetOrder(sp, &order);CHKERRQ(ierr); 633 ierr = PetscDualSpaceSetOrder(*spNew, order);CHKERRQ(ierr); 634 ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr); 635 ierr = PetscDualSpaceSetNumComponents(*spNew, Nc);CHKERRQ(ierr); 636 ierr = PetscDualSpaceLagrangeGetContinuity(sp, &cont);CHKERRQ(ierr); 637 ierr = PetscDualSpaceLagrangeSetContinuity(*spNew, cont);CHKERRQ(ierr); 638 ierr = PetscDualSpaceLagrangeGetTensor(sp, &tensor);CHKERRQ(ierr); 639 ierr = PetscDualSpaceLagrangeSetTensor(*spNew, tensor);CHKERRQ(ierr); 640 PetscFunctionReturn(0); 641 } 642 643 PetscErrorCode PetscDualSpaceSetFromOptions_Lagrange(PetscOptionItems *PetscOptionsObject,PetscDualSpace sp) 644 { 645 PetscBool continuous, tensor, flg; 646 PetscErrorCode ierr; 647 648 PetscFunctionBegin; 649 ierr = PetscDualSpaceLagrangeGetContinuity(sp, &continuous);CHKERRQ(ierr); 650 ierr = PetscDualSpaceLagrangeGetTensor(sp, &tensor);CHKERRQ(ierr); 651 ierr = PetscOptionsHead(PetscOptionsObject,"PetscDualSpace Lagrange Options");CHKERRQ(ierr); 652 ierr = PetscOptionsBool("-petscdualspace_lagrange_continuity", "Flag for continuous element", "PetscDualSpaceLagrangeSetContinuity", continuous, &continuous, &flg);CHKERRQ(ierr); 653 if (flg) {ierr = PetscDualSpaceLagrangeSetContinuity(sp, continuous);CHKERRQ(ierr);} 654 ierr = PetscOptionsBool("-petscdualspace_lagrange_tensor", "Flag for tensor dual space", "PetscDualSpaceLagrangeSetContinuity", tensor, &tensor, &flg);CHKERRQ(ierr); 655 if (flg) {ierr = PetscDualSpaceLagrangeSetTensor(sp, tensor);CHKERRQ(ierr);} 656 ierr = PetscOptionsTail();CHKERRQ(ierr); 657 PetscFunctionReturn(0); 658 } 659 660 PetscErrorCode PetscDualSpaceGetDimension_Lagrange(PetscDualSpace sp, PetscInt *dim) 661 { 662 DM K; 663 const PetscInt *numDof; 664 PetscInt spatialDim, Nc, size = 0, d; 665 PetscErrorCode ierr; 666 667 PetscFunctionBegin; 668 ierr = PetscDualSpaceGetDM(sp, &K);CHKERRQ(ierr); 669 ierr = PetscDualSpaceGetNumDof(sp, &numDof);CHKERRQ(ierr); 670 ierr = DMGetDimension(K, &spatialDim);CHKERRQ(ierr); 671 ierr = DMPlexGetHeightStratum(K, 0, NULL, &Nc);CHKERRQ(ierr); 672 if (Nc == 1) {ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, sp->order, dim);CHKERRQ(ierr); PetscFunctionReturn(0);} 673 for (d = 0; d <= spatialDim; ++d) { 674 PetscInt pStart, pEnd; 675 676 ierr = DMPlexGetDepthStratum(K, d, &pStart, &pEnd);CHKERRQ(ierr); 677 size += (pEnd-pStart)*numDof[d]; 678 } 679 *dim = size; 680 PetscFunctionReturn(0); 681 } 682 683 PetscErrorCode PetscDualSpaceGetNumDof_Lagrange(PetscDualSpace sp, const PetscInt **numDof) 684 { 685 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 686 687 PetscFunctionBegin; 688 *numDof = lag->numDof; 689 PetscFunctionReturn(0); 690 } 691 692 static PetscErrorCode PetscDualSpaceLagrangeGetContinuity_Lagrange(PetscDualSpace sp, PetscBool *continuous) 693 { 694 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 695 696 PetscFunctionBegin; 697 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 698 PetscValidPointer(continuous, 2); 699 *continuous = lag->continuous; 700 PetscFunctionReturn(0); 701 } 702 703 static PetscErrorCode PetscDualSpaceLagrangeSetContinuity_Lagrange(PetscDualSpace sp, PetscBool continuous) 704 { 705 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 706 707 PetscFunctionBegin; 708 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 709 lag->continuous = continuous; 710 PetscFunctionReturn(0); 711 } 712 713 /*@ 714 PetscDualSpaceLagrangeGetContinuity - Retrieves the flag for element continuity 715 716 Not Collective 717 718 Input Parameter: 719 . sp - the PetscDualSpace 720 721 Output Parameter: 722 . continuous - flag for element continuity 723 724 Level: intermediate 725 726 .keywords: PetscDualSpace, Lagrange, continuous, discontinuous 727 .seealso: PetscDualSpaceLagrangeSetContinuity() 728 @*/ 729 PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace sp, PetscBool *continuous) 730 { 731 PetscErrorCode ierr; 732 733 PetscFunctionBegin; 734 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 735 PetscValidPointer(continuous, 2); 736 ierr = PetscTryMethod(sp, "PetscDualSpaceLagrangeGetContinuity_C", (PetscDualSpace,PetscBool*),(sp,continuous));CHKERRQ(ierr); 737 PetscFunctionReturn(0); 738 } 739 740 /*@ 741 PetscDualSpaceLagrangeSetContinuity - Indicate whether the element is continuous 742 743 Logically Collective on PetscDualSpace 744 745 Input Parameters: 746 + sp - the PetscDualSpace 747 - continuous - flag for element continuity 748 749 Options Database: 750 . -petscdualspace_lagrange_continuity <bool> 751 752 Level: intermediate 753 754 .keywords: PetscDualSpace, Lagrange, continuous, discontinuous 755 .seealso: PetscDualSpaceLagrangeGetContinuity() 756 @*/ 757 PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace sp, PetscBool continuous) 758 { 759 PetscErrorCode ierr; 760 761 PetscFunctionBegin; 762 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 763 PetscValidLogicalCollectiveBool(sp, continuous, 2); 764 ierr = PetscTryMethod(sp, "PetscDualSpaceLagrangeSetContinuity_C", (PetscDualSpace,PetscBool),(sp,continuous));CHKERRQ(ierr); 765 PetscFunctionReturn(0); 766 } 767 768 PetscErrorCode PetscDualSpaceGetHeightSubspace_Lagrange(PetscDualSpace sp, PetscInt height, PetscDualSpace *bdsp) 769 { 770 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 771 PetscErrorCode ierr; 772 773 PetscFunctionBegin; 774 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 775 PetscValidPointer(bdsp,2); 776 if (height == 0) { 777 *bdsp = sp; 778 } 779 else { 780 DM dm; 781 PetscInt dim; 782 783 ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 784 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 785 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);} 786 if (height <= lag->height) { 787 *bdsp = lag->subspaces[height-1]; 788 } 789 else { 790 *bdsp = NULL; 791 } 792 } 793 PetscFunctionReturn(0); 794 } 795 796 PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp) 797 { 798 PetscFunctionBegin; 799 sp->ops->setfromoptions = PetscDualSpaceSetFromOptions_Lagrange; 800 sp->ops->setup = PetscDualSpaceSetUp_Lagrange; 801 sp->ops->view = PetscDualSpaceView_Lagrange; 802 sp->ops->destroy = PetscDualSpaceDestroy_Lagrange; 803 sp->ops->duplicate = PetscDualSpaceDuplicate_Lagrange; 804 sp->ops->getdimension = PetscDualSpaceGetDimension_Lagrange; 805 sp->ops->getnumdof = PetscDualSpaceGetNumDof_Lagrange; 806 sp->ops->getheightsubspace = PetscDualSpaceGetHeightSubspace_Lagrange; 807 sp->ops->getsymmetries = PetscDualSpaceGetSymmetries_Lagrange; 808 sp->ops->apply = PetscDualSpaceApplyDefault; 809 sp->ops->applyall = PetscDualSpaceApplyAllDefault; 810 sp->ops->createallpoints = PetscDualSpaceCreateAllPointsDefault; 811 PetscFunctionReturn(0); 812 } 813 814 /*MC 815 PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals 816 817 Level: intermediate 818 819 .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType() 820 M*/ 821 822 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp) 823 { 824 PetscDualSpace_Lag *lag; 825 PetscErrorCode ierr; 826 827 PetscFunctionBegin; 828 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 829 ierr = PetscNewLog(sp,&lag);CHKERRQ(ierr); 830 sp->data = lag; 831 832 lag->numDof = NULL; 833 lag->simplexCell = PETSC_TRUE; 834 lag->tensorSpace = PETSC_FALSE; 835 lag->continuous = PETSC_TRUE; 836 837 ierr = PetscDualSpaceInitialize_Lagrange(sp);CHKERRQ(ierr); 838 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", PetscDualSpaceLagrangeGetContinuity_Lagrange);CHKERRQ(ierr); 839 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", PetscDualSpaceLagrangeSetContinuity_Lagrange);CHKERRQ(ierr); 840 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", PetscDualSpaceLagrangeGetTensor_Lagrange);CHKERRQ(ierr); 841 ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", PetscDualSpaceLagrangeSetTensor_Lagrange);CHKERRQ(ierr); 842 PetscFunctionReturn(0); 843 } 844 845