1 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <../src/sys/utils/hash.h> 3 #include <petsc-private/isimpl.h> 4 #include <petsc-private/petscfeimpl.h> 5 #include <petscsf.h> 6 #include <petscds.h> 7 8 /** hierarchy routines */ 9 10 #undef __FUNCT__ 11 #define __FUNCT__ "DMPlexSetReferenceTree" 12 /*@ 13 DMPlexSetReferenceTree - set the reference tree for hierarchically non-conforming meshes. 14 15 Not collective 16 17 Input Parameters: 18 + dm - The DMPlex object 19 - ref - The reference tree DMPlex object 20 21 Level: intermediate 22 23 .seealso: DMPlexGetReferenceTree(), DMPlexCreateDefaultReferenceTree() 24 @*/ 25 PetscErrorCode DMPlexSetReferenceTree(DM dm, DM ref) 26 { 27 DM_Plex *mesh = (DM_Plex *)dm->data; 28 PetscErrorCode ierr; 29 30 PetscFunctionBegin; 31 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 32 PetscValidHeaderSpecific(ref, DM_CLASSID, 2); 33 ierr = PetscObjectReference((PetscObject)ref);CHKERRQ(ierr); 34 ierr = DMDestroy(&mesh->referenceTree);CHKERRQ(ierr); 35 mesh->referenceTree = ref; 36 PetscFunctionReturn(0); 37 } 38 39 #undef __FUNCT__ 40 #define __FUNCT__ "DMPlexGetReferenceTree" 41 /*@ 42 DMPlexGetReferenceTree - get the reference tree for hierarchically non-conforming meshes. 43 44 Not collective 45 46 Input Parameters: 47 . dm - The DMPlex object 48 49 Output Parameters 50 . ref - The reference tree DMPlex object 51 52 Level: intermediate 53 54 .seealso: DMPlexSetReferenceTree(), DMPlexCreateDefaultReferenceTree() 55 @*/ 56 PetscErrorCode DMPlexGetReferenceTree(DM dm, DM *ref) 57 { 58 DM_Plex *mesh = (DM_Plex *)dm->data; 59 60 PetscFunctionBegin; 61 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 62 PetscValidPointer(ref,2); 63 *ref = mesh->referenceTree; 64 PetscFunctionReturn(0); 65 } 66 67 #undef __FUNCT__ 68 #define __FUNCT__ "DMPlexReferenceTreeGetChildSymmetry_Default" 69 static PetscErrorCode DMPlexReferenceTreeGetChildSymmetry_Default(DM dm, PetscInt parent, PetscInt parentOrientA, PetscInt childOrientA, PetscInt childA, PetscInt parentOrientB, PetscInt *childOrientB, PetscInt *childB) 70 { 71 PetscInt coneSize, dStart, dEnd, dim, ABswap, oAvert, oBvert, ABswapVert; 72 PetscErrorCode ierr; 73 74 PetscFunctionBegin; 75 if (parentOrientA == parentOrientB) { 76 if (childOrientB) *childOrientB = childOrientA; 77 if (childB) *childB = childA; 78 PetscFunctionReturn(0); 79 } 80 for (dim = 0; dim < 3; dim++) { 81 ierr = DMPlexGetDepthStratum(dm,dim,&dStart,&dEnd);CHKERRQ(ierr); 82 if (parent >= dStart && parent <= dEnd) { 83 break; 84 } 85 } 86 if (dim > 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot perform child symmetry for %d-cells",dim); 87 if (!dim) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"A vertex has no children"); 88 if (childA < dStart || childA >= dEnd) { 89 /* this is a lower-dimensional child: bootstrap */ 90 PetscInt size, i, sA = -1, sB, sOrientB, sConeSize; 91 const PetscInt *supp, *coneA, *coneB, *oA, *oB; 92 93 ierr = DMPlexGetSupportSize(dm,childA,&size);CHKERRQ(ierr); 94 ierr = DMPlexGetSupport(dm,childA,&supp);CHKERRQ(ierr); 95 96 /* find a point sA in supp(childA) that has the same parent */ 97 for (i = 0; i < size; i++) { 98 PetscInt sParent; 99 100 sA = supp[i]; 101 if (sA == parent) continue; 102 ierr = DMPlexGetTreeParent(dm,sA,&sParent,NULL);CHKERRQ(ierr); 103 if (sParent == parent) { 104 break; 105 } 106 } 107 if (i == size) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"could not find support in children"); 108 /* find out which point sB is in an equivalent position to sA under 109 * parentOrientB */ 110 ierr = DMPlexReferenceTreeGetChildSymmetry_Default(dm,parent,parentOrientA,0,sA,parentOrientB,&sOrientB,&sB);CHKERRQ(ierr); 111 ierr = DMPlexGetConeSize(dm,sA,&sConeSize);CHKERRQ(ierr); 112 ierr = DMPlexGetCone(dm,sA,&coneA);CHKERRQ(ierr); 113 ierr = DMPlexGetCone(dm,sB,&coneB);CHKERRQ(ierr); 114 ierr = DMPlexGetConeOrientation(dm,sA,&oA);CHKERRQ(ierr); 115 ierr = DMPlexGetConeOrientation(dm,sB,&oB);CHKERRQ(ierr); 116 /* step through the cone of sA in natural order */ 117 for (i = 0; i < sConeSize; i++) { 118 if (coneA[i] == childA) { 119 /* if childA is at position i in coneA, 120 * then we want the point that is at sOrientB*i in coneB */ 121 PetscInt j = (sOrientB >= 0) ? ((sOrientB + i) % sConeSize) : ((sConeSize -(sOrientB+1) - i) % sConeSize); 122 if (childB) *childB = coneB[j]; 123 if (childOrientB) { 124 PetscInt oBtrue; 125 126 ierr = DMPlexGetConeSize(dm,childA,&coneSize);CHKERRQ(ierr); 127 /* compose sOrientB and oB[j] */ 128 if (coneSize != 0 && coneSize != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected a vertex or an edge"); 129 /* we may have to flip an edge */ 130 oBtrue = coneSize ? ((sOrientB >= 0) ? oB[j] : -(oB[j] + 2)) : 0; 131 ABswap = DihedralSwap(coneSize,oA[i],oBtrue);CHKERRQ(ierr); 132 *childOrientB = DihedralCompose(coneSize,childOrientA,ABswap); 133 } 134 break; 135 } 136 } 137 if (i == sConeSize) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"support cone mismatch"); 138 PetscFunctionReturn(0); 139 } 140 /* get the cone size and symmetry swap */ 141 ierr = DMPlexGetConeSize(dm,parent,&coneSize);CHKERRQ(ierr); 142 ABswap = DihedralSwap(coneSize, parentOrientA, parentOrientB); 143 if (dim == 2) { 144 /* orientations refer to cones: we want them to refer to vertices: 145 * if it's a rotation, they are the same, but if the order is reversed, a 146 * permutation that puts side i first does *not* put vertex i first */ 147 oAvert = (parentOrientA >= 0) ? parentOrientA : -((-parentOrientA % coneSize) + 1); 148 oBvert = (parentOrientB >= 0) ? parentOrientB : -((-parentOrientB % coneSize) + 1); 149 ABswapVert = DihedralSwap(coneSize, oAvert, oBvert); 150 } 151 else { 152 oAvert = parentOrientA; 153 oBvert = parentOrientB; 154 ABswapVert = ABswap; 155 } 156 if (childB) { 157 /* assume that each child corresponds to a vertex, in the same order */ 158 PetscInt p, posA = -1, numChildren, i; 159 const PetscInt *children; 160 161 /* count which position the child is in */ 162 ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr); 163 for (i = 0; i < numChildren; i++) { 164 p = children[i]; 165 if (p == childA) { 166 posA = i; 167 break; 168 } 169 } 170 if (posA >= coneSize) { 171 /* this is the triangle in the middle of a uniformly refined triangle: it 172 * is invariant */ 173 if (dim != 2 || posA != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Expected a middle triangle, got something else"); 174 *childB = childA; 175 } 176 else { 177 /* figure out position B by applying ABswapVert */ 178 PetscInt posB; 179 180 posB = (ABswapVert >= 0) ? ((ABswapVert + posA) % coneSize) : ((coneSize -(ABswapVert + 1) - posA) % coneSize); 181 if (childB) *childB = children[posB]; 182 } 183 } 184 if (childOrientB) *childOrientB = DihedralCompose(coneSize,childOrientA,ABswap); 185 PetscFunctionReturn(0); 186 } 187 188 #undef __FUNCT__ 189 #define __FUNCT__ "DMPlexReferenceTreeGetChildSymmetry" 190 /*@ 191 DMPlexReferenceTreeGetChildSymmetry - Given a reference tree, transform a childid and orientation from one parent frame to another 192 193 Input Parameters: 194 + dm - the reference tree DMPlex object 195 . parent - the parent point 196 . parentOrientA - the reference orientation for describing the parent 197 . childOrientA - the reference orientation for describing the child 198 . childA - the reference childID for describing the child 199 - parentOrientB - the new orientation for describing the parent 200 201 Output Parameters: 202 + childOrientB - if not NULL, set to the new oreintation for describing the child 203 . childB - if not NULL, the new childID for describing the child 204 205 Level: developer 206 207 .seealso: DMPlexGetReferenceTree(), DMPlexSetReferenceTree(), DMPlexSetTree() 208 @*/ 209 PetscErrorCode DMPlexReferenceTreeGetChildSymmetry(DM dm, PetscInt parent, PetscInt parentOrientA, PetscInt childOrientA, PetscInt childA, PetscInt parentOrientB, PetscInt *childOrientB, PetscInt *childB) 210 { 211 DM_Plex *mesh = (DM_Plex *)dm->data; 212 PetscErrorCode ierr; 213 214 PetscFunctionBegin; 215 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 216 if (!mesh->getchildsymmetry) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DMPlexReferenceTreeGetChildSymmetry not implemented"); 217 ierr = mesh->getchildsymmetry(dm,parent,parentOrientA,childOrientA,childA,parentOrientB,childOrientB,childB);CHKERRQ(ierr); 218 PetscFunctionReturn(0); 219 } 220 221 static PetscErrorCode DMPlexSetTree_Internal(DM,PetscSection,PetscInt*,PetscInt*,PetscBool); 222 223 #undef __FUNCT__ 224 #define __FUNCT__ "DMPlexCreateDefaultReferenceTree" 225 /*@ 226 DMPlexCreateDefaultReferenceTree - create a reference tree for isotropic hierarchical mesh refinement. 227 228 Collective on comm 229 230 Input Parameters: 231 + comm - the MPI communicator 232 . dim - the spatial dimension 233 - simplex - Flag for simplex, otherwise use a tensor-product cell 234 235 Output Parameters: 236 . ref - the reference tree DMPlex object 237 238 Level: intermediate 239 240 .keywords: reference cell 241 .seealso: DMPlexSetReferenceTree(), DMPlexGetReferenceTree() 242 @*/ 243 PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *ref) 244 { 245 DM_Plex *mesh; 246 DM K, Kref; 247 PetscInt p, pStart, pEnd, pRefStart, pRefEnd, d, offset, parentSize, *parents, *childIDs; 248 PetscInt *permvals, *unionCones, *coneSizes, *unionOrientations, numUnionPoints, *numDimPoints, numCones, numVerts; 249 DMLabel identity, identityRef; 250 PetscSection unionSection, unionConeSection, parentSection; 251 PetscScalar *unionCoords; 252 IS perm; 253 PetscErrorCode ierr; 254 255 PetscFunctionBegin; 256 /* create a reference element */ 257 ierr = DMPlexCreateReferenceCell(comm, dim, simplex, &K);CHKERRQ(ierr); 258 ierr = DMPlexCreateLabel(K, "identity");CHKERRQ(ierr); 259 ierr = DMPlexGetLabel(K, "identity", &identity);CHKERRQ(ierr); 260 ierr = DMPlexGetChart(K, &pStart, &pEnd);CHKERRQ(ierr); 261 for (p = pStart; p < pEnd; p++) { 262 ierr = DMLabelSetValue(identity, p, p);CHKERRQ(ierr); 263 } 264 /* refine it */ 265 ierr = DMRefine(K,comm,&Kref);CHKERRQ(ierr); 266 267 /* the reference tree is the union of these two, without duplicating 268 * points that appear in both */ 269 ierr = DMPlexGetLabel(Kref, "identity", &identityRef);CHKERRQ(ierr); 270 ierr = DMPlexGetChart(Kref, &pRefStart, &pRefEnd);CHKERRQ(ierr); 271 ierr = PetscSectionCreate(comm, &unionSection);CHKERRQ(ierr); 272 ierr = PetscSectionSetChart(unionSection, 0, (pEnd - pStart) + (pRefEnd - pRefStart));CHKERRQ(ierr); 273 /* count points that will go in the union */ 274 for (p = pStart; p < pEnd; p++) { 275 ierr = PetscSectionSetDof(unionSection, p - pStart, 1);CHKERRQ(ierr); 276 } 277 for (p = pRefStart; p < pRefEnd; p++) { 278 PetscInt q, qSize; 279 ierr = DMLabelGetValue(identityRef, p, &q);CHKERRQ(ierr); 280 ierr = DMLabelGetStratumSize(identityRef, q, &qSize);CHKERRQ(ierr); 281 if (qSize > 1) { 282 ierr = PetscSectionSetDof(unionSection, p - pRefStart + (pEnd - pStart), 1);CHKERRQ(ierr); 283 } 284 } 285 ierr = PetscMalloc1((pEnd - pStart) + (pRefEnd - pRefStart),&permvals);CHKERRQ(ierr); 286 offset = 0; 287 /* stratify points in the union by topological dimension */ 288 for (d = 0; d <= dim; d++) { 289 PetscInt cStart, cEnd, c; 290 291 ierr = DMPlexGetHeightStratum(K, d, &cStart, &cEnd);CHKERRQ(ierr); 292 for (c = cStart; c < cEnd; c++) { 293 permvals[offset++] = c; 294 } 295 296 ierr = DMPlexGetHeightStratum(Kref, d, &cStart, &cEnd);CHKERRQ(ierr); 297 for (c = cStart; c < cEnd; c++) { 298 permvals[offset++] = c + (pEnd - pStart); 299 } 300 } 301 ierr = ISCreateGeneral(comm, (pEnd - pStart) + (pRefEnd - pRefStart), permvals, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr); 302 ierr = PetscSectionSetPermutation(unionSection,perm);CHKERRQ(ierr); 303 ierr = PetscSectionSetUp(unionSection);CHKERRQ(ierr); 304 ierr = PetscSectionGetStorageSize(unionSection,&numUnionPoints);CHKERRQ(ierr); 305 ierr = PetscMalloc2(numUnionPoints,&coneSizes,dim+1,&numDimPoints);CHKERRQ(ierr); 306 /* count dimension points */ 307 for (d = 0; d <= dim; d++) { 308 PetscInt cStart, cOff, cOff2; 309 ierr = DMPlexGetHeightStratum(K,d,&cStart,NULL);CHKERRQ(ierr); 310 ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff);CHKERRQ(ierr); 311 if (d < dim) { 312 ierr = DMPlexGetHeightStratum(K,d+1,&cStart,NULL);CHKERRQ(ierr); 313 ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff2);CHKERRQ(ierr); 314 } 315 else { 316 cOff2 = numUnionPoints; 317 } 318 numDimPoints[dim - d] = cOff2 - cOff; 319 } 320 ierr = PetscSectionCreate(comm, &unionConeSection);CHKERRQ(ierr); 321 ierr = PetscSectionSetChart(unionConeSection, 0, numUnionPoints);CHKERRQ(ierr); 322 /* count the cones in the union */ 323 for (p = pStart; p < pEnd; p++) { 324 PetscInt dof, uOff; 325 326 ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr); 327 ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr); 328 ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr); 329 coneSizes[uOff] = dof; 330 } 331 for (p = pRefStart; p < pRefEnd; p++) { 332 PetscInt dof, uDof, uOff; 333 334 ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr); 335 ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 336 ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 337 if (uDof) { 338 ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr); 339 coneSizes[uOff] = dof; 340 } 341 } 342 ierr = PetscSectionSetUp(unionConeSection);CHKERRQ(ierr); 343 ierr = PetscSectionGetStorageSize(unionConeSection,&numCones);CHKERRQ(ierr); 344 ierr = PetscMalloc2(numCones,&unionCones,numCones,&unionOrientations);CHKERRQ(ierr); 345 /* write the cones in the union */ 346 for (p = pStart; p < pEnd; p++) { 347 PetscInt dof, uOff, c, cOff; 348 const PetscInt *cone, *orientation; 349 350 ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr); 351 ierr = DMPlexGetCone(K, p, &cone);CHKERRQ(ierr); 352 ierr = DMPlexGetConeOrientation(K, p, &orientation);CHKERRQ(ierr); 353 ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr); 354 ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr); 355 for (c = 0; c < dof; c++) { 356 PetscInt e, eOff; 357 e = cone[c]; 358 ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr); 359 unionCones[cOff + c] = eOff; 360 unionOrientations[cOff + c] = orientation[c]; 361 } 362 } 363 for (p = pRefStart; p < pRefEnd; p++) { 364 PetscInt dof, uDof, uOff, c, cOff; 365 const PetscInt *cone, *orientation; 366 367 ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr); 368 ierr = DMPlexGetCone(Kref, p, &cone);CHKERRQ(ierr); 369 ierr = DMPlexGetConeOrientation(Kref, p, &orientation);CHKERRQ(ierr); 370 ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 371 ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 372 if (uDof) { 373 ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr); 374 for (c = 0; c < dof; c++) { 375 PetscInt e, eOff, eDof; 376 377 e = cone[c]; 378 ierr = PetscSectionGetDof(unionSection, e - pRefStart + (pEnd - pStart),&eDof);CHKERRQ(ierr); 379 if (eDof) { 380 ierr = PetscSectionGetOffset(unionSection, e - pRefStart + (pEnd - pStart), &eOff);CHKERRQ(ierr); 381 } 382 else { 383 ierr = DMLabelGetValue(identityRef, e, &e);CHKERRQ(ierr); 384 ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr); 385 } 386 unionCones[cOff + c] = eOff; 387 unionOrientations[cOff + c] = orientation[c]; 388 } 389 } 390 } 391 /* get the coordinates */ 392 { 393 PetscInt vStart, vEnd, vRefStart, vRefEnd, v, vDof, vOff; 394 PetscSection KcoordsSec, KrefCoordsSec; 395 Vec KcoordsVec, KrefCoordsVec; 396 PetscScalar *Kcoords; 397 398 DMGetCoordinateSection(K, &KcoordsSec);CHKERRQ(ierr); 399 DMGetCoordinatesLocal(K, &KcoordsVec);CHKERRQ(ierr); 400 DMGetCoordinateSection(Kref, &KrefCoordsSec);CHKERRQ(ierr); 401 DMGetCoordinatesLocal(Kref, &KrefCoordsVec);CHKERRQ(ierr); 402 403 numVerts = numDimPoints[0]; 404 ierr = PetscMalloc1(numVerts * dim,&unionCoords);CHKERRQ(ierr); 405 ierr = DMPlexGetDepthStratum(K,0,&vStart,&vEnd);CHKERRQ(ierr); 406 407 offset = 0; 408 for (v = vStart; v < vEnd; v++) { 409 ierr = PetscSectionGetOffset(unionSection,v - pStart,&vOff);CHKERRQ(ierr); 410 ierr = VecGetValuesSection(KcoordsVec, KcoordsSec, v, &Kcoords);CHKERRQ(ierr); 411 for (d = 0; d < dim; d++) { 412 unionCoords[offset * dim + d] = Kcoords[d]; 413 } 414 offset++; 415 } 416 ierr = DMPlexGetDepthStratum(Kref,0,&vRefStart,&vRefEnd);CHKERRQ(ierr); 417 for (v = vRefStart; v < vRefEnd; v++) { 418 ierr = PetscSectionGetDof(unionSection,v - pRefStart + (pEnd - pStart),&vDof);CHKERRQ(ierr); 419 ierr = PetscSectionGetOffset(unionSection,v - pRefStart + (pEnd - pStart),&vOff);CHKERRQ(ierr); 420 ierr = VecGetValuesSection(KrefCoordsVec, KrefCoordsSec, v, &Kcoords);CHKERRQ(ierr); 421 if (vDof) { 422 for (d = 0; d < dim; d++) { 423 unionCoords[offset * dim + d] = Kcoords[d]; 424 } 425 offset++; 426 } 427 } 428 } 429 ierr = DMCreate(comm,ref);CHKERRQ(ierr); 430 ierr = DMSetType(*ref,DMPLEX);CHKERRQ(ierr); 431 ierr = DMPlexSetDimension(*ref,dim);CHKERRQ(ierr); 432 ierr = DMPlexCreateFromDAG(*ref,dim,numDimPoints,coneSizes,unionCones,unionOrientations,unionCoords);CHKERRQ(ierr); 433 /* set the tree */ 434 ierr = PetscSectionCreate(comm,&parentSection);CHKERRQ(ierr); 435 ierr = PetscSectionSetChart(parentSection,0,numUnionPoints);CHKERRQ(ierr); 436 for (p = pRefStart; p < pRefEnd; p++) { 437 PetscInt uDof, uOff; 438 439 ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 440 ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 441 if (uDof) { 442 PetscSectionSetDof(parentSection,uOff,1);CHKERRQ(ierr); 443 } 444 } 445 ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr); 446 ierr = PetscSectionGetStorageSize(parentSection,&parentSize);CHKERRQ(ierr); 447 ierr = PetscMalloc2(parentSize,&parents,parentSize,&childIDs);CHKERRQ(ierr); 448 for (p = pRefStart; p < pRefEnd; p++) { 449 PetscInt uDof, uOff; 450 451 ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 452 ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 453 if (uDof) { 454 PetscInt pOff, parent, parentU; 455 PetscSectionGetOffset(parentSection,uOff,&pOff);CHKERRQ(ierr); 456 DMLabelGetValue(identityRef,p,&parent);CHKERRQ(ierr); 457 ierr = PetscSectionGetOffset(unionSection, parent - pStart,&parentU);CHKERRQ(ierr); 458 parents[pOff] = parentU; 459 childIDs[pOff] = uOff; 460 } 461 } 462 ierr = DMPlexSetTree_Internal(*ref,parentSection,parents,childIDs,PETSC_TRUE);CHKERRQ(ierr); 463 mesh = (DM_Plex *) (*ref)->data; 464 mesh->getchildsymmetry = DMPlexReferenceTreeGetChildSymmetry_Default; 465 ierr = PetscSectionDestroy(&parentSection);CHKERRQ(ierr); 466 ierr = PetscFree2(parents,childIDs);CHKERRQ(ierr); 467 468 /* clean up */ 469 ierr = PetscSectionDestroy(&unionSection);CHKERRQ(ierr); 470 ierr = PetscSectionDestroy(&unionConeSection);CHKERRQ(ierr); 471 ierr = ISDestroy(&perm);CHKERRQ(ierr); 472 ierr = PetscFree(unionCoords);CHKERRQ(ierr); 473 ierr = PetscFree2(unionCones,unionOrientations);CHKERRQ(ierr); 474 ierr = PetscFree2(coneSizes,numDimPoints);CHKERRQ(ierr); 475 ierr = DMDestroy(&K);CHKERRQ(ierr); 476 ierr = DMDestroy(&Kref);CHKERRQ(ierr); 477 PetscFunctionReturn(0); 478 } 479 480 481 #undef __FUNCT__ 482 #define __FUNCT__ "DMPlexTreeSymmetrize" 483 static PetscErrorCode DMPlexTreeSymmetrize(DM dm) 484 { 485 DM_Plex *mesh = (DM_Plex *)dm->data; 486 PetscSection childSec, pSec; 487 PetscInt p, pSize, cSize, parMax = PETSC_MIN_INT, parMin = PETSC_MAX_INT; 488 PetscInt *offsets, *children, pStart, pEnd; 489 PetscErrorCode ierr; 490 491 PetscFunctionBegin; 492 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 493 ierr = PetscSectionDestroy(&mesh->childSection);CHKERRQ(ierr); 494 ierr = PetscFree(mesh->children);CHKERRQ(ierr); 495 pSec = mesh->parentSection; 496 if (!pSec) PetscFunctionReturn(0); 497 ierr = PetscSectionGetStorageSize(pSec,&pSize);CHKERRQ(ierr); 498 for (p = 0; p < pSize; p++) { 499 PetscInt par = mesh->parents[p]; 500 501 parMax = PetscMax(parMax,par+1); 502 parMin = PetscMin(parMin,par); 503 } 504 if (parMin > parMax) { 505 parMin = -1; 506 parMax = -1; 507 } 508 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)pSec),&childSec);CHKERRQ(ierr); 509 ierr = PetscSectionSetChart(childSec,parMin,parMax);CHKERRQ(ierr); 510 for (p = 0; p < pSize; p++) { 511 PetscInt par = mesh->parents[p]; 512 513 ierr = PetscSectionAddDof(childSec,par,1);CHKERRQ(ierr); 514 } 515 ierr = PetscSectionSetUp(childSec);CHKERRQ(ierr); 516 ierr = PetscSectionGetStorageSize(childSec,&cSize);CHKERRQ(ierr); 517 ierr = PetscMalloc1(cSize,&children);CHKERRQ(ierr); 518 ierr = PetscCalloc1(parMax-parMin,&offsets);CHKERRQ(ierr); 519 ierr = PetscSectionGetChart(pSec,&pStart,&pEnd);CHKERRQ(ierr); 520 for (p = pStart; p < pEnd; p++) { 521 PetscInt dof, off, i; 522 523 ierr = PetscSectionGetDof(pSec,p,&dof);CHKERRQ(ierr); 524 ierr = PetscSectionGetOffset(pSec,p,&off);CHKERRQ(ierr); 525 for (i = 0; i < dof; i++) { 526 PetscInt par = mesh->parents[off + i], cOff; 527 528 ierr = PetscSectionGetOffset(childSec,par,&cOff);CHKERRQ(ierr); 529 children[cOff + offsets[par-parMin]++] = p; 530 } 531 } 532 mesh->childSection = childSec; 533 mesh->children = children; 534 ierr = PetscFree(offsets);CHKERRQ(ierr); 535 PetscFunctionReturn(0); 536 } 537 538 #undef __FUNCT__ 539 #define __FUNCT__ "AnchorsFlatten" 540 static PetscErrorCode AnchorsFlatten (PetscSection section, IS is, PetscSection *sectionNew, IS *isNew) 541 { 542 PetscInt pStart, pEnd, size, sizeNew, i, p, *valsNew = NULL; 543 const PetscInt *vals; 544 PetscSection secNew; 545 PetscBool anyNew, globalAnyNew; 546 PetscBool compress; 547 PetscErrorCode ierr; 548 549 PetscFunctionBegin; 550 ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr); 551 ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr); 552 ierr = ISGetIndices(is,&vals);CHKERRQ(ierr); 553 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secNew);CHKERRQ(ierr); 554 ierr = PetscSectionSetChart(secNew,pStart,pEnd);CHKERRQ(ierr); 555 for (i = 0; i < size; i++) { 556 PetscInt dof; 557 558 p = vals[i]; 559 if (p < pStart || p >= pEnd) continue; 560 ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 561 if (dof) break; 562 } 563 if (i == size) { 564 ierr = PetscSectionSetUp(secNew);CHKERRQ(ierr); 565 anyNew = PETSC_FALSE; 566 compress = PETSC_FALSE; 567 sizeNew = 0; 568 } 569 else { 570 anyNew = PETSC_TRUE; 571 for (p = pStart; p < pEnd; p++) { 572 PetscInt dof, off; 573 574 ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 575 ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 576 for (i = 0; i < dof; i++) { 577 PetscInt q = vals[off + i], qDof = 0; 578 579 if (q >= pStart && q < pEnd) { 580 ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr); 581 } 582 if (qDof) { 583 ierr = PetscSectionAddDof(secNew, p, qDof);CHKERRQ(ierr); 584 } 585 else { 586 ierr = PetscSectionAddDof(secNew, p, 1);CHKERRQ(ierr); 587 } 588 } 589 } 590 ierr = PetscSectionSetUp(secNew);CHKERRQ(ierr); 591 ierr = PetscSectionGetStorageSize(secNew,&sizeNew);CHKERRQ(ierr); 592 ierr = PetscMalloc1(sizeNew,&valsNew);CHKERRQ(ierr); 593 compress = PETSC_FALSE; 594 for (p = pStart; p < pEnd; p++) { 595 PetscInt dof, off, count, offNew, dofNew; 596 597 ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 598 ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 599 ierr = PetscSectionGetDof(secNew, p, &dofNew);CHKERRQ(ierr); 600 ierr = PetscSectionGetOffset(secNew, p, &offNew);CHKERRQ(ierr); 601 count = 0; 602 for (i = 0; i < dof; i++) { 603 PetscInt q = vals[off + i], qDof = 0, qOff = 0, j; 604 605 if (q >= pStart && q < pEnd) { 606 ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr); 607 ierr = PetscSectionGetOffset(section, q, &qOff);CHKERRQ(ierr); 608 } 609 if (qDof) { 610 PetscInt oldCount = count; 611 612 for (j = 0; j < qDof; j++) { 613 PetscInt k, r = vals[qOff + j]; 614 615 for (k = 0; k < oldCount; k++) { 616 if (valsNew[offNew + k] == r) { 617 break; 618 } 619 } 620 if (k == oldCount) { 621 valsNew[offNew + count++] = r; 622 } 623 } 624 } 625 else { 626 PetscInt k, oldCount = count; 627 628 for (k = 0; k < oldCount; k++) { 629 if (valsNew[offNew + k] == q) { 630 break; 631 } 632 } 633 if (k == oldCount) { 634 valsNew[offNew + count++] = q; 635 } 636 } 637 } 638 if (count < dofNew) { 639 ierr = PetscSectionSetDof(secNew, p, count);CHKERRQ(ierr); 640 compress = PETSC_TRUE; 641 } 642 } 643 } 644 ierr = ISRestoreIndices(is,&vals);CHKERRQ(ierr); 645 ierr = MPI_Allreduce(&anyNew,&globalAnyNew,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr); 646 if (!globalAnyNew) { 647 ierr = PetscSectionDestroy(&secNew);CHKERRQ(ierr); 648 *sectionNew = NULL; 649 *isNew = NULL; 650 } 651 else { 652 PetscBool globalCompress; 653 654 ierr = MPI_Allreduce(&compress,&globalCompress,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr); 655 if (compress) { 656 PetscSection secComp; 657 PetscInt *valsComp = NULL; 658 659 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secComp);CHKERRQ(ierr); 660 ierr = PetscSectionSetChart(secComp,pStart,pEnd);CHKERRQ(ierr); 661 for (p = pStart; p < pEnd; p++) { 662 PetscInt dof; 663 664 ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr); 665 ierr = PetscSectionSetDof(secComp, p, dof);CHKERRQ(ierr); 666 } 667 ierr = PetscSectionSetUp(secComp);CHKERRQ(ierr); 668 ierr = PetscSectionGetStorageSize(secComp,&sizeNew);CHKERRQ(ierr); 669 ierr = PetscMalloc1(sizeNew,&valsComp);CHKERRQ(ierr); 670 for (p = pStart; p < pEnd; p++) { 671 PetscInt dof, off, offNew, j; 672 673 ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr); 674 ierr = PetscSectionGetOffset(secNew, p, &off);CHKERRQ(ierr); 675 ierr = PetscSectionGetOffset(secComp, p, &offNew);CHKERRQ(ierr); 676 for (j = 0; j < dof; j++) { 677 valsComp[offNew + j] = valsNew[off + j]; 678 } 679 } 680 ierr = PetscSectionDestroy(&secNew);CHKERRQ(ierr); 681 secNew = secComp; 682 ierr = PetscFree(valsNew);CHKERRQ(ierr); 683 valsNew = valsComp; 684 } 685 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),sizeNew,valsNew,PETSC_OWN_POINTER,isNew);CHKERRQ(ierr); 686 } 687 PetscFunctionReturn(0); 688 } 689 690 #undef __FUNCT__ 691 #define __FUNCT__ "DMPlexComputeConstraints_Tree" 692 static PetscErrorCode DMPlexComputeConstraints_Tree(DM dm) 693 { 694 PetscInt p, pStart, pEnd, *anchors, size; 695 PetscInt aMin = PETSC_MAX_INT, aMax = PETSC_MIN_INT; 696 PetscSection aSec; 697 DMLabel canonLabel; 698 IS aIS; 699 PetscErrorCode ierr; 700 701 PetscFunctionBegin; 702 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 703 ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 704 ierr = DMPlexGetLabel(dm,"canonical",&canonLabel);CHKERRQ(ierr); 705 for (p = pStart; p < pEnd; p++) { 706 PetscInt parent; 707 708 if (canonLabel) { 709 PetscInt canon; 710 711 ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr); 712 if (p != canon) continue; 713 } 714 ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 715 if (parent != p) { 716 aMin = PetscMin(aMin,p); 717 aMax = PetscMax(aMax,p+1); 718 } 719 } 720 if (aMin > aMax) { 721 aMin = -1; 722 aMax = -1; 723 } 724 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&aSec);CHKERRQ(ierr); 725 ierr = PetscSectionSetChart(aSec,aMin,aMax);CHKERRQ(ierr); 726 for (p = aMin; p < aMax; p++) { 727 PetscInt parent, ancestor = p; 728 729 if (canonLabel) { 730 PetscInt canon; 731 732 ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr); 733 if (p != canon) continue; 734 } 735 ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 736 while (parent != ancestor) { 737 ancestor = parent; 738 ierr = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr); 739 } 740 if (ancestor != p) { 741 PetscInt closureSize, *closure = NULL; 742 743 ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 744 ierr = PetscSectionSetDof(aSec,p,closureSize);CHKERRQ(ierr); 745 ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 746 } 747 } 748 ierr = PetscSectionSetUp(aSec);CHKERRQ(ierr); 749 ierr = PetscSectionGetStorageSize(aSec,&size);CHKERRQ(ierr); 750 ierr = PetscMalloc1(size,&anchors);CHKERRQ(ierr); 751 for (p = aMin; p < aMax; p++) { 752 PetscInt parent, ancestor = p; 753 754 if (canonLabel) { 755 PetscInt canon; 756 757 ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr); 758 if (p != canon) continue; 759 } 760 ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 761 while (parent != ancestor) { 762 ancestor = parent; 763 ierr = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr); 764 } 765 if (ancestor != p) { 766 PetscInt j, closureSize, *closure = NULL, aOff; 767 768 ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr); 769 770 ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 771 for (j = 0; j < closureSize; j++) { 772 anchors[aOff + j] = closure[2*j]; 773 } 774 ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 775 } 776 } 777 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm),size,anchors,PETSC_OWN_POINTER,&aIS);CHKERRQ(ierr); 778 { 779 PetscSection aSecNew = aSec; 780 IS aISNew = aIS; 781 782 ierr = PetscObjectReference((PetscObject)aSec);CHKERRQ(ierr); 783 ierr = PetscObjectReference((PetscObject)aIS);CHKERRQ(ierr); 784 while (aSecNew) { 785 ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr); 786 ierr = ISDestroy(&aIS);CHKERRQ(ierr); 787 aSec = aSecNew; 788 aIS = aISNew; 789 aSecNew = NULL; 790 aISNew = NULL; 791 ierr = AnchorsFlatten(aSec,aIS,&aSecNew,&aISNew);CHKERRQ(ierr); 792 } 793 } 794 ierr = DMPlexSetConstraints(dm,aSec,aIS);CHKERRQ(ierr); 795 ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr); 796 ierr = ISDestroy(&aIS);CHKERRQ(ierr); 797 PetscFunctionReturn(0); 798 } 799 800 #undef __FUNCT__ 801 #define __FUNCT__ "DMPlexSetTree_Internal" 802 static PetscErrorCode DMPlexSetTree_Internal(DM dm, PetscSection parentSection, PetscInt *parents, PetscInt *childIDs, PetscBool computeCanonical) 803 { 804 DM_Plex *mesh = (DM_Plex *)dm->data; 805 DM refTree; 806 PetscInt size; 807 PetscErrorCode ierr; 808 809 PetscFunctionBegin; 810 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 811 PetscValidHeaderSpecific(parentSection, PETSC_SECTION_CLASSID, 2); 812 ierr = PetscObjectReference((PetscObject)parentSection);CHKERRQ(ierr); 813 ierr = PetscSectionDestroy(&mesh->parentSection);CHKERRQ(ierr); 814 mesh->parentSection = parentSection; 815 ierr = PetscSectionGetStorageSize(parentSection,&size);CHKERRQ(ierr); 816 if (parents != mesh->parents) { 817 ierr = PetscFree(mesh->parents);CHKERRQ(ierr); 818 ierr = PetscMalloc1(size,&mesh->parents);CHKERRQ(ierr); 819 ierr = PetscMemcpy(mesh->parents, parents, size * sizeof(*parents));CHKERRQ(ierr); 820 } 821 if (childIDs != mesh->childIDs) { 822 ierr = PetscFree(mesh->childIDs);CHKERRQ(ierr); 823 ierr = PetscMalloc1(size,&mesh->childIDs);CHKERRQ(ierr); 824 ierr = PetscMemcpy(mesh->childIDs, childIDs, size * sizeof(*childIDs));CHKERRQ(ierr); 825 } 826 ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr); 827 if (refTree) { 828 DMLabel canonLabel; 829 830 ierr = DMPlexGetLabel(refTree,"canonical",&canonLabel);CHKERRQ(ierr); 831 if (canonLabel) { 832 PetscInt i; 833 834 for (i = 0; i < size; i++) { 835 PetscInt canon; 836 ierr = DMLabelGetValue(canonLabel, mesh->childIDs[i], &canon);CHKERRQ(ierr); 837 if (canon >= 0) { 838 mesh->childIDs[i] = canon; 839 } 840 } 841 } 842 } 843 ierr = DMPlexTreeSymmetrize(dm);CHKERRQ(ierr); 844 if (computeCanonical) { 845 PetscInt d, dim; 846 847 /* add the canonical label */ 848 ierr = DMPlexGetDimension(dm,&dim);CHKERRQ(ierr); 849 ierr = DMPlexCreateLabel(dm,"canonical");CHKERRQ(ierr); 850 for (d = 0; d <= dim; d++) { 851 PetscInt p, dStart, dEnd, canon = -1, cNumChildren; 852 const PetscInt *cChildren; 853 854 ierr = DMPlexGetDepthStratum(dm,d,&dStart,&dEnd);CHKERRQ(ierr); 855 for (p = dStart; p < dEnd; p++) { 856 ierr = DMPlexGetTreeChildren(dm,p,&cNumChildren,&cChildren);CHKERRQ(ierr); 857 if (cNumChildren) { 858 canon = p; 859 break; 860 } 861 } 862 if (canon == -1) continue; 863 for (p = dStart; p < dEnd; p++) { 864 PetscInt numChildren, i; 865 const PetscInt *children; 866 867 ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr); 868 if (numChildren) { 869 if (numChildren != cNumChildren) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"All parent points in a stratum should have the same number of children: %d != %d", numChildren, cNumChildren); 870 ierr = DMPlexSetLabelValue(dm,"canonical",p,canon);CHKERRQ(ierr); 871 for (i = 0; i < numChildren; i++) { 872 ierr = DMPlexSetLabelValue(dm,"canonical",children[i],cChildren[i]);CHKERRQ(ierr); 873 } 874 } 875 } 876 } 877 } 878 ierr = DMPlexComputeConstraints_Tree(dm);CHKERRQ(ierr); 879 PetscFunctionReturn(0); 880 } 881 882 #undef __FUNCT__ 883 #define __FUNCT__ "DMPlexSetTree" 884 /*@ 885 DMPlexSetTree - set the tree that describes the hierarchy of non-conforming mesh points. This routine also creates 886 the point-to-point constraints determined by the tree: a point is constained to the points in the closure of its 887 tree root. 888 889 Collective on dm 890 891 Input Parameters: 892 + dm - the DMPlex object 893 . parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section 894 offset indexes the parent and childID list; the reference count of parentSection is incremented 895 . parents - a list of the point parents; copied, can be destroyed 896 - childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then 897 the child corresponds to the point in the reference tree with index childIDs; copied, can be destroyed 898 899 Level: intermediate 900 901 .seealso: DMPlexGetTree(), DMPlexSetReferenceTree(), DMPlexSetConstraints(), DMPlexGetTreeParent(), DMPlexGetTreeChildren() 902 @*/ 903 PetscErrorCode DMPlexSetTree(DM dm, PetscSection parentSection, PetscInt parents[], PetscInt childIDs[]) 904 { 905 PetscErrorCode ierr; 906 907 PetscFunctionBegin; 908 ierr = DMPlexSetTree_Internal(dm,parentSection,parents,childIDs,PETSC_FALSE);CHKERRQ(ierr); 909 PetscFunctionReturn(0); 910 } 911 912 #undef __FUNCT__ 913 #define __FUNCT__ "DMPlexGetTree" 914 /*@ 915 DMPlexGetTree - get the tree that describes the hierarchy of non-conforming mesh points. 916 Collective on dm 917 918 Input Parameters: 919 . dm - the DMPlex object 920 921 Output Parameters: 922 + parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section 923 offset indexes the parent and childID list 924 . parents - a list of the point parents 925 . childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then 926 the child corresponds to the point in the reference tree with index childID 927 . childSection - the inverse of the parent section 928 - children - a list of the point children 929 930 Level: intermediate 931 932 .seealso: DMPlexSetTree(), DMPlexSetReferenceTree(), DMPlexSetConstraints(), DMPlexGetTreeParent(), DMPlexGetTreeChildren() 933 @*/ 934 PetscErrorCode DMPlexGetTree(DM dm, PetscSection *parentSection, PetscInt *parents[], PetscInt *childIDs[], PetscSection *childSection, PetscInt *children[]) 935 { 936 DM_Plex *mesh = (DM_Plex *)dm->data; 937 938 PetscFunctionBegin; 939 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 940 if (parentSection) *parentSection = mesh->parentSection; 941 if (parents) *parents = mesh->parents; 942 if (childIDs) *childIDs = mesh->childIDs; 943 if (childSection) *childSection = mesh->childSection; 944 if (children) *children = mesh->children; 945 PetscFunctionReturn(0); 946 } 947 948 #undef __FUNCT__ 949 #define __FUNCT__ "DMPlexGetTreeParent" 950 /*@ 951 DMPlexGetTreeParent - get the parent of a point in the tree describing the point hierarchy (not the Sieve DAG) 952 953 Input Parameters: 954 + dm - the DMPlex object 955 - point - the query point 956 957 Output Parameters: 958 + parent - if not NULL, set to the parent of the point, or the point itself if the point does not have a parent 959 - childID - if not NULL, set to the child ID of the point with respect to its parent, or 0 if the point 960 does not have a parent 961 962 Level: intermediate 963 964 .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeChildren() 965 @*/ 966 PetscErrorCode DMPlexGetTreeParent(DM dm, PetscInt point, PetscInt *parent, PetscInt *childID) 967 { 968 DM_Plex *mesh = (DM_Plex *)dm->data; 969 PetscSection pSec; 970 PetscErrorCode ierr; 971 972 PetscFunctionBegin; 973 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 974 pSec = mesh->parentSection; 975 if (pSec && point >= pSec->pStart && point < pSec->pEnd) { 976 PetscInt dof; 977 978 ierr = PetscSectionGetDof (pSec, point, &dof);CHKERRQ(ierr); 979 if (dof) { 980 PetscInt off; 981 982 ierr = PetscSectionGetOffset (pSec, point, &off);CHKERRQ(ierr); 983 if (parent) *parent = mesh->parents[off]; 984 if (childID) *childID = mesh->childIDs[off]; 985 PetscFunctionReturn(0); 986 } 987 } 988 if (parent) { 989 *parent = point; 990 } 991 if (childID) { 992 *childID = 0; 993 } 994 PetscFunctionReturn(0); 995 } 996 997 #undef __FUNCT__ 998 #define __FUNCT__ "DMPlexGetTreeChildren" 999 /*@C 1000 DMPlexGetTreeChildren - get the children of a point in the tree describing the point hierarchy (not the Sieve DAG) 1001 1002 Input Parameters: 1003 + dm - the DMPlex object 1004 - point - the query point 1005 1006 Output Parameters: 1007 + numChildren - if not NULL, set to the number of children 1008 - children - if not NULL, set to a list children, or set to NULL if the point has no children 1009 1010 Level: intermediate 1011 1012 Fortran Notes: 1013 Since it returns an array, this routine is only available in Fortran 90, and you must 1014 include petsc.h90 in your code. 1015 1016 .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeParent() 1017 @*/ 1018 PetscErrorCode DMPlexGetTreeChildren(DM dm, PetscInt point, PetscInt *numChildren, const PetscInt *children[]) 1019 { 1020 DM_Plex *mesh = (DM_Plex *)dm->data; 1021 PetscSection childSec; 1022 PetscInt dof = 0; 1023 PetscErrorCode ierr; 1024 1025 PetscFunctionBegin; 1026 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1027 childSec = mesh->childSection; 1028 if (childSec && point >= childSec->pStart && point < childSec->pEnd) { 1029 ierr = PetscSectionGetDof (childSec, point, &dof);CHKERRQ(ierr); 1030 } 1031 if (numChildren) *numChildren = dof; 1032 if (children) { 1033 if (dof) { 1034 PetscInt off; 1035 1036 ierr = PetscSectionGetOffset (childSec, point, &off);CHKERRQ(ierr); 1037 *children = &mesh->children[off]; 1038 } 1039 else { 1040 *children = NULL; 1041 } 1042 } 1043 PetscFunctionReturn(0); 1044 } 1045 1046 #undef __FUNCT__ 1047 #define __FUNCT__ "DMPlexComputeConstraintMatrix_ReferenceTree" 1048 PetscErrorCode DMPlexComputeConstraintMatrix_ReferenceTree(DM dm) 1049 { 1050 PetscDS ds; 1051 PetscInt spdim; 1052 PetscInt numFields, f, c, cStart, cEnd, pStart, pEnd, conStart, conEnd; 1053 const PetscInt *anchors; 1054 PetscSection section, cSec, aSec; 1055 Mat cMat; 1056 PetscReal *v0, *v0parent, *vtmp, *J, *Jparent, *invJparent, detJ, detJparent; 1057 IS aIS; 1058 PetscErrorCode ierr; 1059 1060 PetscFunctionBegin; 1061 ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 1062 ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 1063 ierr = PetscDSGetNumFields(ds,&numFields);CHKERRQ(ierr); 1064 ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 1065 ierr = DMGetDefaultSection(dm,§ion);CHKERRQ(ierr); 1066 ierr = DMPlexGetConstraintSection(dm,&cSec);CHKERRQ(ierr); 1067 ierr = DMPlexGetConstraints(dm,&aSec,&aIS);CHKERRQ(ierr); 1068 ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr); 1069 ierr = PetscSectionGetChart(cSec,&conStart,&conEnd);CHKERRQ(ierr); 1070 ierr = DMPlexGetConstraintMatrix(dm,&cMat);CHKERRQ(ierr); 1071 ierr = DMPlexGetDimension(dm,&spdim);CHKERRQ(ierr); 1072 ierr = PetscMalloc6(spdim,&v0,spdim,&v0parent,spdim,&vtmp,spdim*spdim,&J,spdim*spdim,&Jparent,spdim*spdim,&invJparent);CHKERRQ(ierr); 1073 1074 for (f = 0; f < numFields; f++) { 1075 PetscFE fe; 1076 PetscDualSpace space; 1077 PetscInt i, j, k, nPoints, offset; 1078 PetscInt fSize, fComp; 1079 PetscScalar *B = NULL; 1080 PetscReal *weights, *pointsRef, *pointsReal; 1081 Mat Amat, Bmat, Xmat; 1082 1083 ierr = PetscDSGetDiscretization(ds,f,(PetscObject *)(&fe));CHKERRQ(ierr); 1084 ierr = PetscFEGetDualSpace(fe,&space);CHKERRQ(ierr); 1085 ierr = PetscDualSpaceGetDimension(space,&fSize);CHKERRQ(ierr); 1086 ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr); 1087 ierr = MatCreate(PETSC_COMM_SELF,&Amat);CHKERRQ(ierr); 1088 ierr = MatSetSizes(Amat,fSize,fSize,fSize,fSize);CHKERRQ(ierr); 1089 ierr = MatSetType(Amat,MATSEQDENSE);CHKERRQ(ierr); 1090 ierr = MatSetUp(Amat);CHKERRQ(ierr); 1091 ierr = MatDuplicate(Amat,MAT_DO_NOT_COPY_VALUES,&Bmat);CHKERRQ(ierr); 1092 ierr = MatDuplicate(Amat,MAT_DO_NOT_COPY_VALUES,&Xmat);CHKERRQ(ierr); 1093 nPoints = 0; 1094 for (i = 0; i < fSize; i++) { 1095 PetscInt qPoints; 1096 PetscQuadrature quad; 1097 1098 ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 1099 ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr); 1100 nPoints += qPoints; 1101 } 1102 ierr = PetscMalloc3(nPoints,&weights,spdim*nPoints,&pointsRef,spdim*nPoints,&pointsReal);CHKERRQ(ierr); 1103 offset = 0; 1104 for (i = 0; i < fSize; i++) { 1105 PetscInt qPoints; 1106 const PetscReal *p, *w; 1107 PetscQuadrature quad; 1108 1109 ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 1110 ierr = PetscQuadratureGetData(quad,NULL,&qPoints,&p,&w);CHKERRQ(ierr); 1111 ierr = PetscMemcpy(weights+offset,w,qPoints*sizeof(*w));CHKERRQ(ierr); 1112 ierr = PetscMemcpy(pointsRef+spdim*offset,p,spdim*qPoints*sizeof(*p));CHKERRQ(ierr); 1113 offset += qPoints; 1114 } 1115 ierr = PetscFEGetTabulation(fe,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr); 1116 offset = 0; 1117 for (i = 0; i < fSize; i++) { 1118 PetscInt qPoints; 1119 PetscQuadrature quad; 1120 1121 ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 1122 ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr); 1123 for (j = 0; j < fSize; j++) { 1124 PetscScalar val = 0.; 1125 1126 for (k = 0; k < qPoints; k++) { 1127 val += B[((offset + k) * fSize + j) * fComp] * weights[k]; 1128 } 1129 ierr = MatSetValue(Amat,i,j,val,INSERT_VALUES);CHKERRQ(ierr); 1130 } 1131 offset += qPoints; 1132 } 1133 ierr = PetscFERestoreTabulation(fe,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr); 1134 ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1135 ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1136 ierr = MatLUFactor(Amat,NULL,NULL,NULL);CHKERRQ(ierr); 1137 for (c = cStart; c < cEnd; c++) { 1138 PetscInt parent; 1139 PetscInt closureSize, closureSizeP, *closure = NULL, *closureP = NULL; 1140 PetscInt *childOffsets, *parentOffsets; 1141 1142 ierr = DMPlexGetTreeParent(dm,c,&parent,NULL);CHKERRQ(ierr); 1143 if (parent == c) continue; 1144 ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1145 for (i = 0; i < closureSize; i++) { 1146 PetscInt p = closure[2*i]; 1147 PetscInt conDof; 1148 1149 if (p < conStart || p >= conEnd) continue; 1150 if (numFields > 1) { 1151 ierr = PetscSectionGetFieldDof(cSec,p,f,&conDof);CHKERRQ(ierr); 1152 } 1153 else { 1154 ierr = PetscSectionGetDof(cSec,p,&conDof);CHKERRQ(ierr); 1155 } 1156 if (conDof) break; 1157 } 1158 if (i == closureSize) { 1159 ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1160 continue; 1161 } 1162 1163 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr); 1164 ierr = DMPlexComputeCellGeometry(dm, parent, v0parent, Jparent, invJparent, &detJparent);CHKERRQ(ierr); 1165 for (i = 0; i < nPoints; i++) { 1166 CoordinatesRefToReal(spdim, spdim, v0, J, &pointsRef[i*spdim],vtmp);CHKERRQ(ierr); 1167 CoordinatesRealToRef(spdim, spdim, v0parent, invJparent, vtmp, &pointsReal[i*spdim]);CHKERRQ(ierr); 1168 } 1169 ierr = PetscFEGetTabulation(fe,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr); 1170 offset = 0; 1171 for (i = 0; i < fSize; i++) { 1172 PetscInt qPoints; 1173 PetscQuadrature quad; 1174 1175 ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 1176 ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr); 1177 for (j = 0; j < fSize; j++) { 1178 PetscScalar val = 0.; 1179 1180 for (k = 0; k < qPoints; k++) { 1181 val += B[((offset + k) * fSize + j) * fComp] * weights[k]; 1182 } 1183 MatSetValue(Bmat,i,j,val,INSERT_VALUES);CHKERRQ(ierr); 1184 } 1185 offset += qPoints; 1186 } 1187 ierr = PetscFERestoreTabulation(fe,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr); 1188 ierr = MatAssemblyBegin(Bmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1189 ierr = MatAssemblyEnd(Bmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1190 ierr = MatMatSolve(Amat,Bmat,Xmat);CHKERRQ(ierr); 1191 ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSizeP,&closureP);CHKERRQ(ierr); 1192 ierr = PetscMalloc2(closureSize+1,&childOffsets,closureSizeP+1,&parentOffsets);CHKERRQ(ierr); 1193 childOffsets[0] = 0; 1194 for (i = 0; i < closureSize; i++) { 1195 PetscInt p = closure[2*i]; 1196 PetscInt dof; 1197 1198 if (numFields > 1) { 1199 ierr = PetscSectionGetFieldDof(section,p,f,&dof);CHKERRQ(ierr); 1200 } 1201 else { 1202 ierr = PetscSectionGetDof(section,p,&dof);CHKERRQ(ierr); 1203 } 1204 childOffsets[i+1]=childOffsets[i]+dof / fComp; 1205 } 1206 parentOffsets[0] = 0; 1207 for (i = 0; i < closureSizeP; i++) { 1208 PetscInt p = closureP[2*i]; 1209 PetscInt dof; 1210 1211 if (numFields > 1) { 1212 ierr = PetscSectionGetFieldDof(section,p,f,&dof);CHKERRQ(ierr); 1213 } 1214 else { 1215 ierr = PetscSectionGetDof(section,p,&dof);CHKERRQ(ierr); 1216 } 1217 parentOffsets[i+1]=parentOffsets[i]+dof / fComp; 1218 } 1219 for (i = 0; i < closureSize; i++) { 1220 PetscInt conDof, conOff, aDof, aOff; 1221 PetscInt p = closure[2*i]; 1222 PetscInt o = closure[2*i+1]; 1223 1224 if (p < conStart || p >= conEnd) continue; 1225 if (numFields > 1) { 1226 ierr = PetscSectionGetFieldDof(cSec,p,f,&conDof);CHKERRQ(ierr); 1227 ierr = PetscSectionGetFieldOffset(cSec,p,f,&conOff);CHKERRQ(ierr); 1228 } 1229 else { 1230 ierr = PetscSectionGetDof(cSec,p,&conDof);CHKERRQ(ierr); 1231 ierr = PetscSectionGetOffset(cSec,p,&conOff);CHKERRQ(ierr); 1232 } 1233 if (!conDof) continue; 1234 ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr); 1235 ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr); 1236 for (k = 0; k < aDof; k++) { 1237 PetscInt a = anchors[aOff + k]; 1238 PetscInt aSecDof, aSecOff; 1239 1240 if (numFields > 1) { 1241 ierr = PetscSectionGetFieldDof(section,a,f,&aSecDof);CHKERRQ(ierr); 1242 ierr = PetscSectionGetFieldOffset(section,a,f,&aSecOff);CHKERRQ(ierr); 1243 } 1244 else { 1245 ierr = PetscSectionGetDof(section,a,&aSecDof);CHKERRQ(ierr); 1246 ierr = PetscSectionGetOffset(section,a,&aSecOff);CHKERRQ(ierr); 1247 } 1248 if (!aSecDof) continue; 1249 1250 for (j = 0; j < closureSizeP; j++) { 1251 PetscInt q = closureP[2*j]; 1252 PetscInt oq = closureP[2*j+1]; 1253 1254 if (q == a) { 1255 PetscInt r, s, t; 1256 1257 for (r = childOffsets[i]; r < childOffsets[i+1]; r++) { 1258 for (s = parentOffsets[j]; s < parentOffsets[j+1]; s++) { 1259 PetscScalar val; 1260 PetscInt insertCol, insertRow; 1261 1262 ierr = MatGetValue(Xmat,r,s,&val);CHKERRQ(ierr); 1263 if (o >= 0) { 1264 insertRow = conOff + fComp * (r - childOffsets[i]); 1265 } 1266 else { 1267 insertRow = conOff + fComp * (childOffsets[i + 1] - 1 - r); 1268 } 1269 if (oq >= 0) { 1270 insertCol = aSecOff + fComp * (s - parentOffsets[j]); 1271 } 1272 else { 1273 insertCol = aSecOff + fComp * (parentOffsets[j + 1] - 1 - s); 1274 } 1275 for (t = 0; t < fComp; t++) { 1276 ierr = MatSetValue(cMat,insertRow + t,insertCol + t,val,INSERT_VALUES);CHKERRQ(ierr); 1277 } 1278 } 1279 } 1280 } 1281 } 1282 } 1283 } 1284 ierr = PetscFree2(childOffsets,parentOffsets);CHKERRQ(ierr); 1285 ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1286 ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSizeP,&closureP);CHKERRQ(ierr); 1287 } 1288 ierr = MatDestroy(&Amat);CHKERRQ(ierr); 1289 ierr = MatDestroy(&Bmat);CHKERRQ(ierr); 1290 ierr = MatDestroy(&Xmat);CHKERRQ(ierr); 1291 ierr = PetscFree3(weights,pointsRef,pointsReal);CHKERRQ(ierr); 1292 } 1293 ierr = MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1294 ierr = MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1295 ierr = PetscFree6(v0,v0parent,vtmp,J,Jparent,invJparent);CHKERRQ(ierr); 1296 ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr); 1297 1298 PetscFunctionReturn(0); 1299 } 1300 1301 #undef __FUNCT__ 1302 #define __FUNCT__ "DMPlexComputeConstraintMatrix_Tree" 1303 PetscErrorCode DMPlexComputeConstraintMatrix_Tree(DM dm) 1304 { 1305 DM refTree; 1306 PetscDS ds; 1307 Mat refCmat, cMat; 1308 PetscInt numFields, f, pRefStart, pRefEnd, p, *rows, *cols, maxDof, maxAnDof, *perm, *iperm, pStart, pEnd, conStart, conEnd, **refPointFieldN; 1309 PetscScalar ***refPointFieldMats, *pointWork; 1310 PetscSection refConSec, conSec, refAnSec, anSec, section, refSection; 1311 IS refAnIS, anIS; 1312 const PetscInt *refAnchors, *anchors; 1313 PetscErrorCode ierr; 1314 1315 PetscFunctionBegin; 1316 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1317 ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 1318 ierr = PetscDSGetNumFields(ds,&numFields);CHKERRQ(ierr); 1319 ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr); 1320 ierr = DMSetDS(refTree,ds);CHKERRQ(ierr); 1321 ierr = DMPlexComputeConstraintMatrix_ReferenceTree(refTree);CHKERRQ(ierr); 1322 ierr = DMPlexGetConstraintMatrix(refTree,&refCmat);CHKERRQ(ierr); 1323 ierr = DMPlexGetChart(refTree,&pRefStart,&pRefEnd);CHKERRQ(ierr); 1324 ierr = DMPlexGetConstraints(refTree,&refAnSec,&refAnIS);CHKERRQ(ierr); 1325 ierr = DMPlexGetConstraints(dm,&anSec,&anIS);CHKERRQ(ierr); 1326 ierr = ISGetIndices(refAnIS,&refAnchors);CHKERRQ(ierr); 1327 ierr = ISGetIndices(anIS,&anchors);CHKERRQ(ierr); 1328 ierr = DMPlexGetConstraintSection(refTree,&refConSec);CHKERRQ(ierr); 1329 ierr = DMPlexGetConstraintSection(dm,&conSec);CHKERRQ(ierr); 1330 ierr = DMGetDefaultSection(refTree,&refSection);CHKERRQ(ierr); 1331 ierr = DMGetDefaultSection(dm,§ion);CHKERRQ(ierr); 1332 ierr = PetscSectionGetChart(refConSec,&pRefStart,&pRefEnd);CHKERRQ(ierr); 1333 ierr = PetscSectionGetChart(conSec,&conStart,&conEnd);CHKERRQ(ierr); 1334 ierr = PetscMalloc1(pRefEnd-pRefStart,&refPointFieldMats);CHKERRQ(ierr); 1335 ierr = PetscMalloc1(pRefEnd-pRefStart,&refPointFieldN);CHKERRQ(ierr); 1336 ierr = PetscSectionGetMaxDof(refConSec,&maxDof);CHKERRQ(ierr); 1337 ierr = PetscSectionGetMaxDof(refAnSec,&maxAnDof);CHKERRQ(ierr); 1338 ierr = PetscMalloc1(maxDof,&rows);CHKERRQ(ierr); 1339 ierr = PetscMalloc1(maxDof*maxAnDof,&cols);CHKERRQ(ierr); 1340 ierr = PetscMalloc1(maxDof*maxDof*maxAnDof,&pointWork);CHKERRQ(ierr); 1341 1342 /* step 1: get submats for every constrained point in the reference tree */ 1343 for (p = pRefStart; p < pRefEnd; p++) { 1344 PetscInt parent, closureSize, *closure = NULL, pDof; 1345 1346 ierr = DMPlexGetTreeParent(refTree,p,&parent,NULL);CHKERRQ(ierr); 1347 ierr = PetscSectionGetDof(refConSec,p,&pDof);CHKERRQ(ierr); 1348 if (!pDof || parent == p) continue; 1349 1350 ierr = PetscMalloc1(numFields,&refPointFieldMats[p-pRefStart]);CHKERRQ(ierr); 1351 ierr = PetscCalloc1(numFields,&refPointFieldN[p-pRefStart]);CHKERRQ(ierr); 1352 ierr = DMPlexGetTransitiveClosure(refTree,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1353 for (f = 0; f < numFields; f++) { 1354 PetscInt cDof, cOff, numCols, r, i, fComp; 1355 PetscFE fe; 1356 1357 ierr = PetscDSGetDiscretization(ds,f,(PetscObject *) &fe);CHKERRQ(ierr); 1358 ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr); 1359 1360 if (numFields > 1) { 1361 ierr = PetscSectionGetFieldDof(refConSec,p,f,&cDof);CHKERRQ(ierr); 1362 ierr = PetscSectionGetFieldOffset(refConSec,p,f,&cOff);CHKERRQ(ierr); 1363 } 1364 else { 1365 ierr = PetscSectionGetDof(refConSec,p,&cDof);CHKERRQ(ierr); 1366 ierr = PetscSectionGetOffset(refConSec,p,&cOff);CHKERRQ(ierr); 1367 } 1368 1369 if (!cDof) continue; 1370 for (r = 0; r < cDof; r++) { 1371 rows[r] = cOff + r; 1372 } 1373 numCols = 0; 1374 for (i = 0; i < closureSize; i++) { 1375 PetscInt q = closure[2*i]; 1376 PetscInt o = closure[2*i+1]; 1377 PetscInt aDof, aOff, j; 1378 1379 if (numFields > 1) { 1380 ierr = PetscSectionGetFieldDof(refSection,q,f,&aDof);CHKERRQ(ierr); 1381 ierr = PetscSectionGetFieldOffset(refSection,q,f,&aOff);CHKERRQ(ierr); 1382 } 1383 else { 1384 ierr = PetscSectionGetDof(refSection,q,&aDof);CHKERRQ(ierr); 1385 ierr = PetscSectionGetOffset(refSection,q,&aOff);CHKERRQ(ierr); 1386 } 1387 1388 for (j = 0; j < aDof; j++) { 1389 PetscInt node = (o >= 0) ? (j / fComp) : ((aDof - 1 - j) / fComp); 1390 PetscInt comp = (j % fComp); 1391 1392 cols[numCols++] = aOff + node * fComp + comp; 1393 } 1394 } 1395 refPointFieldN[p-pRefStart][f] = numCols; 1396 ierr = PetscMalloc1(cDof*numCols,&refPointFieldMats[p-pRefStart][f]);CHKERRQ(ierr); 1397 ierr = MatGetValues(refCmat,cDof,rows,numCols,cols,refPointFieldMats[p-pRefStart][f]);CHKERRQ(ierr); 1398 } 1399 ierr = DMPlexRestoreTransitiveClosure(refTree,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1400 } 1401 1402 /* step 2: compute the preorder */ 1403 ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 1404 ierr = PetscMalloc2(pEnd-pStart,&perm,pEnd-pStart,&iperm);CHKERRQ(ierr); 1405 for (p = pStart; p < pEnd; p++) { 1406 perm[p - pStart] = p; 1407 iperm[p - pStart] = p-pStart; 1408 } 1409 for (p = 0; p < pEnd - pStart;) { 1410 PetscInt point = perm[p]; 1411 PetscInt parent; 1412 1413 ierr = DMPlexGetTreeParent(dm,point,&parent,NULL);CHKERRQ(ierr); 1414 if (parent == point) { 1415 p++; 1416 } 1417 else { 1418 PetscInt size, closureSize, *closure = NULL, i; 1419 1420 ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1421 for (i = 0; i < closureSize; i++) { 1422 PetscInt q = closure[2*i]; 1423 if (iperm[q-pStart] > iperm[point-pStart]) { 1424 /* swap */ 1425 perm[p] = q; 1426 perm[iperm[q-pStart]] = point; 1427 iperm[point-pStart] = iperm[q-pStart]; 1428 iperm[q-pStart] = p; 1429 break; 1430 } 1431 } 1432 size = closureSize; 1433 ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1434 if (i == size) { 1435 p++; 1436 } 1437 } 1438 } 1439 1440 /* step 3: fill the constraint matrix */ 1441 /* we are going to use a preorder progressive fill strategy. Mat doesn't 1442 * allow progressive fill without assembly, so we are going to set up the 1443 * values outside of the Mat first. 1444 */ 1445 { 1446 PetscInt nRows, row, nnz; 1447 PetscBool done; 1448 const PetscInt *ia, *ja; 1449 PetscScalar *vals; 1450 1451 ierr = DMPlexGetConstraintMatrix(dm,&cMat);CHKERRQ(ierr); 1452 ierr = MatGetRowIJ(cMat,0,PETSC_FALSE,PETSC_FALSE,&nRows,&ia,&ja,&done);CHKERRQ(ierr); 1453 if (!done) SETERRQ(PetscObjectComm((PetscObject)cMat),PETSC_ERR_PLIB,"Could not get RowIJ of constraint matrix"); 1454 nnz = ia[nRows]; 1455 /* malloc and then zero rows right before we fill them: this way valgrind 1456 * can tell if we are doing progressive fill in the wrong order */ 1457 ierr = PetscMalloc1(nnz,&vals);CHKERRQ(ierr); 1458 for (p = 0; p < pEnd - pStart; p++) { 1459 PetscInt parent, childid, closureSize, *closure = NULL; 1460 PetscInt point = perm[p], pointDof; 1461 1462 ierr = DMPlexGetTreeParent(dm,point,&parent,&childid);CHKERRQ(ierr); 1463 if ((point < conStart) || (point >= conEnd) || (parent == point)) continue; 1464 ierr = PetscSectionGetDof(conSec,point,&pointDof);CHKERRQ(ierr); 1465 if (!pointDof) continue; 1466 ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1467 for (f = 0; f < numFields; f++) { 1468 PetscInt cDof, cOff, numCols, numFillCols, i, r, fComp, matOffset, offset; 1469 PetscScalar *pointMat; 1470 PetscFE fe; 1471 1472 ierr = PetscDSGetDiscretization(ds,f,(PetscObject *) &fe);CHKERRQ(ierr); 1473 ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr); 1474 1475 if (numFields > 1) { 1476 ierr = PetscSectionGetFieldDof(conSec,point,f,&cDof);CHKERRQ(ierr); 1477 ierr = PetscSectionGetFieldOffset(conSec,point,f,&cOff);CHKERRQ(ierr); 1478 } 1479 else { 1480 ierr = PetscSectionGetDof(conSec,point,&cDof);CHKERRQ(ierr); 1481 ierr = PetscSectionGetOffset(conSec,point,&cOff);CHKERRQ(ierr); 1482 } 1483 if (!cDof) continue; 1484 1485 /* make sure that every row for this point is the same size */ 1486 #if defined(PETSC_USE_DEBUG) 1487 for (r = 0; r < cDof; r++) { 1488 if (cDof > 1 && r) { 1489 if ((ia[cOff+r+1]-ia[cOff+r]) != (ia[cOff+r]-ia[cOff+r-1])) { 1490 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Two point rows have different nnz: %d vs. %d", (ia[rows[r]+1]-ia[rows[r]]), (ia[rows[r]]-ia[rows[r]-1])); 1491 } 1492 } 1493 } 1494 #endif 1495 /* zero rows */ 1496 for (i = ia[cOff] ; i< ia[cOff+cDof];i++) { 1497 vals[i] = 0.; 1498 } 1499 matOffset = ia[cOff]; 1500 numFillCols = ia[cOff+1] - matOffset; 1501 pointMat = refPointFieldMats[childid-pRefStart][f]; 1502 numCols = refPointFieldN[childid-pRefStart][f]; 1503 offset = 0; 1504 for (i = 0; i < closureSize; i++) { 1505 PetscInt q = closure[2*i]; 1506 PetscInt o = closure[2*i+1]; 1507 PetscInt aDof, aOff, j, k, qConDof, qConOff; 1508 1509 qConDof = qConOff = 0; 1510 if (numFields > 1) { 1511 ierr = PetscSectionGetFieldDof(section,q,f,&aDof);CHKERRQ(ierr); 1512 ierr = PetscSectionGetFieldOffset(section,q,f,&aOff);CHKERRQ(ierr); 1513 if (q >= conStart && q < conEnd) { 1514 ierr = PetscSectionGetFieldDof(conSec,q,f,&qConDof);CHKERRQ(ierr); 1515 ierr = PetscSectionGetFieldOffset(conSec,q,f,&qConOff);CHKERRQ(ierr); 1516 } 1517 } 1518 else { 1519 ierr = PetscSectionGetDof(section,q,&aDof);CHKERRQ(ierr); 1520 ierr = PetscSectionGetOffset(section,q,&aOff);CHKERRQ(ierr); 1521 if (q >= conStart && q < conEnd) { 1522 ierr = PetscSectionGetDof(conSec,q,&qConDof);CHKERRQ(ierr); 1523 ierr = PetscSectionGetOffset(conSec,q,&qConOff);CHKERRQ(ierr); 1524 } 1525 } 1526 if (!aDof) continue; 1527 if (qConDof) { 1528 /* this point has anchors: its rows of the matrix should already 1529 * be filled, thanks to preordering */ 1530 /* first multiply into pointWork, then set in matrix */ 1531 PetscInt aMatOffset = ia[qConOff]; 1532 PetscInt aNumFillCols = ia[qConOff + 1] - aMatOffset; 1533 for (r = 0; r < cDof; r++) { 1534 for (j = 0; j < aNumFillCols; j++) { 1535 PetscScalar inVal = 0; 1536 for (k = 0; k < aDof; k++) { 1537 PetscInt node = (o >= 0) ? (k / fComp) : ((aDof - 1 - k) / fComp); 1538 PetscInt comp = (k % fComp); 1539 PetscInt col = node * fComp + comp; 1540 1541 inVal += pointMat[r * numCols + offset + col] * vals[aMatOffset + aNumFillCols * k + j]; 1542 } 1543 pointWork[r * aNumFillCols + j] = inVal; 1544 } 1545 } 1546 /* assume that the columns are sorted, spend less time searching */ 1547 for (j = 0, k = 0; j < aNumFillCols; j++) { 1548 PetscInt col = ja[aMatOffset + j]; 1549 for (;k < numFillCols; k++) { 1550 if (ja[matOffset + k] == col) { 1551 break; 1552 } 1553 } 1554 if (k == numFillCols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"No nonzero space for (%d, %d)", cOff, col); 1555 for (r = 0; r < cDof; r++) { 1556 vals[matOffset + numFillCols * r + k] = pointWork[r * aNumFillCols + j]; 1557 } 1558 } 1559 } 1560 else { 1561 /* find where to put this portion of pointMat into the matrix */ 1562 for (k = 0; k < numFillCols; k++) { 1563 if (ja[matOffset + k] == aOff) { 1564 break; 1565 } 1566 } 1567 if (k == numFillCols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"No nonzero space for (%d, %d)", cOff, aOff); 1568 for (j = 0; j < aDof; j++) { 1569 PetscInt node = (o >= 0) ? (j / fComp) : ((aDof - 1 - j) / fComp); 1570 PetscInt comp = (j % fComp); 1571 PetscInt col = node * fComp + comp; 1572 for (r = 0; r < cDof; r++) { 1573 vals[matOffset + numFillCols * r + k + col] += pointMat[r * numCols + offset + col]; 1574 } 1575 } 1576 } 1577 offset += aDof; 1578 } 1579 } 1580 ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1581 } 1582 for (row = 0; row < nRows; row++) { 1583 ierr = MatSetValues(cMat,1,&row,ia[row+1]-ia[row],&ja[ia[row]],&vals[ia[row]],INSERT_VALUES);CHKERRQ(ierr); 1584 } 1585 ierr = MatRestoreRowIJ(cMat,0,PETSC_FALSE,PETSC_FALSE,&nRows,&ia,&ja,&done);CHKERRQ(ierr); 1586 if (!done) SETERRQ(PetscObjectComm((PetscObject)cMat),PETSC_ERR_PLIB,"Could not restore RowIJ of constraint matrix"); 1587 ierr = MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1588 ierr = MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1589 ierr = PetscFree(vals);CHKERRQ(ierr); 1590 } 1591 1592 /* clean up */ 1593 ierr = ISRestoreIndices(refAnIS,&refAnchors);CHKERRQ(ierr); 1594 ierr = ISRestoreIndices(anIS,&anchors);CHKERRQ(ierr); 1595 ierr = PetscFree2(perm,iperm);CHKERRQ(ierr); 1596 ierr = PetscFree(rows);CHKERRQ(ierr); 1597 ierr = PetscFree(cols);CHKERRQ(ierr); 1598 ierr = PetscFree(pointWork);CHKERRQ(ierr); 1599 for (p = pRefStart; p < pRefEnd; p++) { 1600 PetscInt parent, pDof; 1601 1602 ierr = DMPlexGetTreeParent(refTree,p,&parent,NULL);CHKERRQ(ierr); 1603 ierr = PetscSectionGetDof(refConSec,p,&pDof);CHKERRQ(ierr); 1604 if (!pDof || parent == p) continue; 1605 1606 for (f = 0; f < numFields; f++) { 1607 PetscInt cDof; 1608 1609 if (numFields > 1) { 1610 ierr = PetscSectionGetFieldDof(refConSec,p,f,&cDof);CHKERRQ(ierr); 1611 } 1612 else { 1613 ierr = PetscSectionGetDof(refConSec,p,&cDof);CHKERRQ(ierr); 1614 } 1615 1616 if (!cDof) continue; 1617 ierr = PetscFree(refPointFieldMats[p - pRefStart][f]);CHKERRQ(ierr); 1618 } 1619 ierr = PetscFree(refPointFieldMats[p - pRefStart]);CHKERRQ(ierr); 1620 ierr = PetscFree(refPointFieldN[p - pRefStart]);CHKERRQ(ierr); 1621 } 1622 ierr = PetscFree(refPointFieldMats);CHKERRQ(ierr); 1623 ierr = PetscFree(refPointFieldN);CHKERRQ(ierr); 1624 PetscFunctionReturn(0); 1625 } 1626 1627