1d6a7ad0dSToby Isaac #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2d6a7ad0dSToby Isaac #include <../src/sys/utils/hash.h> 3d6a7ad0dSToby Isaac #include <petsc-private/isimpl.h> 40c37af3bSToby Isaac #include <petsc-private/petscfeimpl.h> 5d6a7ad0dSToby Isaac #include <petscsf.h> 60c37af3bSToby Isaac #include <petscds.h> 7d6a7ad0dSToby Isaac 8d6a7ad0dSToby Isaac /** hierarchy routines */ 9d6a7ad0dSToby Isaac 10d6a7ad0dSToby Isaac #undef __FUNCT__ 11d6a7ad0dSToby Isaac #define __FUNCT__ "DMPlexSetReferenceTree" 12d6a7ad0dSToby Isaac /*@ 13d6a7ad0dSToby Isaac DMPlexSetReferenceTree - set the reference tree for hierarchically non-conforming meshes. 14d6a7ad0dSToby Isaac 15d6a7ad0dSToby Isaac Not collective 16d6a7ad0dSToby Isaac 17d6a7ad0dSToby Isaac Input Parameters: 18d6a7ad0dSToby Isaac + dm - The DMPlex object 19d6a7ad0dSToby Isaac - ref - The reference tree DMPlex object 20d6a7ad0dSToby Isaac 210b7167a0SToby Isaac Level: intermediate 22d6a7ad0dSToby Isaac 23da43764aSToby Isaac .seealso: DMPlexGetReferenceTree(), DMPlexCreateDefaultReferenceTree() 24d6a7ad0dSToby Isaac @*/ 25d6a7ad0dSToby Isaac PetscErrorCode DMPlexSetReferenceTree(DM dm, DM ref) 26d6a7ad0dSToby Isaac { 27d6a7ad0dSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 28d6a7ad0dSToby Isaac PetscErrorCode ierr; 29d6a7ad0dSToby Isaac 30d6a7ad0dSToby Isaac PetscFunctionBegin; 31d6a7ad0dSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 32d6a7ad0dSToby Isaac PetscValidHeaderSpecific(ref, DM_CLASSID, 2); 33d6a7ad0dSToby Isaac ierr = PetscObjectReference((PetscObject)ref);CHKERRQ(ierr); 34d6a7ad0dSToby Isaac ierr = DMDestroy(&mesh->referenceTree);CHKERRQ(ierr); 35d6a7ad0dSToby Isaac mesh->referenceTree = ref; 36d6a7ad0dSToby Isaac PetscFunctionReturn(0); 37d6a7ad0dSToby Isaac } 38d6a7ad0dSToby Isaac 39d6a7ad0dSToby Isaac #undef __FUNCT__ 40d6a7ad0dSToby Isaac #define __FUNCT__ "DMPlexGetReferenceTree" 41d6a7ad0dSToby Isaac /*@ 42d6a7ad0dSToby Isaac DMPlexGetReferenceTree - get the reference tree for hierarchically non-conforming meshes. 43d6a7ad0dSToby Isaac 44d6a7ad0dSToby Isaac Not collective 45d6a7ad0dSToby Isaac 46d6a7ad0dSToby Isaac Input Parameters: 47d6a7ad0dSToby Isaac . dm - The DMPlex object 48d6a7ad0dSToby Isaac 49d6a7ad0dSToby Isaac Output Parameters 50d6a7ad0dSToby Isaac . ref - The reference tree DMPlex object 51d6a7ad0dSToby Isaac 520b7167a0SToby Isaac Level: intermediate 53d6a7ad0dSToby Isaac 54da43764aSToby Isaac .seealso: DMPlexSetReferenceTree(), DMPlexCreateDefaultReferenceTree() 55d6a7ad0dSToby Isaac @*/ 56d6a7ad0dSToby Isaac PetscErrorCode DMPlexGetReferenceTree(DM dm, DM *ref) 57d6a7ad0dSToby Isaac { 58d6a7ad0dSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 59d6a7ad0dSToby Isaac 60d6a7ad0dSToby Isaac PetscFunctionBegin; 61d6a7ad0dSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 62d6a7ad0dSToby Isaac PetscValidPointer(ref,2); 63d6a7ad0dSToby Isaac *ref = mesh->referenceTree; 64d6a7ad0dSToby Isaac PetscFunctionReturn(0); 65d6a7ad0dSToby Isaac } 66d6a7ad0dSToby Isaac 67*f9f063d4SToby Isaac static PetscErrorCode DMPlexSetTree_Internal(DM,PetscSection,PetscInt*,PetscInt*,PetscBool); 68*f9f063d4SToby Isaac 69da43764aSToby Isaac #undef __FUNCT__ 70da43764aSToby Isaac #define __FUNCT__ "DMPlexCreateDefaultReferenceTree" 71da43764aSToby Isaac /*@ 72da43764aSToby Isaac DMPlexCreateDefaultReferenceTree - create a reference tree for isotropic hierarchical mesh refinement. 73da43764aSToby Isaac 74da43764aSToby Isaac Collective on comm 75da43764aSToby Isaac 76da43764aSToby Isaac Input Parameters: 77da43764aSToby Isaac + comm - the MPI communicator 78da43764aSToby Isaac . dim - the spatial dimension 79da43764aSToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell 80da43764aSToby Isaac 81da43764aSToby Isaac Output Parameters: 82da43764aSToby Isaac . ref - the reference tree DMPlex object 83da43764aSToby Isaac 840b7167a0SToby Isaac Level: intermediate 85da43764aSToby Isaac 86da43764aSToby Isaac .keywords: reference cell 87da43764aSToby Isaac .seealso: DMPlexSetReferenceTree(), DMPlexGetReferenceTree() 88da43764aSToby Isaac @*/ 89da43764aSToby Isaac PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *ref) 90da43764aSToby Isaac { 91da43764aSToby Isaac DM K, Kref; 9210f7e118SToby Isaac PetscInt p, pStart, pEnd, pRefStart, pRefEnd, d, offset, parentSize, *parents, *childIDs; 93da43764aSToby Isaac PetscInt *permvals, *unionCones, *coneSizes, *unionOrientations, numUnionPoints, *numDimPoints, numCones, numVerts; 94da43764aSToby Isaac DMLabel identity, identityRef; 9510f7e118SToby Isaac PetscSection unionSection, unionConeSection, parentSection; 96da43764aSToby Isaac PetscScalar *unionCoords; 97da43764aSToby Isaac IS perm; 98da43764aSToby Isaac PetscErrorCode ierr; 99da43764aSToby Isaac 100da43764aSToby Isaac PetscFunctionBegin; 101da43764aSToby Isaac /* create a reference element */ 102da43764aSToby Isaac ierr = DMPlexCreateReferenceCell(comm, dim, simplex, &K);CHKERRQ(ierr); 103da43764aSToby Isaac ierr = DMPlexCreateLabel(K, "identity");CHKERRQ(ierr); 104da43764aSToby Isaac ierr = DMPlexGetLabel(K, "identity", &identity);CHKERRQ(ierr); 105da43764aSToby Isaac ierr = DMPlexGetChart(K, &pStart, &pEnd);CHKERRQ(ierr); 106da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 107da43764aSToby Isaac ierr = DMLabelSetValue(identity, p, p);CHKERRQ(ierr); 108da43764aSToby Isaac } 109da43764aSToby Isaac /* refine it */ 110da43764aSToby Isaac ierr = DMRefine(K,comm,&Kref);CHKERRQ(ierr); 111da43764aSToby Isaac 112da43764aSToby Isaac /* the reference tree is the union of these two, without duplicating 113da43764aSToby Isaac * points that appear in both */ 114da43764aSToby Isaac ierr = DMPlexGetLabel(Kref, "identity", &identityRef);CHKERRQ(ierr); 115da43764aSToby Isaac ierr = DMPlexGetChart(Kref, &pRefStart, &pRefEnd);CHKERRQ(ierr); 116da43764aSToby Isaac ierr = PetscSectionCreate(comm, &unionSection);CHKERRQ(ierr); 117da43764aSToby Isaac ierr = PetscSectionSetChart(unionSection, 0, (pEnd - pStart) + (pRefEnd - pRefStart));CHKERRQ(ierr); 118da43764aSToby Isaac /* count points that will go in the union */ 119da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 120da43764aSToby Isaac ierr = PetscSectionSetDof(unionSection, p - pStart, 1);CHKERRQ(ierr); 121da43764aSToby Isaac } 122da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 123da43764aSToby Isaac PetscInt q, qSize; 124da43764aSToby Isaac ierr = DMLabelGetValue(identityRef, p, &q);CHKERRQ(ierr); 125da43764aSToby Isaac ierr = DMLabelGetStratumSize(identityRef, q, &qSize);CHKERRQ(ierr); 126da43764aSToby Isaac if (qSize > 1) { 127da43764aSToby Isaac ierr = PetscSectionSetDof(unionSection, p - pRefStart + (pEnd - pStart), 1);CHKERRQ(ierr); 128da43764aSToby Isaac } 129da43764aSToby Isaac } 130da43764aSToby Isaac ierr = PetscMalloc1((pEnd - pStart) + (pRefEnd - pRefStart),&permvals);CHKERRQ(ierr); 131da43764aSToby Isaac offset = 0; 132da43764aSToby Isaac /* stratify points in the union by topological dimension */ 133da43764aSToby Isaac for (d = 0; d <= dim; d++) { 134da43764aSToby Isaac PetscInt cStart, cEnd, c; 135da43764aSToby Isaac 136da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K, d, &cStart, &cEnd);CHKERRQ(ierr); 137da43764aSToby Isaac for (c = cStart; c < cEnd; c++) { 138da43764aSToby Isaac permvals[offset++] = c; 139da43764aSToby Isaac } 140da43764aSToby Isaac 141da43764aSToby Isaac ierr = DMPlexGetHeightStratum(Kref, d, &cStart, &cEnd);CHKERRQ(ierr); 142da43764aSToby Isaac for (c = cStart; c < cEnd; c++) { 143da43764aSToby Isaac permvals[offset++] = c + (pEnd - pStart); 144da43764aSToby Isaac } 145da43764aSToby Isaac } 146da43764aSToby Isaac ierr = ISCreateGeneral(comm, (pEnd - pStart) + (pRefEnd - pRefStart), permvals, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr); 147da43764aSToby Isaac ierr = PetscSectionSetPermutation(unionSection,perm);CHKERRQ(ierr); 148da43764aSToby Isaac ierr = PetscSectionSetUp(unionSection);CHKERRQ(ierr); 149da43764aSToby Isaac ierr = PetscSectionGetStorageSize(unionSection,&numUnionPoints);CHKERRQ(ierr); 150da43764aSToby Isaac ierr = PetscMalloc2(numUnionPoints,&coneSizes,dim+1,&numDimPoints);CHKERRQ(ierr); 151da43764aSToby Isaac /* count dimension points */ 152da43764aSToby Isaac for (d = 0; d <= dim; d++) { 153da43764aSToby Isaac PetscInt cStart, cOff, cOff2; 154da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K,d,&cStart,NULL);CHKERRQ(ierr); 155da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff);CHKERRQ(ierr); 156da43764aSToby Isaac if (d < dim) { 157da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K,d+1,&cStart,NULL);CHKERRQ(ierr); 158da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff2);CHKERRQ(ierr); 159da43764aSToby Isaac } 160da43764aSToby Isaac else { 161da43764aSToby Isaac cOff2 = numUnionPoints; 162da43764aSToby Isaac } 163da43764aSToby Isaac numDimPoints[dim - d] = cOff2 - cOff; 164da43764aSToby Isaac } 165da43764aSToby Isaac ierr = PetscSectionCreate(comm, &unionConeSection);CHKERRQ(ierr); 166da43764aSToby Isaac ierr = PetscSectionSetChart(unionConeSection, 0, numUnionPoints);CHKERRQ(ierr); 167da43764aSToby Isaac /* count the cones in the union */ 168da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 169da43764aSToby Isaac PetscInt dof, uOff; 170da43764aSToby Isaac 171da43764aSToby Isaac ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr); 172da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr); 173da43764aSToby Isaac ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr); 174da43764aSToby Isaac coneSizes[uOff] = dof; 175da43764aSToby Isaac } 176da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 177da43764aSToby Isaac PetscInt dof, uDof, uOff; 178da43764aSToby Isaac 179da43764aSToby Isaac ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr); 180da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 181da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 182da43764aSToby Isaac if (uDof) { 183da43764aSToby Isaac ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr); 184da43764aSToby Isaac coneSizes[uOff] = dof; 185da43764aSToby Isaac } 186da43764aSToby Isaac } 187da43764aSToby Isaac ierr = PetscSectionSetUp(unionConeSection);CHKERRQ(ierr); 188da43764aSToby Isaac ierr = PetscSectionGetStorageSize(unionConeSection,&numCones);CHKERRQ(ierr); 189da43764aSToby Isaac ierr = PetscMalloc2(numCones,&unionCones,numCones,&unionOrientations);CHKERRQ(ierr); 190da43764aSToby Isaac /* write the cones in the union */ 191da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 192da43764aSToby Isaac PetscInt dof, uOff, c, cOff; 193da43764aSToby Isaac const PetscInt *cone, *orientation; 194da43764aSToby Isaac 195da43764aSToby Isaac ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr); 196da43764aSToby Isaac ierr = DMPlexGetCone(K, p, &cone);CHKERRQ(ierr); 197da43764aSToby Isaac ierr = DMPlexGetConeOrientation(K, p, &orientation);CHKERRQ(ierr); 198da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr); 199da43764aSToby Isaac ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr); 200da43764aSToby Isaac for (c = 0; c < dof; c++) { 201da43764aSToby Isaac PetscInt e, eOff; 202da43764aSToby Isaac e = cone[c]; 203da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr); 204da43764aSToby Isaac unionCones[cOff + c] = eOff; 205da43764aSToby Isaac unionOrientations[cOff + c] = orientation[c]; 206da43764aSToby Isaac } 207da43764aSToby Isaac } 208da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 209da43764aSToby Isaac PetscInt dof, uDof, uOff, c, cOff; 210da43764aSToby Isaac const PetscInt *cone, *orientation; 211da43764aSToby Isaac 212da43764aSToby Isaac ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr); 213da43764aSToby Isaac ierr = DMPlexGetCone(Kref, p, &cone);CHKERRQ(ierr); 214da43764aSToby Isaac ierr = DMPlexGetConeOrientation(Kref, p, &orientation);CHKERRQ(ierr); 215da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 216da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 217da43764aSToby Isaac if (uDof) { 218da43764aSToby Isaac ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr); 219da43764aSToby Isaac for (c = 0; c < dof; c++) { 220da43764aSToby Isaac PetscInt e, eOff, eDof; 221da43764aSToby Isaac 222da43764aSToby Isaac e = cone[c]; 223da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, e - pRefStart + (pEnd - pStart),&eDof);CHKERRQ(ierr); 224da43764aSToby Isaac if (eDof) { 225da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pRefStart + (pEnd - pStart), &eOff);CHKERRQ(ierr); 226da43764aSToby Isaac } 227da43764aSToby Isaac else { 228da43764aSToby Isaac ierr = DMLabelGetValue(identityRef, e, &e);CHKERRQ(ierr); 229da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr); 230da43764aSToby Isaac } 231da43764aSToby Isaac unionCones[cOff + c] = eOff; 232da43764aSToby Isaac unionOrientations[cOff + c] = orientation[c]; 233da43764aSToby Isaac } 234da43764aSToby Isaac } 235da43764aSToby Isaac } 236da43764aSToby Isaac /* get the coordinates */ 237da43764aSToby Isaac { 238da43764aSToby Isaac PetscInt vStart, vEnd, vRefStart, vRefEnd, v, vDof, vOff; 239da43764aSToby Isaac PetscSection KcoordsSec, KrefCoordsSec; 240da43764aSToby Isaac Vec KcoordsVec, KrefCoordsVec; 241da43764aSToby Isaac PetscScalar *Kcoords; 242da43764aSToby Isaac 243da43764aSToby Isaac DMGetCoordinateSection(K, &KcoordsSec);CHKERRQ(ierr); 244da43764aSToby Isaac DMGetCoordinatesLocal(K, &KcoordsVec);CHKERRQ(ierr); 245da43764aSToby Isaac DMGetCoordinateSection(Kref, &KrefCoordsSec);CHKERRQ(ierr); 246da43764aSToby Isaac DMGetCoordinatesLocal(Kref, &KrefCoordsVec);CHKERRQ(ierr); 247da43764aSToby Isaac 248da43764aSToby Isaac numVerts = numDimPoints[0]; 249da43764aSToby Isaac ierr = PetscMalloc1(numVerts * dim,&unionCoords);CHKERRQ(ierr); 250da43764aSToby Isaac ierr = DMPlexGetDepthStratum(K,0,&vStart,&vEnd);CHKERRQ(ierr); 251da43764aSToby Isaac 252da43764aSToby Isaac offset = 0; 253da43764aSToby Isaac for (v = vStart; v < vEnd; v++) { 254da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,v - pStart,&vOff);CHKERRQ(ierr); 255da43764aSToby Isaac ierr = VecGetValuesSection(KcoordsVec, KcoordsSec, v, &Kcoords);CHKERRQ(ierr); 256da43764aSToby Isaac for (d = 0; d < dim; d++) { 257da43764aSToby Isaac unionCoords[offset * dim + d] = Kcoords[d]; 258da43764aSToby Isaac } 259da43764aSToby Isaac offset++; 260da43764aSToby Isaac } 261da43764aSToby Isaac ierr = DMPlexGetDepthStratum(Kref,0,&vRefStart,&vRefEnd);CHKERRQ(ierr); 262da43764aSToby Isaac for (v = vRefStart; v < vRefEnd; v++) { 263da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection,v - pRefStart + (pEnd - pStart),&vDof);CHKERRQ(ierr); 264da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,v - pRefStart + (pEnd - pStart),&vOff);CHKERRQ(ierr); 265da43764aSToby Isaac ierr = VecGetValuesSection(KrefCoordsVec, KrefCoordsSec, v, &Kcoords);CHKERRQ(ierr); 266da43764aSToby Isaac if (vDof) { 267da43764aSToby Isaac for (d = 0; d < dim; d++) { 268da43764aSToby Isaac unionCoords[offset * dim + d] = Kcoords[d]; 269da43764aSToby Isaac } 270da43764aSToby Isaac offset++; 271da43764aSToby Isaac } 272da43764aSToby Isaac } 273da43764aSToby Isaac } 274da43764aSToby Isaac ierr = DMCreate(comm,ref);CHKERRQ(ierr); 275da43764aSToby Isaac ierr = DMSetType(*ref,DMPLEX);CHKERRQ(ierr); 276da43764aSToby Isaac ierr = DMPlexSetDimension(*ref,dim);CHKERRQ(ierr); 277da43764aSToby Isaac ierr = DMPlexCreateFromDAG(*ref,dim,numDimPoints,coneSizes,unionCones,unionOrientations,unionCoords);CHKERRQ(ierr); 27810f7e118SToby Isaac /* set the tree */ 27910f7e118SToby Isaac ierr = PetscSectionCreate(comm,&parentSection);CHKERRQ(ierr); 28010f7e118SToby Isaac ierr = PetscSectionSetChart(parentSection,0,numUnionPoints);CHKERRQ(ierr); 28110f7e118SToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 28210f7e118SToby Isaac PetscInt uDof, uOff; 28310f7e118SToby Isaac 28410f7e118SToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 28510f7e118SToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 28610f7e118SToby Isaac if (uDof) { 28710f7e118SToby Isaac PetscSectionSetDof(parentSection,uOff,1);CHKERRQ(ierr); 28810f7e118SToby Isaac } 28910f7e118SToby Isaac } 29010f7e118SToby Isaac ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr); 29110f7e118SToby Isaac ierr = PetscSectionGetStorageSize(parentSection,&parentSize);CHKERRQ(ierr); 29210f7e118SToby Isaac ierr = PetscMalloc2(parentSize,&parents,parentSize,&childIDs);CHKERRQ(ierr); 29310f7e118SToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 29410f7e118SToby Isaac PetscInt uDof, uOff; 29510f7e118SToby Isaac 29610f7e118SToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 29710f7e118SToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 29810f7e118SToby Isaac if (uDof) { 29910f7e118SToby Isaac PetscInt pOff, parent, parentU; 30010f7e118SToby Isaac PetscSectionGetOffset(parentSection,uOff,&pOff);CHKERRQ(ierr); 30110f7e118SToby Isaac DMLabelGetValue(identityRef,p,&parent);CHKERRQ(ierr); 30210f7e118SToby Isaac ierr = PetscSectionGetOffset(unionSection, parent - pStart,&parentU);CHKERRQ(ierr); 30310f7e118SToby Isaac parents[pOff] = parentU; 30410f7e118SToby Isaac childIDs[pOff] = uOff; 30510f7e118SToby Isaac } 30610f7e118SToby Isaac } 307*f9f063d4SToby Isaac ierr = DMPlexSetTree_Internal(*ref,parentSection,parents,childIDs,PETSC_TRUE);CHKERRQ(ierr); 30810f7e118SToby Isaac ierr = PetscSectionDestroy(&parentSection);CHKERRQ(ierr); 30910f7e118SToby Isaac ierr = PetscFree2(parents,childIDs);CHKERRQ(ierr); 31010f7e118SToby Isaac 311da43764aSToby Isaac /* clean up */ 312da43764aSToby Isaac ierr = PetscSectionDestroy(&unionSection);CHKERRQ(ierr); 313da43764aSToby Isaac ierr = PetscSectionDestroy(&unionConeSection);CHKERRQ(ierr); 314da43764aSToby Isaac ierr = ISDestroy(&perm);CHKERRQ(ierr); 315da43764aSToby Isaac ierr = PetscFree(unionCoords);CHKERRQ(ierr); 316da43764aSToby Isaac ierr = PetscFree2(unionCones,unionOrientations);CHKERRQ(ierr); 317da43764aSToby Isaac ierr = PetscFree2(coneSizes,numDimPoints);CHKERRQ(ierr); 318da43764aSToby Isaac ierr = DMDestroy(&K);CHKERRQ(ierr); 319da43764aSToby Isaac ierr = DMDestroy(&Kref);CHKERRQ(ierr); 320da43764aSToby Isaac PetscFunctionReturn(0); 321da43764aSToby Isaac } 322da43764aSToby Isaac 323*f9f063d4SToby Isaac 324d961a43aSToby Isaac #undef __FUNCT__ 325878b19aaSToby Isaac #define __FUNCT__ "DMPlexTreeSymmetrize" 326878b19aaSToby Isaac static PetscErrorCode DMPlexTreeSymmetrize(DM dm) 327878b19aaSToby Isaac { 328878b19aaSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 329878b19aaSToby Isaac PetscSection childSec, pSec; 330878b19aaSToby Isaac PetscInt p, pSize, cSize, parMax = PETSC_MIN_INT, parMin = PETSC_MAX_INT; 331878b19aaSToby Isaac PetscInt *offsets, *children, pStart, pEnd; 332878b19aaSToby Isaac PetscErrorCode ierr; 333878b19aaSToby Isaac 334878b19aaSToby Isaac PetscFunctionBegin; 335878b19aaSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 336878b19aaSToby Isaac ierr = PetscSectionDestroy(&mesh->childSection);CHKERRQ(ierr); 337878b19aaSToby Isaac ierr = PetscFree(mesh->children);CHKERRQ(ierr); 338878b19aaSToby Isaac pSec = mesh->parentSection; 339878b19aaSToby Isaac if (!pSec) PetscFunctionReturn(0); 340878b19aaSToby Isaac ierr = PetscSectionGetStorageSize(pSec,&pSize);CHKERRQ(ierr); 341878b19aaSToby Isaac for (p = 0; p < pSize; p++) { 342878b19aaSToby Isaac PetscInt par = mesh->parents[p]; 343878b19aaSToby Isaac 344878b19aaSToby Isaac parMax = PetscMax(parMax,par+1); 345878b19aaSToby Isaac parMin = PetscMin(parMin,par); 346878b19aaSToby Isaac } 347878b19aaSToby Isaac if (parMin > parMax) { 348878b19aaSToby Isaac parMin = -1; 349878b19aaSToby Isaac parMax = -1; 350878b19aaSToby Isaac } 351878b19aaSToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)pSec),&childSec);CHKERRQ(ierr); 352878b19aaSToby Isaac ierr = PetscSectionSetChart(childSec,parMin,parMax);CHKERRQ(ierr); 353878b19aaSToby Isaac for (p = 0; p < pSize; p++) { 354878b19aaSToby Isaac PetscInt par = mesh->parents[p]; 355878b19aaSToby Isaac 356878b19aaSToby Isaac ierr = PetscSectionAddDof(childSec,par,1);CHKERRQ(ierr); 357878b19aaSToby Isaac } 358878b19aaSToby Isaac ierr = PetscSectionSetUp(childSec);CHKERRQ(ierr); 359878b19aaSToby Isaac ierr = PetscSectionGetStorageSize(childSec,&cSize);CHKERRQ(ierr); 360878b19aaSToby Isaac ierr = PetscMalloc1(cSize,&children);CHKERRQ(ierr); 361878b19aaSToby Isaac ierr = PetscCalloc1(parMax-parMin,&offsets);CHKERRQ(ierr); 362878b19aaSToby Isaac ierr = PetscSectionGetChart(pSec,&pStart,&pEnd);CHKERRQ(ierr); 363878b19aaSToby Isaac for (p = pStart; p < pEnd; p++) { 364878b19aaSToby Isaac PetscInt dof, off, i; 365878b19aaSToby Isaac 366878b19aaSToby Isaac ierr = PetscSectionGetDof(pSec,p,&dof);CHKERRQ(ierr); 367878b19aaSToby Isaac ierr = PetscSectionGetOffset(pSec,p,&off);CHKERRQ(ierr); 368878b19aaSToby Isaac for (i = 0; i < dof; i++) { 369878b19aaSToby Isaac PetscInt par = mesh->parents[off + i], cOff; 370878b19aaSToby Isaac 371878b19aaSToby Isaac ierr = PetscSectionGetOffset(childSec,par,&cOff);CHKERRQ(ierr); 372878b19aaSToby Isaac children[cOff + offsets[par-parMin]++] = p; 373878b19aaSToby Isaac } 374878b19aaSToby Isaac } 375878b19aaSToby Isaac mesh->childSection = childSec; 376878b19aaSToby Isaac mesh->children = children; 377878b19aaSToby Isaac ierr = PetscFree(offsets);CHKERRQ(ierr); 378878b19aaSToby Isaac PetscFunctionReturn(0); 379878b19aaSToby Isaac } 380878b19aaSToby Isaac 381878b19aaSToby Isaac #undef __FUNCT__ 3826dd5a8c8SToby Isaac #define __FUNCT__ "AnchorsFlatten" 3836dd5a8c8SToby Isaac static PetscErrorCode AnchorsFlatten (PetscSection section, IS is, PetscSection *sectionNew, IS *isNew) 3846dd5a8c8SToby Isaac { 3856dd5a8c8SToby Isaac PetscInt pStart, pEnd, size, sizeNew, i, p, *valsNew = NULL; 3866dd5a8c8SToby Isaac const PetscInt *vals; 3876dd5a8c8SToby Isaac PetscSection secNew; 3886dd5a8c8SToby Isaac PetscBool anyNew, globalAnyNew; 3896dd5a8c8SToby Isaac PetscBool compress; 3906dd5a8c8SToby Isaac PetscErrorCode ierr; 3916dd5a8c8SToby Isaac 3926dd5a8c8SToby Isaac PetscFunctionBegin; 3936dd5a8c8SToby Isaac ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr); 3946dd5a8c8SToby Isaac ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr); 3956dd5a8c8SToby Isaac ierr = ISGetIndices(is,&vals);CHKERRQ(ierr); 3966dd5a8c8SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secNew);CHKERRQ(ierr); 3976dd5a8c8SToby Isaac ierr = PetscSectionSetChart(secNew,pStart,pEnd);CHKERRQ(ierr); 3986dd5a8c8SToby Isaac for (i = 0; i < size; i++) { 3996dd5a8c8SToby Isaac PetscInt dof; 4006dd5a8c8SToby Isaac 4016dd5a8c8SToby Isaac p = vals[i]; 4026dd5a8c8SToby Isaac if (p < pStart || p >= pEnd) continue; 4036dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 4046dd5a8c8SToby Isaac if (dof) break; 4056dd5a8c8SToby Isaac } 4066dd5a8c8SToby Isaac if (i == size) { 4076dd5a8c8SToby Isaac ierr = PetscSectionSetUp(secNew);CHKERRQ(ierr); 4086dd5a8c8SToby Isaac anyNew = PETSC_FALSE; 4096dd5a8c8SToby Isaac compress = PETSC_FALSE; 4106dd5a8c8SToby Isaac sizeNew = 0; 4116dd5a8c8SToby Isaac } 4126dd5a8c8SToby Isaac else { 4136dd5a8c8SToby Isaac anyNew = PETSC_TRUE; 4146dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 4156dd5a8c8SToby Isaac PetscInt dof, off; 4166dd5a8c8SToby Isaac 4176dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 4186dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 4196dd5a8c8SToby Isaac for (i = 0; i < dof; i++) { 4206dd5a8c8SToby Isaac PetscInt q = vals[off + i], qDof = 0; 4216dd5a8c8SToby Isaac 4226dd5a8c8SToby Isaac if (q >= pStart && q < pEnd) { 4236dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr); 4246dd5a8c8SToby Isaac } 4256dd5a8c8SToby Isaac if (qDof) { 4266dd5a8c8SToby Isaac ierr = PetscSectionAddDof(secNew, p, qDof);CHKERRQ(ierr); 4276dd5a8c8SToby Isaac } 4286dd5a8c8SToby Isaac else { 4296dd5a8c8SToby Isaac ierr = PetscSectionAddDof(secNew, p, 1);CHKERRQ(ierr); 4306dd5a8c8SToby Isaac } 4316dd5a8c8SToby Isaac } 4326dd5a8c8SToby Isaac } 4336dd5a8c8SToby Isaac ierr = PetscSectionSetUp(secNew);CHKERRQ(ierr); 4346dd5a8c8SToby Isaac ierr = PetscSectionGetStorageSize(secNew,&sizeNew);CHKERRQ(ierr); 4356dd5a8c8SToby Isaac ierr = PetscMalloc1(sizeNew,&valsNew);CHKERRQ(ierr); 4366dd5a8c8SToby Isaac compress = PETSC_FALSE; 4376dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 4386dd5a8c8SToby Isaac PetscInt dof, off, count, offNew, dofNew; 4396dd5a8c8SToby Isaac 4406dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 4416dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 4426dd5a8c8SToby Isaac ierr = PetscSectionGetDof(secNew, p, &dofNew);CHKERRQ(ierr); 4436dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(secNew, p, &offNew);CHKERRQ(ierr); 4446dd5a8c8SToby Isaac count = 0; 4456dd5a8c8SToby Isaac for (i = 0; i < dof; i++) { 4466dd5a8c8SToby Isaac PetscInt q = vals[off + i], qDof = 0, qOff = 0, j; 4476dd5a8c8SToby Isaac 4486dd5a8c8SToby Isaac if (q >= pStart && q < pEnd) { 4496dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr); 4506dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(section, q, &qOff);CHKERRQ(ierr); 4516dd5a8c8SToby Isaac } 4526dd5a8c8SToby Isaac if (qDof) { 4536dd5a8c8SToby Isaac PetscInt oldCount = count; 4546dd5a8c8SToby Isaac 4556dd5a8c8SToby Isaac for (j = 0; j < qDof; j++) { 4566dd5a8c8SToby Isaac PetscInt k, r = vals[qOff + j]; 4576dd5a8c8SToby Isaac 4586dd5a8c8SToby Isaac for (k = 0; k < oldCount; k++) { 4596dd5a8c8SToby Isaac if (valsNew[offNew + k] == r) { 4606dd5a8c8SToby Isaac break; 4616dd5a8c8SToby Isaac } 4626dd5a8c8SToby Isaac } 4636dd5a8c8SToby Isaac if (k == oldCount) { 4646dd5a8c8SToby Isaac valsNew[offNew + count++] = r; 4656dd5a8c8SToby Isaac } 4666dd5a8c8SToby Isaac } 4676dd5a8c8SToby Isaac } 4686dd5a8c8SToby Isaac else { 4696dd5a8c8SToby Isaac PetscInt k, oldCount = count; 4706dd5a8c8SToby Isaac 4716dd5a8c8SToby Isaac for (k = 0; k < oldCount; k++) { 4726dd5a8c8SToby Isaac if (valsNew[offNew + k] == q) { 4736dd5a8c8SToby Isaac break; 4746dd5a8c8SToby Isaac } 4756dd5a8c8SToby Isaac } 4766dd5a8c8SToby Isaac if (k == oldCount) { 4776dd5a8c8SToby Isaac valsNew[offNew + count++] = q; 4786dd5a8c8SToby Isaac } 4796dd5a8c8SToby Isaac } 4806dd5a8c8SToby Isaac } 4816dd5a8c8SToby Isaac if (count < dofNew) { 4826dd5a8c8SToby Isaac ierr = PetscSectionSetDof(secNew, p, count);CHKERRQ(ierr); 4836dd5a8c8SToby Isaac compress = PETSC_TRUE; 4846dd5a8c8SToby Isaac } 4856dd5a8c8SToby Isaac } 4866dd5a8c8SToby Isaac } 4876dd5a8c8SToby Isaac ierr = ISRestoreIndices(is,&vals);CHKERRQ(ierr); 4886dd5a8c8SToby Isaac ierr = MPI_Allreduce(&anyNew,&globalAnyNew,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr); 4896dd5a8c8SToby Isaac if (!globalAnyNew) { 4906dd5a8c8SToby Isaac ierr = PetscSectionDestroy(&secNew);CHKERRQ(ierr); 4916dd5a8c8SToby Isaac *sectionNew = NULL; 4926dd5a8c8SToby Isaac *isNew = NULL; 4936dd5a8c8SToby Isaac } 4946dd5a8c8SToby Isaac else { 4956dd5a8c8SToby Isaac PetscBool globalCompress; 4966dd5a8c8SToby Isaac 4976dd5a8c8SToby Isaac ierr = MPI_Allreduce(&compress,&globalCompress,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr); 4986dd5a8c8SToby Isaac if (compress) { 4996dd5a8c8SToby Isaac PetscSection secComp; 5006dd5a8c8SToby Isaac PetscInt *valsComp = NULL; 5016dd5a8c8SToby Isaac 5026dd5a8c8SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secComp);CHKERRQ(ierr); 5036dd5a8c8SToby Isaac ierr = PetscSectionSetChart(secComp,pStart,pEnd);CHKERRQ(ierr); 5046dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 5056dd5a8c8SToby Isaac PetscInt dof; 5066dd5a8c8SToby Isaac 5076dd5a8c8SToby Isaac ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr); 5086dd5a8c8SToby Isaac ierr = PetscSectionSetDof(secComp, p, dof);CHKERRQ(ierr); 5096dd5a8c8SToby Isaac } 5106dd5a8c8SToby Isaac ierr = PetscSectionSetUp(secComp);CHKERRQ(ierr); 5116dd5a8c8SToby Isaac ierr = PetscSectionGetStorageSize(secComp,&sizeNew);CHKERRQ(ierr); 5126dd5a8c8SToby Isaac ierr = PetscMalloc1(sizeNew,&valsComp);CHKERRQ(ierr); 5136dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 5146dd5a8c8SToby Isaac PetscInt dof, off, offNew, j; 5156dd5a8c8SToby Isaac 5166dd5a8c8SToby Isaac ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr); 5176dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(secNew, p, &off);CHKERRQ(ierr); 5186dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(secComp, p, &offNew);CHKERRQ(ierr); 5196dd5a8c8SToby Isaac for (j = 0; j < dof; j++) { 5206dd5a8c8SToby Isaac valsComp[offNew + j] = valsNew[off + j]; 5216dd5a8c8SToby Isaac } 5226dd5a8c8SToby Isaac } 5236dd5a8c8SToby Isaac ierr = PetscSectionDestroy(&secNew);CHKERRQ(ierr); 5246dd5a8c8SToby Isaac secNew = secComp; 5256dd5a8c8SToby Isaac ierr = PetscFree(valsNew);CHKERRQ(ierr); 5266dd5a8c8SToby Isaac valsNew = valsComp; 5276dd5a8c8SToby Isaac } 5286dd5a8c8SToby Isaac ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),sizeNew,valsNew,PETSC_OWN_POINTER,isNew);CHKERRQ(ierr); 5296dd5a8c8SToby Isaac } 5306dd5a8c8SToby Isaac PetscFunctionReturn(0); 5316dd5a8c8SToby Isaac } 5326dd5a8c8SToby Isaac 5336dd5a8c8SToby Isaac #undef __FUNCT__ 53466af876cSToby Isaac #define __FUNCT__ "DMPlexComputeConstraints_Tree" 53566af876cSToby Isaac static PetscErrorCode DMPlexComputeConstraints_Tree(DM dm) 53666af876cSToby Isaac { 53766af876cSToby Isaac PetscInt p, pStart, pEnd, *anchors, size; 53866af876cSToby Isaac PetscInt aMin = PETSC_MAX_INT, aMax = PETSC_MIN_INT; 53966af876cSToby Isaac PetscSection aSec; 540*f9f063d4SToby Isaac DMLabel canonLabel; 54166af876cSToby Isaac IS aIS; 54266af876cSToby Isaac PetscErrorCode ierr; 54366af876cSToby Isaac 54466af876cSToby Isaac PetscFunctionBegin; 54566af876cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 54666af876cSToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 547*f9f063d4SToby Isaac ierr = DMPlexGetLabel(dm,"canonical",&canonLabel);CHKERRQ(ierr); 54866af876cSToby Isaac for (p = pStart; p < pEnd; p++) { 54966af876cSToby Isaac PetscInt parent; 55066af876cSToby Isaac 551*f9f063d4SToby Isaac if (canonLabel) { 552*f9f063d4SToby Isaac PetscInt canon; 553*f9f063d4SToby Isaac 554*f9f063d4SToby Isaac ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr); 555*f9f063d4SToby Isaac if (p != canon) continue; 556*f9f063d4SToby Isaac } 55766af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 55866af876cSToby Isaac if (parent != p) { 55966af876cSToby Isaac aMin = PetscMin(aMin,p); 56066af876cSToby Isaac aMax = PetscMax(aMax,p+1); 56166af876cSToby Isaac } 56266af876cSToby Isaac } 56366af876cSToby Isaac if (aMin > aMax) { 56466af876cSToby Isaac aMin = -1; 56566af876cSToby Isaac aMax = -1; 56666af876cSToby Isaac } 56766af876cSToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&aSec);CHKERRQ(ierr); 56866af876cSToby Isaac ierr = PetscSectionSetChart(aSec,aMin,aMax);CHKERRQ(ierr); 56966af876cSToby Isaac for (p = aMin; p < aMax; p++) { 57066af876cSToby Isaac PetscInt parent, ancestor = p; 57166af876cSToby Isaac 572*f9f063d4SToby Isaac if (canonLabel) { 573*f9f063d4SToby Isaac PetscInt canon; 574*f9f063d4SToby Isaac 575*f9f063d4SToby Isaac ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr); 576*f9f063d4SToby Isaac if (p != canon) continue; 577*f9f063d4SToby Isaac } 57866af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 57966af876cSToby Isaac while (parent != ancestor) { 58066af876cSToby Isaac ancestor = parent; 58166af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr); 58266af876cSToby Isaac } 58366af876cSToby Isaac if (ancestor != p) { 58466af876cSToby Isaac PetscInt closureSize, *closure = NULL; 58566af876cSToby Isaac 58666af876cSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 58766af876cSToby Isaac ierr = PetscSectionSetDof(aSec,p,closureSize);CHKERRQ(ierr); 58866af876cSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 58966af876cSToby Isaac } 59066af876cSToby Isaac } 59166af876cSToby Isaac ierr = PetscSectionSetUp(aSec);CHKERRQ(ierr); 59266af876cSToby Isaac ierr = PetscSectionGetStorageSize(aSec,&size);CHKERRQ(ierr); 59366af876cSToby Isaac ierr = PetscMalloc1(size,&anchors);CHKERRQ(ierr); 59466af876cSToby Isaac for (p = aMin; p < aMax; p++) { 59566af876cSToby Isaac PetscInt parent, ancestor = p; 59666af876cSToby Isaac 597*f9f063d4SToby Isaac if (canonLabel) { 598*f9f063d4SToby Isaac PetscInt canon; 599*f9f063d4SToby Isaac 600*f9f063d4SToby Isaac ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr); 601*f9f063d4SToby Isaac if (p != canon) continue; 602*f9f063d4SToby Isaac } 60366af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 60466af876cSToby Isaac while (parent != ancestor) { 60566af876cSToby Isaac ancestor = parent; 60666af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr); 60766af876cSToby Isaac } 60866af876cSToby Isaac if (ancestor != p) { 60966af876cSToby Isaac PetscInt j, closureSize, *closure = NULL, aOff; 61066af876cSToby Isaac 61166af876cSToby Isaac ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr); 61266af876cSToby Isaac 61366af876cSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 61466af876cSToby Isaac for (j = 0; j < closureSize; j++) { 61566af876cSToby Isaac anchors[aOff + j] = closure[2*j]; 61666af876cSToby Isaac } 61766af876cSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 61866af876cSToby Isaac } 61966af876cSToby Isaac } 62066af876cSToby Isaac ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm),size,anchors,PETSC_OWN_POINTER,&aIS);CHKERRQ(ierr); 6216dd5a8c8SToby Isaac { 6226dd5a8c8SToby Isaac PetscSection aSecNew = aSec; 6236dd5a8c8SToby Isaac IS aISNew = aIS; 6246dd5a8c8SToby Isaac 6256dd5a8c8SToby Isaac ierr = PetscObjectReference((PetscObject)aSec);CHKERRQ(ierr); 6266dd5a8c8SToby Isaac ierr = PetscObjectReference((PetscObject)aIS);CHKERRQ(ierr); 6276dd5a8c8SToby Isaac while (aSecNew) { 6286dd5a8c8SToby Isaac ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr); 6296dd5a8c8SToby Isaac ierr = ISDestroy(&aIS);CHKERRQ(ierr); 6306dd5a8c8SToby Isaac aSec = aSecNew; 6316dd5a8c8SToby Isaac aIS = aISNew; 6326dd5a8c8SToby Isaac aSecNew = NULL; 6336dd5a8c8SToby Isaac aISNew = NULL; 6346dd5a8c8SToby Isaac ierr = AnchorsFlatten(aSec,aIS,&aSecNew,&aISNew);CHKERRQ(ierr); 6356dd5a8c8SToby Isaac } 6366dd5a8c8SToby Isaac } 63766af876cSToby Isaac ierr = DMPlexSetConstraints(dm,aSec,aIS);CHKERRQ(ierr); 63866af876cSToby Isaac ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr); 63966af876cSToby Isaac ierr = ISDestroy(&aIS);CHKERRQ(ierr); 64066af876cSToby Isaac PetscFunctionReturn(0); 64166af876cSToby Isaac } 64266af876cSToby Isaac 64366af876cSToby Isaac #undef __FUNCT__ 644*f9f063d4SToby Isaac #define __FUNCT__ "DMPlexSetTree_Internal" 645*f9f063d4SToby Isaac static PetscErrorCode DMPlexSetTree_Internal(DM dm, PetscSection parentSection, PetscInt *parents, PetscInt *childIDs, PetscBool computeCanonical) 646*f9f063d4SToby Isaac { 647*f9f063d4SToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 648*f9f063d4SToby Isaac DM refTree; 649*f9f063d4SToby Isaac PetscInt size; 650*f9f063d4SToby Isaac PetscErrorCode ierr; 651*f9f063d4SToby Isaac 652*f9f063d4SToby Isaac PetscFunctionBegin; 653*f9f063d4SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 654*f9f063d4SToby Isaac PetscValidHeaderSpecific(parentSection, PETSC_SECTION_CLASSID, 2); 655*f9f063d4SToby Isaac ierr = PetscObjectReference((PetscObject)parentSection);CHKERRQ(ierr); 656*f9f063d4SToby Isaac ierr = PetscSectionDestroy(&mesh->parentSection);CHKERRQ(ierr); 657*f9f063d4SToby Isaac mesh->parentSection = parentSection; 658*f9f063d4SToby Isaac ierr = PetscSectionGetStorageSize(parentSection,&size);CHKERRQ(ierr); 659*f9f063d4SToby Isaac if (parents != mesh->parents) { 660*f9f063d4SToby Isaac ierr = PetscFree(mesh->parents);CHKERRQ(ierr); 661*f9f063d4SToby Isaac ierr = PetscMalloc1(size,&mesh->parents);CHKERRQ(ierr); 662*f9f063d4SToby Isaac ierr = PetscMemcpy(mesh->parents, parents, size * sizeof(*parents));CHKERRQ(ierr); 663*f9f063d4SToby Isaac } 664*f9f063d4SToby Isaac if (childIDs != mesh->childIDs) { 665*f9f063d4SToby Isaac ierr = PetscFree(mesh->childIDs);CHKERRQ(ierr); 666*f9f063d4SToby Isaac ierr = PetscMalloc1(size,&mesh->childIDs);CHKERRQ(ierr); 667*f9f063d4SToby Isaac ierr = PetscMemcpy(mesh->childIDs, childIDs, size * sizeof(*childIDs));CHKERRQ(ierr); 668*f9f063d4SToby Isaac } 669*f9f063d4SToby Isaac ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr); 670*f9f063d4SToby Isaac if (refTree) { 671*f9f063d4SToby Isaac DMLabel canonLabel; 672*f9f063d4SToby Isaac 673*f9f063d4SToby Isaac ierr = DMPlexGetLabel(refTree,"canonical",&canonLabel);CHKERRQ(ierr); 674*f9f063d4SToby Isaac if (canonLabel) { 675*f9f063d4SToby Isaac PetscInt i; 676*f9f063d4SToby Isaac 677*f9f063d4SToby Isaac for (i = 0; i < size; i++) { 678*f9f063d4SToby Isaac PetscInt canon; 679*f9f063d4SToby Isaac ierr = DMLabelGetValue(canonLabel, mesh->childIDs[i], &canon);CHKERRQ(ierr); 680*f9f063d4SToby Isaac if (canon >= 0) { 681*f9f063d4SToby Isaac mesh->childIDs[i] = canon; 682*f9f063d4SToby Isaac } 683*f9f063d4SToby Isaac } 684*f9f063d4SToby Isaac } 685*f9f063d4SToby Isaac } 686*f9f063d4SToby Isaac ierr = DMPlexTreeSymmetrize(dm);CHKERRQ(ierr); 687*f9f063d4SToby Isaac if (computeCanonical) { 688*f9f063d4SToby Isaac PetscInt d, dim; 689*f9f063d4SToby Isaac 690*f9f063d4SToby Isaac /* add the canonical label */ 691*f9f063d4SToby Isaac ierr = DMPlexGetDimension(dm,&dim);CHKERRQ(ierr); 692*f9f063d4SToby Isaac ierr = DMPlexCreateLabel(dm,"canonical");CHKERRQ(ierr); 693*f9f063d4SToby Isaac for (d = 0; d <= dim; d++) { 694*f9f063d4SToby Isaac PetscInt p, dStart, dEnd, canon = -1, cNumChildren; 695*f9f063d4SToby Isaac const PetscInt *cChildren; 696*f9f063d4SToby Isaac 697*f9f063d4SToby Isaac ierr = DMPlexGetDepthStratum(dm,d,&dStart,&dEnd);CHKERRQ(ierr); 698*f9f063d4SToby Isaac for (p = dStart; p < dEnd; p++) { 699*f9f063d4SToby Isaac ierr = DMPlexGetTreeChildren(dm,p,&cNumChildren,&cChildren);CHKERRQ(ierr); 700*f9f063d4SToby Isaac if (cNumChildren) { 701*f9f063d4SToby Isaac canon = p; 702*f9f063d4SToby Isaac break; 703*f9f063d4SToby Isaac } 704*f9f063d4SToby Isaac } 705*f9f063d4SToby Isaac if (canon == -1) continue; 706*f9f063d4SToby Isaac for (p = dStart; p < dEnd; p++) { 707*f9f063d4SToby Isaac PetscInt numChildren, i; 708*f9f063d4SToby Isaac const PetscInt *children; 709*f9f063d4SToby Isaac 710*f9f063d4SToby Isaac ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr); 711*f9f063d4SToby Isaac if (numChildren) { 712*f9f063d4SToby Isaac 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); 713*f9f063d4SToby Isaac ierr = DMPlexSetLabelValue(dm,"canonical",p,canon);CHKERRQ(ierr); 714*f9f063d4SToby Isaac for (i = 0; i < numChildren; i++) { 715*f9f063d4SToby Isaac ierr = DMPlexSetLabelValue(dm,"canonical",children[i],cChildren[i]);CHKERRQ(ierr); 716*f9f063d4SToby Isaac } 717*f9f063d4SToby Isaac } 718*f9f063d4SToby Isaac } 719*f9f063d4SToby Isaac } 720*f9f063d4SToby Isaac } 721*f9f063d4SToby Isaac ierr = DMPlexComputeConstraints_Tree(dm);CHKERRQ(ierr); 722*f9f063d4SToby Isaac PetscFunctionReturn(0); 723*f9f063d4SToby Isaac } 724*f9f063d4SToby Isaac 725*f9f063d4SToby Isaac #undef __FUNCT__ 7260b7167a0SToby Isaac #define __FUNCT__ "DMPlexSetTree" 7270b7167a0SToby Isaac /*@ 7280b7167a0SToby Isaac DMPlexSetTree - set the tree that describes the hierarchy of non-conforming mesh points. This routine also creates 7290b7167a0SToby Isaac the point-to-point constraints determined by the tree: a point is constained to the points in the closure of its 7300b7167a0SToby Isaac tree root. 7310b7167a0SToby Isaac 7320b7167a0SToby Isaac Collective on dm 7330b7167a0SToby Isaac 7340b7167a0SToby Isaac Input Parameters: 7350b7167a0SToby Isaac + dm - the DMPlex object 7360b7167a0SToby Isaac . parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section 7370b7167a0SToby Isaac offset indexes the parent and childID list; the reference count of parentSection is incremented 7380b7167a0SToby Isaac . parents - a list of the point parents; copied, can be destroyed 7390b7167a0SToby Isaac - childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then 7400b7167a0SToby Isaac the child corresponds to the point in the reference tree with index childIDs; copied, can be destroyed 7410b7167a0SToby Isaac 7420b7167a0SToby Isaac Level: intermediate 7430b7167a0SToby Isaac 744b2f41788SToby Isaac .seealso: DMPlexGetTree(), DMPlexSetReferenceTree(), DMPlexSetConstraints(), DMPlexGetTreeParent(), DMPlexGetTreeChildren() 7450b7167a0SToby Isaac @*/ 746b2f41788SToby Isaac PetscErrorCode DMPlexSetTree(DM dm, PetscSection parentSection, PetscInt parents[], PetscInt childIDs[]) 7470b7167a0SToby Isaac { 7480b7167a0SToby Isaac PetscErrorCode ierr; 7490b7167a0SToby Isaac 7500b7167a0SToby Isaac PetscFunctionBegin; 751*f9f063d4SToby Isaac ierr = DMPlexSetTree_Internal(dm,parentSection,parents,childIDs,PETSC_FALSE);CHKERRQ(ierr); 7520b7167a0SToby Isaac PetscFunctionReturn(0); 7530b7167a0SToby Isaac } 7540b7167a0SToby Isaac 7550b7167a0SToby Isaac #undef __FUNCT__ 756b2f41788SToby Isaac #define __FUNCT__ "DMPlexGetTree" 757b2f41788SToby Isaac /*@ 758b2f41788SToby Isaac DMPlexGetTree - get the tree that describes the hierarchy of non-conforming mesh points. 759b2f41788SToby Isaac Collective on dm 760b2f41788SToby Isaac 761b2f41788SToby Isaac Input Parameters: 762b2f41788SToby Isaac . dm - the DMPlex object 763b2f41788SToby Isaac 764b2f41788SToby Isaac Output Parameters: 765b2f41788SToby Isaac + parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section 766b2f41788SToby Isaac offset indexes the parent and childID list 767b2f41788SToby Isaac . parents - a list of the point parents 768b2f41788SToby Isaac . childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then 769b2f41788SToby Isaac the child corresponds to the point in the reference tree with index childID 770b2f41788SToby Isaac . childSection - the inverse of the parent section 771b2f41788SToby Isaac - children - a list of the point children 772b2f41788SToby Isaac 773b2f41788SToby Isaac Level: intermediate 774b2f41788SToby Isaac 775b2f41788SToby Isaac .seealso: DMPlexSetTree(), DMPlexSetReferenceTree(), DMPlexSetConstraints(), DMPlexGetTreeParent(), DMPlexGetTreeChildren() 776b2f41788SToby Isaac @*/ 777b2f41788SToby Isaac PetscErrorCode DMPlexGetTree(DM dm, PetscSection *parentSection, PetscInt *parents[], PetscInt *childIDs[], PetscSection *childSection, PetscInt *children[]) 778b2f41788SToby Isaac { 779b2f41788SToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 780b2f41788SToby Isaac 781b2f41788SToby Isaac PetscFunctionBegin; 782b2f41788SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 783b2f41788SToby Isaac if (parentSection) *parentSection = mesh->parentSection; 784b2f41788SToby Isaac if (parents) *parents = mesh->parents; 785b2f41788SToby Isaac if (childIDs) *childIDs = mesh->childIDs; 786b2f41788SToby Isaac if (childSection) *childSection = mesh->childSection; 787b2f41788SToby Isaac if (children) *children = mesh->children; 788b2f41788SToby Isaac PetscFunctionReturn(0); 789b2f41788SToby Isaac } 790b2f41788SToby Isaac 791b2f41788SToby Isaac #undef __FUNCT__ 792d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeParent" 793d961a43aSToby Isaac /*@ 794d961a43aSToby Isaac DMPlexGetTreeParent - get the parent of a point in the tree describing the point hierarchy (not the Sieve DAG) 795d961a43aSToby Isaac 796d961a43aSToby Isaac Input Parameters: 797d961a43aSToby Isaac + dm - the DMPlex object 798d961a43aSToby Isaac - point - the query point 799d961a43aSToby Isaac 800d961a43aSToby Isaac Output Parameters: 801d961a43aSToby Isaac + parent - if not NULL, set to the parent of the point, or the point itself if the point does not have a parent 802d961a43aSToby Isaac - childID - if not NULL, set to the child ID of the point with respect to its parent, or 0 if the point 803d961a43aSToby Isaac does not have a parent 804d961a43aSToby Isaac 805d961a43aSToby Isaac Level: intermediate 806d961a43aSToby Isaac 807d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeChildren() 808d961a43aSToby Isaac @*/ 809d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeParent(DM dm, PetscInt point, PetscInt *parent, PetscInt *childID) 810d961a43aSToby Isaac { 811d961a43aSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 812d961a43aSToby Isaac PetscSection pSec; 813d961a43aSToby Isaac PetscErrorCode ierr; 814d961a43aSToby Isaac 815d961a43aSToby Isaac PetscFunctionBegin; 816d961a43aSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 817d961a43aSToby Isaac pSec = mesh->parentSection; 818d961a43aSToby Isaac if (pSec && point >= pSec->pStart && point < pSec->pEnd) { 819d961a43aSToby Isaac PetscInt dof; 820d961a43aSToby Isaac 821d961a43aSToby Isaac ierr = PetscSectionGetDof (pSec, point, &dof);CHKERRQ(ierr); 822d961a43aSToby Isaac if (dof) { 823d961a43aSToby Isaac PetscInt off; 824d961a43aSToby Isaac 825d961a43aSToby Isaac ierr = PetscSectionGetOffset (pSec, point, &off);CHKERRQ(ierr); 826d961a43aSToby Isaac if (parent) *parent = mesh->parents[off]; 827d961a43aSToby Isaac if (childID) *childID = mesh->childIDs[off]; 828d961a43aSToby Isaac PetscFunctionReturn(0); 829d961a43aSToby Isaac } 830d961a43aSToby Isaac } 831d961a43aSToby Isaac if (parent) { 832d961a43aSToby Isaac *parent = point; 833d961a43aSToby Isaac } 834d961a43aSToby Isaac if (childID) { 835d961a43aSToby Isaac *childID = 0; 836d961a43aSToby Isaac } 837d961a43aSToby Isaac PetscFunctionReturn(0); 838d961a43aSToby Isaac } 839d961a43aSToby Isaac 840d961a43aSToby Isaac #undef __FUNCT__ 841d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeChildren" 842d961a43aSToby Isaac /*@C 843d961a43aSToby Isaac DMPlexGetTreeChildren - get the children of a point in the tree describing the point hierarchy (not the Sieve DAG) 844d961a43aSToby Isaac 845d961a43aSToby Isaac Input Parameters: 846d961a43aSToby Isaac + dm - the DMPlex object 847d961a43aSToby Isaac - point - the query point 848d961a43aSToby Isaac 849d961a43aSToby Isaac Output Parameters: 850d961a43aSToby Isaac + numChildren - if not NULL, set to the number of children 851d961a43aSToby Isaac - children - if not NULL, set to a list children, or set to NULL if the point has no children 852d961a43aSToby Isaac 853d961a43aSToby Isaac Level: intermediate 854d961a43aSToby Isaac 855d961a43aSToby Isaac Fortran Notes: 856d961a43aSToby Isaac Since it returns an array, this routine is only available in Fortran 90, and you must 857d961a43aSToby Isaac include petsc.h90 in your code. 858d961a43aSToby Isaac 859d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeParent() 860d961a43aSToby Isaac @*/ 861d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeChildren(DM dm, PetscInt point, PetscInt *numChildren, const PetscInt *children[]) 862d961a43aSToby Isaac { 863d961a43aSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 864d961a43aSToby Isaac PetscSection childSec; 865d961a43aSToby Isaac PetscInt dof = 0; 866d961a43aSToby Isaac PetscErrorCode ierr; 867d961a43aSToby Isaac 868d961a43aSToby Isaac PetscFunctionBegin; 869d961a43aSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 870d961a43aSToby Isaac childSec = mesh->childSection; 871d961a43aSToby Isaac if (childSec && point >= childSec->pStart && point < childSec->pEnd) { 872d961a43aSToby Isaac ierr = PetscSectionGetDof (childSec, point, &dof);CHKERRQ(ierr); 873d961a43aSToby Isaac } 874d961a43aSToby Isaac if (numChildren) *numChildren = dof; 875d961a43aSToby Isaac if (children) { 876d961a43aSToby Isaac if (dof) { 877d961a43aSToby Isaac PetscInt off; 878d961a43aSToby Isaac 879d961a43aSToby Isaac ierr = PetscSectionGetOffset (childSec, point, &off);CHKERRQ(ierr); 880d961a43aSToby Isaac *children = &mesh->children[off]; 881d961a43aSToby Isaac } 882d961a43aSToby Isaac else { 883d961a43aSToby Isaac *children = NULL; 884d961a43aSToby Isaac } 885d961a43aSToby Isaac } 886d961a43aSToby Isaac PetscFunctionReturn(0); 887d961a43aSToby Isaac } 8880c37af3bSToby Isaac 8890c37af3bSToby Isaac #undef __FUNCT__ 8900c37af3bSToby Isaac #define __FUNCT__ "DMPlexComputeConstraintMatrix_ReferenceTree" 8910c37af3bSToby Isaac PetscErrorCode DMPlexComputeConstraintMatrix_ReferenceTree(DM dm) 8920c37af3bSToby Isaac { 8930c37af3bSToby Isaac PetscDS ds; 8940c37af3bSToby Isaac PetscInt spdim; 8950c37af3bSToby Isaac PetscInt numFields, f, c, cStart, cEnd, pStart, pEnd, conStart, conEnd; 8960c37af3bSToby Isaac const PetscInt *anchors; 8970c37af3bSToby Isaac PetscSection section, cSec, aSec; 8980c37af3bSToby Isaac Mat cMat; 8990c37af3bSToby Isaac PetscReal *v0, *v0parent, *vtmp, *J, *Jparent, *invJparent, detJ, detJparent; 9000c37af3bSToby Isaac IS aIS; 9010c37af3bSToby Isaac PetscErrorCode ierr; 9020c37af3bSToby Isaac 9030c37af3bSToby Isaac PetscFunctionBegin; 9040c37af3bSToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 9050c37af3bSToby Isaac ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 9060c37af3bSToby Isaac ierr = PetscDSGetNumFields(ds,&numFields);CHKERRQ(ierr); 9070c37af3bSToby Isaac ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 9080c37af3bSToby Isaac ierr = DMGetDefaultSection(dm,§ion);CHKERRQ(ierr); 9090c37af3bSToby Isaac ierr = DMPlexGetConstraintSection(dm,&cSec);CHKERRQ(ierr); 9100c37af3bSToby Isaac ierr = DMPlexGetConstraints(dm,&aSec,&aIS);CHKERRQ(ierr); 9110c37af3bSToby Isaac ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr); 9120c37af3bSToby Isaac ierr = PetscSectionGetChart(cSec,&conStart,&conEnd);CHKERRQ(ierr); 9130c37af3bSToby Isaac ierr = DMPlexGetConstraintMatrix(dm,&cMat);CHKERRQ(ierr); 9140c37af3bSToby Isaac ierr = DMPlexGetDimension(dm,&spdim);CHKERRQ(ierr); 9150c37af3bSToby Isaac ierr = PetscMalloc6(spdim,&v0,spdim,&v0parent,spdim,&vtmp,spdim*spdim,&J,spdim*spdim,&Jparent,spdim*spdim,&invJparent);CHKERRQ(ierr); 9160c37af3bSToby Isaac 9170c37af3bSToby Isaac for (f = 0; f < numFields; f++) { 9180c37af3bSToby Isaac PetscFE fe; 9190c37af3bSToby Isaac PetscDualSpace space; 9200c37af3bSToby Isaac PetscInt i, j, k, nPoints, offset; 9210c37af3bSToby Isaac PetscInt fSize, fComp; 9220c37af3bSToby Isaac PetscScalar *B = NULL; 9230c37af3bSToby Isaac PetscReal *weights, *pointsRef, *pointsReal; 9240c37af3bSToby Isaac Mat Amat, Bmat, Xmat; 9250c37af3bSToby Isaac 9260c37af3bSToby Isaac ierr = PetscDSGetDiscretization(ds,f,(PetscObject *)(&fe));CHKERRQ(ierr); 9270c37af3bSToby Isaac ierr = PetscFEGetDualSpace(fe,&space);CHKERRQ(ierr); 9280c37af3bSToby Isaac ierr = PetscDualSpaceGetDimension(space,&fSize);CHKERRQ(ierr); 9290c37af3bSToby Isaac ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr); 9300c37af3bSToby Isaac ierr = MatCreate(PETSC_COMM_SELF,&Amat);CHKERRQ(ierr); 9310c37af3bSToby Isaac ierr = MatSetSizes(Amat,fSize,fSize,fSize,fSize);CHKERRQ(ierr); 9320c37af3bSToby Isaac ierr = MatSetType(Amat,MATSEQDENSE);CHKERRQ(ierr); 9330c37af3bSToby Isaac ierr = MatSetUp(Amat);CHKERRQ(ierr); 9340c37af3bSToby Isaac ierr = MatDuplicate(Amat,MAT_DO_NOT_COPY_VALUES,&Bmat);CHKERRQ(ierr); 9350c37af3bSToby Isaac ierr = MatDuplicate(Amat,MAT_DO_NOT_COPY_VALUES,&Xmat);CHKERRQ(ierr); 9360c37af3bSToby Isaac nPoints = 0; 9370c37af3bSToby Isaac for (i = 0; i < fSize; i++) { 9380c37af3bSToby Isaac PetscInt qPoints; 9390c37af3bSToby Isaac PetscQuadrature quad; 9400c37af3bSToby Isaac 9410c37af3bSToby Isaac ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 9420c37af3bSToby Isaac ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr); 9430c37af3bSToby Isaac nPoints += qPoints; 9440c37af3bSToby Isaac } 9450c37af3bSToby Isaac ierr = PetscMalloc3(nPoints,&weights,spdim*nPoints,&pointsRef,spdim*nPoints,&pointsReal);CHKERRQ(ierr); 9460c37af3bSToby Isaac offset = 0; 9470c37af3bSToby Isaac for (i = 0; i < fSize; i++) { 9480c37af3bSToby Isaac PetscInt qPoints; 9490c37af3bSToby Isaac const PetscReal *p, *w; 9500c37af3bSToby Isaac PetscQuadrature quad; 9510c37af3bSToby Isaac 9520c37af3bSToby Isaac ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 9530c37af3bSToby Isaac ierr = PetscQuadratureGetData(quad,NULL,&qPoints,&p,&w);CHKERRQ(ierr); 9540c37af3bSToby Isaac ierr = PetscMemcpy(weights+offset,w,qPoints*sizeof(*w));CHKERRQ(ierr); 9550c37af3bSToby Isaac ierr = PetscMemcpy(pointsRef+spdim*offset,p,spdim*qPoints*sizeof(*p));CHKERRQ(ierr); 9560c37af3bSToby Isaac offset += qPoints; 9570c37af3bSToby Isaac } 9580c37af3bSToby Isaac ierr = PetscFEGetTabulation(fe,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr); 9590c37af3bSToby Isaac offset = 0; 9600c37af3bSToby Isaac for (i = 0; i < fSize; i++) { 9610c37af3bSToby Isaac PetscInt qPoints; 9620c37af3bSToby Isaac PetscQuadrature quad; 9630c37af3bSToby Isaac 9640c37af3bSToby Isaac ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 9650c37af3bSToby Isaac ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr); 9660c37af3bSToby Isaac for (j = 0; j < fSize; j++) { 9670c37af3bSToby Isaac PetscScalar val = 0.; 9680c37af3bSToby Isaac 9690c37af3bSToby Isaac for (k = 0; k < qPoints; k++) { 9700c37af3bSToby Isaac val += B[((offset + k) * fSize + j) * fComp] * weights[k]; 9710c37af3bSToby Isaac } 9720c37af3bSToby Isaac ierr = MatSetValue(Amat,i,j,val,INSERT_VALUES);CHKERRQ(ierr); 9730c37af3bSToby Isaac } 9740c37af3bSToby Isaac offset += qPoints; 9750c37af3bSToby Isaac } 9760c37af3bSToby Isaac ierr = PetscFERestoreTabulation(fe,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr); 9770c37af3bSToby Isaac ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9780c37af3bSToby Isaac ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9790c37af3bSToby Isaac ierr = MatLUFactor(Amat,NULL,NULL,NULL);CHKERRQ(ierr); 9800c37af3bSToby Isaac for (c = cStart; c < cEnd; c++) { 9810c37af3bSToby Isaac PetscInt parent; 9820c37af3bSToby Isaac PetscInt closureSize, closureSizeP, *closure = NULL, *closureP = NULL; 9830c37af3bSToby Isaac PetscInt *childOffsets, *parentOffsets; 9840c37af3bSToby Isaac 9850c37af3bSToby Isaac ierr = DMPlexGetTreeParent(dm,c,&parent,NULL);CHKERRQ(ierr); 9860c37af3bSToby Isaac if (parent == c) continue; 9870c37af3bSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 9880c37af3bSToby Isaac for (i = 0; i < closureSize; i++) { 9890c37af3bSToby Isaac PetscInt p = closure[2*i]; 9900c37af3bSToby Isaac PetscInt conDof; 9910c37af3bSToby Isaac 9920c37af3bSToby Isaac if (p < conStart || p >= conEnd) continue; 9930c37af3bSToby Isaac if (numFields > 1) { 9940c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(cSec,p,f,&conDof);CHKERRQ(ierr); 9950c37af3bSToby Isaac } 9960c37af3bSToby Isaac else { 9970c37af3bSToby Isaac ierr = PetscSectionGetDof(cSec,p,&conDof);CHKERRQ(ierr); 9980c37af3bSToby Isaac } 9990c37af3bSToby Isaac if (conDof) break; 10000c37af3bSToby Isaac } 10010c37af3bSToby Isaac if (i == closureSize) { 10020c37af3bSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 10030c37af3bSToby Isaac continue; 10040c37af3bSToby Isaac } 10050c37af3bSToby Isaac 10060c37af3bSToby Isaac ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr); 10070c37af3bSToby Isaac ierr = DMPlexComputeCellGeometry(dm, parent, v0parent, Jparent, invJparent, &detJparent);CHKERRQ(ierr); 10080c37af3bSToby Isaac for (i = 0; i < nPoints; i++) { 10090c37af3bSToby Isaac CoordinatesRefToReal(spdim, spdim, v0, J, &pointsRef[i*spdim],vtmp);CHKERRQ(ierr); 10100c37af3bSToby Isaac CoordinatesRealToRef(spdim, spdim, v0parent, invJparent, vtmp, &pointsReal[i*spdim]);CHKERRQ(ierr); 10110c37af3bSToby Isaac } 10120c37af3bSToby Isaac ierr = PetscFEGetTabulation(fe,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr); 10130c37af3bSToby Isaac offset = 0; 10140c37af3bSToby Isaac for (i = 0; i < fSize; i++) { 10150c37af3bSToby Isaac PetscInt qPoints; 10160c37af3bSToby Isaac PetscQuadrature quad; 10170c37af3bSToby Isaac 10180c37af3bSToby Isaac ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 10190c37af3bSToby Isaac ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr); 10200c37af3bSToby Isaac for (j = 0; j < fSize; j++) { 10210c37af3bSToby Isaac PetscScalar val = 0.; 10220c37af3bSToby Isaac 10230c37af3bSToby Isaac for (k = 0; k < qPoints; k++) { 10240c37af3bSToby Isaac val += B[((offset + k) * fSize + j) * fComp] * weights[k]; 10250c37af3bSToby Isaac } 10260c37af3bSToby Isaac MatSetValue(Bmat,i,j,val,INSERT_VALUES);CHKERRQ(ierr); 10270c37af3bSToby Isaac } 10280c37af3bSToby Isaac offset += qPoints; 10290c37af3bSToby Isaac } 10300c37af3bSToby Isaac ierr = PetscFERestoreTabulation(fe,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr); 10310c37af3bSToby Isaac ierr = MatAssemblyBegin(Bmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10320c37af3bSToby Isaac ierr = MatAssemblyEnd(Bmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10330c37af3bSToby Isaac ierr = MatMatSolve(Amat,Bmat,Xmat);CHKERRQ(ierr); 10340c37af3bSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSizeP,&closureP);CHKERRQ(ierr); 10350c37af3bSToby Isaac ierr = PetscMalloc2(closureSize+1,&childOffsets,closureSizeP+1,&parentOffsets);CHKERRQ(ierr); 10360c37af3bSToby Isaac childOffsets[0] = 0; 10370c37af3bSToby Isaac for (i = 0; i < closureSize; i++) { 10380c37af3bSToby Isaac PetscInt p = closure[2*i]; 10390c37af3bSToby Isaac PetscInt dof; 10400c37af3bSToby Isaac 10410c37af3bSToby Isaac if (numFields > 1) { 10420c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(section,p,f,&dof);CHKERRQ(ierr); 10430c37af3bSToby Isaac } 10440c37af3bSToby Isaac else { 10450c37af3bSToby Isaac ierr = PetscSectionGetDof(section,p,&dof);CHKERRQ(ierr); 10460c37af3bSToby Isaac } 10470c37af3bSToby Isaac childOffsets[i+1]=childOffsets[i]+dof / fComp; 10480c37af3bSToby Isaac } 10490c37af3bSToby Isaac parentOffsets[0] = 0; 10500c37af3bSToby Isaac for (i = 0; i < closureSizeP; i++) { 10510c37af3bSToby Isaac PetscInt p = closureP[2*i]; 10520c37af3bSToby Isaac PetscInt dof; 10530c37af3bSToby Isaac 10540c37af3bSToby Isaac if (numFields > 1) { 10550c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(section,p,f,&dof);CHKERRQ(ierr); 10560c37af3bSToby Isaac } 10570c37af3bSToby Isaac else { 10580c37af3bSToby Isaac ierr = PetscSectionGetDof(section,p,&dof);CHKERRQ(ierr); 10590c37af3bSToby Isaac } 10600c37af3bSToby Isaac parentOffsets[i+1]=parentOffsets[i]+dof / fComp; 10610c37af3bSToby Isaac } 10620c37af3bSToby Isaac for (i = 0; i < closureSize; i++) { 10630c37af3bSToby Isaac PetscInt conDof, conOff, aDof, aOff; 10640c37af3bSToby Isaac PetscInt p = closure[2*i]; 10650c37af3bSToby Isaac PetscInt o = closure[2*i+1]; 10660c37af3bSToby Isaac 10670c37af3bSToby Isaac if (p < conStart || p >= conEnd) continue; 10680c37af3bSToby Isaac if (numFields > 1) { 10690c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(cSec,p,f,&conDof);CHKERRQ(ierr); 10700c37af3bSToby Isaac ierr = PetscSectionGetFieldOffset(cSec,p,f,&conOff);CHKERRQ(ierr); 10710c37af3bSToby Isaac } 10720c37af3bSToby Isaac else { 10730c37af3bSToby Isaac ierr = PetscSectionGetDof(cSec,p,&conDof);CHKERRQ(ierr); 10740c37af3bSToby Isaac ierr = PetscSectionGetOffset(cSec,p,&conOff);CHKERRQ(ierr); 10750c37af3bSToby Isaac } 10760c37af3bSToby Isaac if (!conDof) continue; 10770c37af3bSToby Isaac ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr); 10780c37af3bSToby Isaac ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr); 10790c37af3bSToby Isaac for (k = 0; k < aDof; k++) { 10800c37af3bSToby Isaac PetscInt a = anchors[aOff + k]; 10810c37af3bSToby Isaac PetscInt aSecDof, aSecOff; 10820c37af3bSToby Isaac 10830c37af3bSToby Isaac if (numFields > 1) { 10840c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(section,a,f,&aSecDof);CHKERRQ(ierr); 10850c37af3bSToby Isaac ierr = PetscSectionGetFieldOffset(section,a,f,&aSecOff);CHKERRQ(ierr); 10860c37af3bSToby Isaac } 10870c37af3bSToby Isaac else { 10880c37af3bSToby Isaac ierr = PetscSectionGetDof(section,a,&aSecDof);CHKERRQ(ierr); 10890c37af3bSToby Isaac ierr = PetscSectionGetOffset(section,a,&aSecOff);CHKERRQ(ierr); 10900c37af3bSToby Isaac } 10910c37af3bSToby Isaac if (!aSecDof) continue; 10920c37af3bSToby Isaac 10930c37af3bSToby Isaac for (j = 0; j < closureSizeP; j++) { 10940c37af3bSToby Isaac PetscInt q = closureP[2*j]; 10950c37af3bSToby Isaac PetscInt oq = closureP[2*j+1]; 10960c37af3bSToby Isaac 10970c37af3bSToby Isaac if (q == a) { 10980c37af3bSToby Isaac PetscInt r, s, t; 10990c37af3bSToby Isaac 11000c37af3bSToby Isaac for (r = childOffsets[i]; r < childOffsets[i+1]; r++) { 11010c37af3bSToby Isaac for (s = parentOffsets[j]; s < parentOffsets[j+1]; s++) { 11020c37af3bSToby Isaac PetscScalar val; 11030c37af3bSToby Isaac PetscInt insertCol, insertRow; 11040c37af3bSToby Isaac 11050c37af3bSToby Isaac ierr = MatGetValue(Xmat,r,s,&val);CHKERRQ(ierr); 11060c37af3bSToby Isaac if (o >= 0) { 11070c37af3bSToby Isaac insertRow = conOff + fComp * (r - childOffsets[i]); 11080c37af3bSToby Isaac } 11090c37af3bSToby Isaac else { 11100c37af3bSToby Isaac insertRow = conOff + fComp * (childOffsets[i + 1] - 1 - r); 11110c37af3bSToby Isaac } 11120c37af3bSToby Isaac if (oq >= 0) { 11130c37af3bSToby Isaac insertCol = aSecOff + fComp * (s - parentOffsets[j]); 11140c37af3bSToby Isaac } 11150c37af3bSToby Isaac else { 11160c37af3bSToby Isaac insertCol = aSecOff + fComp * (parentOffsets[j + 1] - 1 - s); 11170c37af3bSToby Isaac } 11180c37af3bSToby Isaac for (t = 0; t < fComp; t++) { 11190c37af3bSToby Isaac ierr = MatSetValue(cMat,insertRow + t,insertCol + t,val,INSERT_VALUES);CHKERRQ(ierr); 11200c37af3bSToby Isaac } 11210c37af3bSToby Isaac } 11220c37af3bSToby Isaac } 11230c37af3bSToby Isaac } 11240c37af3bSToby Isaac } 11250c37af3bSToby Isaac } 11260c37af3bSToby Isaac } 11270c37af3bSToby Isaac ierr = PetscFree2(childOffsets,parentOffsets);CHKERRQ(ierr); 11280c37af3bSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 11290c37af3bSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSizeP,&closureP);CHKERRQ(ierr); 11300c37af3bSToby Isaac } 11310c37af3bSToby Isaac ierr = MatDestroy(&Amat);CHKERRQ(ierr); 11320c37af3bSToby Isaac ierr = MatDestroy(&Bmat);CHKERRQ(ierr); 11330c37af3bSToby Isaac ierr = MatDestroy(&Xmat);CHKERRQ(ierr); 11340c37af3bSToby Isaac ierr = PetscFree3(weights,pointsRef,pointsReal);CHKERRQ(ierr); 11350c37af3bSToby Isaac } 11360c37af3bSToby Isaac ierr = MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11370c37af3bSToby Isaac ierr = MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11380c37af3bSToby Isaac ierr = PetscFree6(v0,v0parent,vtmp,J,Jparent,invJparent);CHKERRQ(ierr); 11390c37af3bSToby Isaac ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr); 11400c37af3bSToby Isaac 11410c37af3bSToby Isaac PetscFunctionReturn(0); 11420c37af3bSToby Isaac } 114395a0b26dSToby Isaac 114495a0b26dSToby Isaac #undef __FUNCT__ 114595a0b26dSToby Isaac #define __FUNCT__ "DMPlexComputeConstraintMatrix_Tree" 114695a0b26dSToby Isaac PetscErrorCode DMPlexComputeConstraintMatrix_Tree(DM dm) 114795a0b26dSToby Isaac { 114895a0b26dSToby Isaac DM refTree; 114995a0b26dSToby Isaac PetscDS ds; 115095a0b26dSToby Isaac Mat refCmat, cMat; 115195a0b26dSToby Isaac PetscInt numFields, f, pRefStart, pRefEnd, p, *rows, *cols, maxDof, maxAnDof, *perm, *iperm, pStart, pEnd, conStart, conEnd, **refPointFieldN; 115295a0b26dSToby Isaac PetscScalar ***refPointFieldMats, *pointWork; 115395a0b26dSToby Isaac PetscSection refConSec, conSec, refAnSec, anSec, section, refSection; 115495a0b26dSToby Isaac IS refAnIS, anIS; 115595a0b26dSToby Isaac const PetscInt *refAnchors, *anchors; 115695a0b26dSToby Isaac PetscErrorCode ierr; 115795a0b26dSToby Isaac 115895a0b26dSToby Isaac PetscFunctionBegin; 115995a0b26dSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 116095a0b26dSToby Isaac ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 116195a0b26dSToby Isaac ierr = PetscDSGetNumFields(ds,&numFields);CHKERRQ(ierr); 116295a0b26dSToby Isaac ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr); 116395a0b26dSToby Isaac ierr = DMSetDS(refTree,ds);CHKERRQ(ierr); 116495a0b26dSToby Isaac ierr = DMPlexComputeConstraintMatrix_ReferenceTree(refTree);CHKERRQ(ierr); 116595a0b26dSToby Isaac ierr = DMPlexGetConstraintMatrix(refTree,&refCmat);CHKERRQ(ierr); 116695a0b26dSToby Isaac ierr = DMPlexGetChart(refTree,&pRefStart,&pRefEnd);CHKERRQ(ierr); 116795a0b26dSToby Isaac ierr = DMPlexGetConstraints(refTree,&refAnSec,&refAnIS);CHKERRQ(ierr); 116895a0b26dSToby Isaac ierr = DMPlexGetConstraints(dm,&anSec,&anIS);CHKERRQ(ierr); 116995a0b26dSToby Isaac ierr = ISGetIndices(refAnIS,&refAnchors);CHKERRQ(ierr); 117095a0b26dSToby Isaac ierr = ISGetIndices(anIS,&anchors);CHKERRQ(ierr); 117195a0b26dSToby Isaac ierr = DMPlexGetConstraintSection(refTree,&refConSec);CHKERRQ(ierr); 117295a0b26dSToby Isaac ierr = DMPlexGetConstraintSection(dm,&conSec);CHKERRQ(ierr); 117395a0b26dSToby Isaac ierr = DMGetDefaultSection(refTree,&refSection);CHKERRQ(ierr); 117495a0b26dSToby Isaac ierr = DMGetDefaultSection(dm,§ion);CHKERRQ(ierr); 117595a0b26dSToby Isaac ierr = PetscSectionGetChart(refConSec,&pRefStart,&pRefEnd);CHKERRQ(ierr); 117695a0b26dSToby Isaac ierr = PetscSectionGetChart(conSec,&conStart,&conEnd);CHKERRQ(ierr); 117795a0b26dSToby Isaac ierr = PetscMalloc1(pRefEnd-pRefStart,&refPointFieldMats);CHKERRQ(ierr); 117895a0b26dSToby Isaac ierr = PetscMalloc1(pRefEnd-pRefStart,&refPointFieldN);CHKERRQ(ierr); 117995a0b26dSToby Isaac ierr = PetscSectionGetMaxDof(refConSec,&maxDof);CHKERRQ(ierr); 118095a0b26dSToby Isaac ierr = PetscSectionGetMaxDof(refAnSec,&maxAnDof);CHKERRQ(ierr); 118195a0b26dSToby Isaac ierr = PetscMalloc1(maxDof,&rows);CHKERRQ(ierr); 118295a0b26dSToby Isaac ierr = PetscMalloc1(maxDof*maxAnDof,&cols);CHKERRQ(ierr); 118395a0b26dSToby Isaac ierr = PetscMalloc1(maxDof*maxDof*maxAnDof,&pointWork);CHKERRQ(ierr); 118495a0b26dSToby Isaac 118595a0b26dSToby Isaac /* step 1: get submats for every constrained point in the reference tree */ 118695a0b26dSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 118795a0b26dSToby Isaac PetscInt parent, closureSize, *closure = NULL, pDof; 118895a0b26dSToby Isaac 118995a0b26dSToby Isaac ierr = DMPlexGetTreeParent(refTree,p,&parent,NULL);CHKERRQ(ierr); 119095a0b26dSToby Isaac ierr = PetscSectionGetDof(refConSec,p,&pDof);CHKERRQ(ierr); 119195a0b26dSToby Isaac if (!pDof || parent == p) continue; 119295a0b26dSToby Isaac 119395a0b26dSToby Isaac ierr = PetscMalloc1(numFields,&refPointFieldMats[p-pRefStart]);CHKERRQ(ierr); 119495a0b26dSToby Isaac ierr = PetscCalloc1(numFields,&refPointFieldN[p-pRefStart]);CHKERRQ(ierr); 119595a0b26dSToby Isaac ierr = DMPlexGetTransitiveClosure(refTree,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 119695a0b26dSToby Isaac for (f = 0; f < numFields; f++) { 119795a0b26dSToby Isaac PetscInt cDof, cOff, numCols, r, i, fComp; 119895a0b26dSToby Isaac PetscFE fe; 119995a0b26dSToby Isaac 120095a0b26dSToby Isaac ierr = PetscDSGetDiscretization(ds,f,(PetscObject *) &fe);CHKERRQ(ierr); 120195a0b26dSToby Isaac ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr); 120295a0b26dSToby Isaac 120395a0b26dSToby Isaac if (numFields > 1) { 120495a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(refConSec,p,f,&cDof);CHKERRQ(ierr); 120595a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(refConSec,p,f,&cOff);CHKERRQ(ierr); 120695a0b26dSToby Isaac } 120795a0b26dSToby Isaac else { 120895a0b26dSToby Isaac ierr = PetscSectionGetDof(refConSec,p,&cDof);CHKERRQ(ierr); 120995a0b26dSToby Isaac ierr = PetscSectionGetOffset(refConSec,p,&cOff);CHKERRQ(ierr); 121095a0b26dSToby Isaac } 121195a0b26dSToby Isaac 121295a0b26dSToby Isaac if (!cDof) continue; 121395a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 121495a0b26dSToby Isaac rows[r] = cOff + r; 121595a0b26dSToby Isaac } 121695a0b26dSToby Isaac numCols = 0; 121795a0b26dSToby Isaac for (i = 0; i < closureSize; i++) { 121895a0b26dSToby Isaac PetscInt q = closure[2*i]; 121995a0b26dSToby Isaac PetscInt o = closure[2*i+1]; 122095a0b26dSToby Isaac PetscInt aDof, aOff, j; 122195a0b26dSToby Isaac 122295a0b26dSToby Isaac if (numFields > 1) { 122395a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(refSection,q,f,&aDof);CHKERRQ(ierr); 122495a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(refSection,q,f,&aOff);CHKERRQ(ierr); 122595a0b26dSToby Isaac } 122695a0b26dSToby Isaac else { 122795a0b26dSToby Isaac ierr = PetscSectionGetDof(refSection,q,&aDof);CHKERRQ(ierr); 122895a0b26dSToby Isaac ierr = PetscSectionGetOffset(refSection,q,&aOff);CHKERRQ(ierr); 122995a0b26dSToby Isaac } 123095a0b26dSToby Isaac 123195a0b26dSToby Isaac for (j = 0; j < aDof; j++) { 123295a0b26dSToby Isaac PetscInt node = (o >= 0) ? (j / fComp) : ((aDof - 1 - j) / fComp); 123395a0b26dSToby Isaac PetscInt comp = (j % fComp); 123495a0b26dSToby Isaac 123595a0b26dSToby Isaac cols[numCols++] = aOff + node * fComp + comp; 123695a0b26dSToby Isaac } 123795a0b26dSToby Isaac } 123895a0b26dSToby Isaac refPointFieldN[p-pRefStart][f] = numCols; 123995a0b26dSToby Isaac ierr = PetscMalloc1(cDof*numCols,&refPointFieldMats[p-pRefStart][f]);CHKERRQ(ierr); 124095a0b26dSToby Isaac ierr = MatGetValues(refCmat,cDof,rows,numCols,cols,refPointFieldMats[p-pRefStart][f]);CHKERRQ(ierr); 124195a0b26dSToby Isaac } 124295a0b26dSToby Isaac ierr = DMPlexRestoreTransitiveClosure(refTree,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 124395a0b26dSToby Isaac } 124495a0b26dSToby Isaac 124595a0b26dSToby Isaac /* step 2: compute the preorder */ 124695a0b26dSToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 124795a0b26dSToby Isaac ierr = PetscMalloc2(pEnd-pStart,&perm,pEnd-pStart,&iperm);CHKERRQ(ierr); 124895a0b26dSToby Isaac for (p = pStart; p < pEnd; p++) { 124995a0b26dSToby Isaac perm[p - pStart] = p; 125095a0b26dSToby Isaac iperm[p - pStart] = p-pStart; 125195a0b26dSToby Isaac } 125295a0b26dSToby Isaac for (p = 0; p < pEnd - pStart;) { 125395a0b26dSToby Isaac PetscInt point = perm[p]; 125495a0b26dSToby Isaac PetscInt parent; 125595a0b26dSToby Isaac 125695a0b26dSToby Isaac ierr = DMPlexGetTreeParent(dm,point,&parent,NULL);CHKERRQ(ierr); 125795a0b26dSToby Isaac if (parent == point) { 125895a0b26dSToby Isaac p++; 125995a0b26dSToby Isaac } 126095a0b26dSToby Isaac else { 126195a0b26dSToby Isaac PetscInt size, closureSize, *closure = NULL, i; 126295a0b26dSToby Isaac 126395a0b26dSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 126495a0b26dSToby Isaac for (i = 0; i < closureSize; i++) { 126595a0b26dSToby Isaac PetscInt q = closure[2*i]; 126695a0b26dSToby Isaac if (iperm[q-pStart] > iperm[point-pStart]) { 126795a0b26dSToby Isaac /* swap */ 126895a0b26dSToby Isaac perm[p] = q; 126995a0b26dSToby Isaac perm[iperm[q-pStart]] = point; 127095a0b26dSToby Isaac iperm[point-pStart] = iperm[q-pStart]; 127195a0b26dSToby Isaac iperm[q-pStart] = p; 127295a0b26dSToby Isaac break; 127395a0b26dSToby Isaac } 127495a0b26dSToby Isaac } 127595a0b26dSToby Isaac size = closureSize; 127695a0b26dSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 127795a0b26dSToby Isaac if (i == size) { 127895a0b26dSToby Isaac p++; 127995a0b26dSToby Isaac } 128095a0b26dSToby Isaac } 128195a0b26dSToby Isaac } 128295a0b26dSToby Isaac 128395a0b26dSToby Isaac /* step 3: fill the constraint matrix */ 128495a0b26dSToby Isaac /* we are going to use a preorder progressive fill strategy. Mat doesn't 128595a0b26dSToby Isaac * allow progressive fill without assembly, so we are going to set up the 128695a0b26dSToby Isaac * values outside of the Mat first. 128795a0b26dSToby Isaac */ 128895a0b26dSToby Isaac { 128995a0b26dSToby Isaac PetscInt nRows, row, nnz; 129095a0b26dSToby Isaac PetscBool done; 129195a0b26dSToby Isaac const PetscInt *ia, *ja; 129295a0b26dSToby Isaac PetscScalar *vals; 129395a0b26dSToby Isaac 129495a0b26dSToby Isaac ierr = DMPlexGetConstraintMatrix(dm,&cMat);CHKERRQ(ierr); 129595a0b26dSToby Isaac ierr = MatGetRowIJ(cMat,0,PETSC_FALSE,PETSC_FALSE,&nRows,&ia,&ja,&done);CHKERRQ(ierr); 129695a0b26dSToby Isaac if (!done) SETERRQ(PetscObjectComm((PetscObject)cMat),PETSC_ERR_PLIB,"Could not get RowIJ of constraint matrix"); 129795a0b26dSToby Isaac nnz = ia[nRows]; 129895a0b26dSToby Isaac /* malloc and then zero rows right before we fill them: this way valgrind 129995a0b26dSToby Isaac * can tell if we are doing progressive fill in the wrong order */ 130095a0b26dSToby Isaac ierr = PetscMalloc1(nnz,&vals);CHKERRQ(ierr); 130195a0b26dSToby Isaac for (p = 0; p < pEnd - pStart; p++) { 130295a0b26dSToby Isaac PetscInt parent, childid, closureSize, *closure = NULL; 130395a0b26dSToby Isaac PetscInt point = perm[p], pointDof; 130495a0b26dSToby Isaac 130595a0b26dSToby Isaac ierr = DMPlexGetTreeParent(dm,point,&parent,&childid);CHKERRQ(ierr); 130695a0b26dSToby Isaac if ((point < conStart) || (point >= conEnd) || (parent == point)) continue; 130795a0b26dSToby Isaac ierr = PetscSectionGetDof(conSec,point,&pointDof);CHKERRQ(ierr); 130895a0b26dSToby Isaac if (!pointDof) continue; 130995a0b26dSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 131095a0b26dSToby Isaac for (f = 0; f < numFields; f++) { 131195a0b26dSToby Isaac PetscInt cDof, cOff, numCols, numFillCols, i, r, fComp, matOffset, offset; 131295a0b26dSToby Isaac PetscScalar *pointMat; 131395a0b26dSToby Isaac PetscFE fe; 131495a0b26dSToby Isaac 131595a0b26dSToby Isaac ierr = PetscDSGetDiscretization(ds,f,(PetscObject *) &fe);CHKERRQ(ierr); 131695a0b26dSToby Isaac ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr); 131795a0b26dSToby Isaac 131895a0b26dSToby Isaac if (numFields > 1) { 131995a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(conSec,point,f,&cDof);CHKERRQ(ierr); 132095a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(conSec,point,f,&cOff);CHKERRQ(ierr); 132195a0b26dSToby Isaac } 132295a0b26dSToby Isaac else { 132395a0b26dSToby Isaac ierr = PetscSectionGetDof(conSec,point,&cDof);CHKERRQ(ierr); 132495a0b26dSToby Isaac ierr = PetscSectionGetOffset(conSec,point,&cOff);CHKERRQ(ierr); 132595a0b26dSToby Isaac } 132695a0b26dSToby Isaac if (!cDof) continue; 132795a0b26dSToby Isaac 132895a0b26dSToby Isaac /* make sure that every row for this point is the same size */ 132995a0b26dSToby Isaac #if defined(PETSC_USE_DEBUG) 133095a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 133195a0b26dSToby Isaac if (cDof > 1 && r) { 133295a0b26dSToby Isaac if ((ia[cOff+r+1]-ia[cOff+r]) != (ia[cOff+r]-ia[cOff+r-1])) { 133395a0b26dSToby Isaac 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])); 133495a0b26dSToby Isaac } 133595a0b26dSToby Isaac } 133695a0b26dSToby Isaac } 133795a0b26dSToby Isaac #endif 133895a0b26dSToby Isaac /* zero rows */ 133995a0b26dSToby Isaac for (i = ia[cOff] ; i< ia[cOff+cDof];i++) { 134095a0b26dSToby Isaac vals[i] = 0.; 134195a0b26dSToby Isaac } 134295a0b26dSToby Isaac matOffset = ia[cOff]; 134395a0b26dSToby Isaac numFillCols = ia[cOff+1] - matOffset; 134495a0b26dSToby Isaac pointMat = refPointFieldMats[childid-pRefStart][f]; 134595a0b26dSToby Isaac numCols = refPointFieldN[childid-pRefStart][f]; 134695a0b26dSToby Isaac offset = 0; 134795a0b26dSToby Isaac for (i = 0; i < closureSize; i++) { 134895a0b26dSToby Isaac PetscInt q = closure[2*i]; 134995a0b26dSToby Isaac PetscInt o = closure[2*i+1]; 135095a0b26dSToby Isaac PetscInt aDof, aOff, j, k, qConDof, qConOff; 135195a0b26dSToby Isaac 135295a0b26dSToby Isaac qConDof = qConOff = 0; 135395a0b26dSToby Isaac if (numFields > 1) { 135495a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(section,q,f,&aDof);CHKERRQ(ierr); 135595a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(section,q,f,&aOff);CHKERRQ(ierr); 135695a0b26dSToby Isaac if (q >= conStart && q < conEnd) { 135795a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(conSec,q,f,&qConDof);CHKERRQ(ierr); 135895a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(conSec,q,f,&qConOff);CHKERRQ(ierr); 135995a0b26dSToby Isaac } 136095a0b26dSToby Isaac } 136195a0b26dSToby Isaac else { 136295a0b26dSToby Isaac ierr = PetscSectionGetDof(section,q,&aDof);CHKERRQ(ierr); 136395a0b26dSToby Isaac ierr = PetscSectionGetOffset(section,q,&aOff);CHKERRQ(ierr); 136495a0b26dSToby Isaac if (q >= conStart && q < conEnd) { 136595a0b26dSToby Isaac ierr = PetscSectionGetDof(conSec,q,&qConDof);CHKERRQ(ierr); 136695a0b26dSToby Isaac ierr = PetscSectionGetOffset(conSec,q,&qConOff);CHKERRQ(ierr); 136795a0b26dSToby Isaac } 136895a0b26dSToby Isaac } 136995a0b26dSToby Isaac if (!aDof) continue; 137095a0b26dSToby Isaac if (qConDof) { 137195a0b26dSToby Isaac /* this point has anchors: its rows of the matrix should already 137295a0b26dSToby Isaac * be filled, thanks to preordering */ 137395a0b26dSToby Isaac /* first multiply into pointWork, then set in matrix */ 137495a0b26dSToby Isaac PetscInt aMatOffset = ia[qConOff]; 137595a0b26dSToby Isaac PetscInt aNumFillCols = ia[qConOff + 1] - aMatOffset; 137695a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 137795a0b26dSToby Isaac for (j = 0; j < aNumFillCols; j++) { 137895a0b26dSToby Isaac PetscScalar inVal = 0; 137995a0b26dSToby Isaac for (k = 0; k < aDof; k++) { 138095a0b26dSToby Isaac PetscInt node = (o >= 0) ? (k / fComp) : ((aDof - 1 - k) / fComp); 138195a0b26dSToby Isaac PetscInt comp = (k % fComp); 138295a0b26dSToby Isaac PetscInt col = node * fComp + comp; 138395a0b26dSToby Isaac 138495a0b26dSToby Isaac inVal += pointMat[r * numCols + offset + col] * vals[aMatOffset + aNumFillCols * k + j]; 138595a0b26dSToby Isaac } 138695a0b26dSToby Isaac pointWork[r * aNumFillCols + j] = inVal; 138795a0b26dSToby Isaac } 138895a0b26dSToby Isaac } 138995a0b26dSToby Isaac /* assume that the columns are sorted, spend less time searching */ 139095a0b26dSToby Isaac for (j = 0, k = 0; j < aNumFillCols; j++) { 139195a0b26dSToby Isaac PetscInt col = ja[aMatOffset + j]; 139295a0b26dSToby Isaac for (;k < numFillCols; k++) { 139395a0b26dSToby Isaac if (ja[matOffset + k] == col) { 139495a0b26dSToby Isaac break; 139595a0b26dSToby Isaac } 139695a0b26dSToby Isaac } 139795a0b26dSToby Isaac if (k == numFillCols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"No nonzero space for (%d, %d)", cOff, col); 139895a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 139995a0b26dSToby Isaac vals[matOffset + numFillCols * r + k] = pointWork[r * aNumFillCols + j]; 140095a0b26dSToby Isaac } 140195a0b26dSToby Isaac } 140295a0b26dSToby Isaac } 140395a0b26dSToby Isaac else { 140495a0b26dSToby Isaac /* find where to put this portion of pointMat into the matrix */ 140595a0b26dSToby Isaac for (k = 0; k < numFillCols; k++) { 140695a0b26dSToby Isaac if (ja[matOffset + k] == aOff) { 140795a0b26dSToby Isaac break; 140895a0b26dSToby Isaac } 140995a0b26dSToby Isaac } 141095a0b26dSToby Isaac if (k == numFillCols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"No nonzero space for (%d, %d)", cOff, aOff); 141195a0b26dSToby Isaac for (j = 0; j < aDof; j++) { 141295a0b26dSToby Isaac PetscInt node = (o >= 0) ? (j / fComp) : ((aDof - 1 - j) / fComp); 141395a0b26dSToby Isaac PetscInt comp = (j % fComp); 141495a0b26dSToby Isaac PetscInt col = node * fComp + comp; 141595a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 141695a0b26dSToby Isaac vals[matOffset + numFillCols * r + k + col] += pointMat[r * numCols + offset + col]; 141795a0b26dSToby Isaac } 141895a0b26dSToby Isaac } 141995a0b26dSToby Isaac } 142095a0b26dSToby Isaac offset += aDof; 142195a0b26dSToby Isaac } 142295a0b26dSToby Isaac } 142395a0b26dSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 142495a0b26dSToby Isaac } 142595a0b26dSToby Isaac for (row = 0; row < nRows; row++) { 142695a0b26dSToby Isaac ierr = MatSetValues(cMat,1,&row,ia[row+1]-ia[row],&ja[ia[row]],&vals[ia[row]],INSERT_VALUES);CHKERRQ(ierr); 142795a0b26dSToby Isaac } 142895a0b26dSToby Isaac ierr = MatRestoreRowIJ(cMat,0,PETSC_FALSE,PETSC_FALSE,&nRows,&ia,&ja,&done);CHKERRQ(ierr); 142995a0b26dSToby Isaac if (!done) SETERRQ(PetscObjectComm((PetscObject)cMat),PETSC_ERR_PLIB,"Could not restore RowIJ of constraint matrix"); 143095a0b26dSToby Isaac ierr = MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 143195a0b26dSToby Isaac ierr = MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 143295a0b26dSToby Isaac ierr = PetscFree(vals);CHKERRQ(ierr); 143395a0b26dSToby Isaac } 143495a0b26dSToby Isaac 143595a0b26dSToby Isaac /* clean up */ 143695a0b26dSToby Isaac ierr = ISRestoreIndices(refAnIS,&refAnchors);CHKERRQ(ierr); 143795a0b26dSToby Isaac ierr = ISRestoreIndices(anIS,&anchors);CHKERRQ(ierr); 143895a0b26dSToby Isaac ierr = PetscFree2(perm,iperm);CHKERRQ(ierr); 143995a0b26dSToby Isaac ierr = PetscFree(rows);CHKERRQ(ierr); 144095a0b26dSToby Isaac ierr = PetscFree(cols);CHKERRQ(ierr); 144195a0b26dSToby Isaac ierr = PetscFree(pointWork);CHKERRQ(ierr); 144295a0b26dSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 144395a0b26dSToby Isaac PetscInt parent, pDof; 144495a0b26dSToby Isaac 144595a0b26dSToby Isaac ierr = DMPlexGetTreeParent(refTree,p,&parent,NULL);CHKERRQ(ierr); 144695a0b26dSToby Isaac ierr = PetscSectionGetDof(refConSec,p,&pDof);CHKERRQ(ierr); 144795a0b26dSToby Isaac if (!pDof || parent == p) continue; 144895a0b26dSToby Isaac 144995a0b26dSToby Isaac for (f = 0; f < numFields; f++) { 145095a0b26dSToby Isaac PetscInt cDof; 145195a0b26dSToby Isaac 145295a0b26dSToby Isaac if (numFields > 1) { 145395a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(refConSec,p,f,&cDof);CHKERRQ(ierr); 145495a0b26dSToby Isaac } 145595a0b26dSToby Isaac else { 145695a0b26dSToby Isaac ierr = PetscSectionGetDof(refConSec,p,&cDof);CHKERRQ(ierr); 145795a0b26dSToby Isaac } 145895a0b26dSToby Isaac 145995a0b26dSToby Isaac if (!cDof) continue; 146095a0b26dSToby Isaac ierr = PetscFree(refPointFieldMats[p - pRefStart][f]);CHKERRQ(ierr); 146195a0b26dSToby Isaac } 146295a0b26dSToby Isaac ierr = PetscFree(refPointFieldMats[p - pRefStart]);CHKERRQ(ierr); 146395a0b26dSToby Isaac ierr = PetscFree(refPointFieldN[p - pRefStart]);CHKERRQ(ierr); 146495a0b26dSToby Isaac } 146595a0b26dSToby Isaac ierr = PetscFree(refPointFieldMats);CHKERRQ(ierr); 146695a0b26dSToby Isaac ierr = PetscFree(refPointFieldN);CHKERRQ(ierr); 146795a0b26dSToby Isaac PetscFunctionReturn(0); 146895a0b26dSToby Isaac } 146995a0b26dSToby Isaac 1470