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 67da43764aSToby Isaac #undef __FUNCT__ 68da43764aSToby Isaac #define __FUNCT__ "DMPlexCreateDefaultReferenceTree" 69da43764aSToby Isaac /*@ 70da43764aSToby Isaac DMPlexCreateDefaultReferenceTree - create a reference tree for isotropic hierarchical mesh refinement. 71da43764aSToby Isaac 72da43764aSToby Isaac Collective on comm 73da43764aSToby Isaac 74da43764aSToby Isaac Input Parameters: 75da43764aSToby Isaac + comm - the MPI communicator 76da43764aSToby Isaac . dim - the spatial dimension 77da43764aSToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell 78da43764aSToby Isaac 79da43764aSToby Isaac Output Parameters: 80da43764aSToby Isaac . ref - the reference tree DMPlex object 81da43764aSToby Isaac 820b7167a0SToby Isaac Level: intermediate 83da43764aSToby Isaac 84da43764aSToby Isaac .keywords: reference cell 85da43764aSToby Isaac .seealso: DMPlexSetReferenceTree(), DMPlexGetReferenceTree() 86da43764aSToby Isaac @*/ 87da43764aSToby Isaac PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *ref) 88da43764aSToby Isaac { 89da43764aSToby Isaac DM K, Kref; 9010f7e118SToby Isaac PetscInt p, pStart, pEnd, pRefStart, pRefEnd, d, offset, parentSize, *parents, *childIDs; 91da43764aSToby Isaac PetscInt *permvals, *unionCones, *coneSizes, *unionOrientations, numUnionPoints, *numDimPoints, numCones, numVerts; 92da43764aSToby Isaac DMLabel identity, identityRef; 9310f7e118SToby Isaac PetscSection unionSection, unionConeSection, parentSection; 94da43764aSToby Isaac PetscScalar *unionCoords; 95da43764aSToby Isaac IS perm; 96da43764aSToby Isaac PetscErrorCode ierr; 97da43764aSToby Isaac 98da43764aSToby Isaac PetscFunctionBegin; 99da43764aSToby Isaac /* create a reference element */ 100da43764aSToby Isaac ierr = DMPlexCreateReferenceCell(comm, dim, simplex, &K);CHKERRQ(ierr); 101da43764aSToby Isaac ierr = DMPlexCreateLabel(K, "identity");CHKERRQ(ierr); 102da43764aSToby Isaac ierr = DMPlexGetLabel(K, "identity", &identity);CHKERRQ(ierr); 103da43764aSToby Isaac ierr = DMPlexGetChart(K, &pStart, &pEnd);CHKERRQ(ierr); 104da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 105da43764aSToby Isaac ierr = DMLabelSetValue(identity, p, p);CHKERRQ(ierr); 106da43764aSToby Isaac } 107da43764aSToby Isaac /* refine it */ 108da43764aSToby Isaac ierr = DMRefine(K,comm,&Kref);CHKERRQ(ierr); 109da43764aSToby Isaac 110da43764aSToby Isaac /* the reference tree is the union of these two, without duplicating 111da43764aSToby Isaac * points that appear in both */ 112da43764aSToby Isaac ierr = DMPlexGetLabel(Kref, "identity", &identityRef);CHKERRQ(ierr); 113da43764aSToby Isaac ierr = DMPlexGetChart(Kref, &pRefStart, &pRefEnd);CHKERRQ(ierr); 114da43764aSToby Isaac ierr = PetscSectionCreate(comm, &unionSection);CHKERRQ(ierr); 115da43764aSToby Isaac ierr = PetscSectionSetChart(unionSection, 0, (pEnd - pStart) + (pRefEnd - pRefStart));CHKERRQ(ierr); 116da43764aSToby Isaac /* count points that will go in the union */ 117da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 118da43764aSToby Isaac ierr = PetscSectionSetDof(unionSection, p - pStart, 1);CHKERRQ(ierr); 119da43764aSToby Isaac } 120da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 121da43764aSToby Isaac PetscInt q, qSize; 122da43764aSToby Isaac ierr = DMLabelGetValue(identityRef, p, &q);CHKERRQ(ierr); 123da43764aSToby Isaac ierr = DMLabelGetStratumSize(identityRef, q, &qSize);CHKERRQ(ierr); 124da43764aSToby Isaac if (qSize > 1) { 125da43764aSToby Isaac ierr = PetscSectionSetDof(unionSection, p - pRefStart + (pEnd - pStart), 1);CHKERRQ(ierr); 126da43764aSToby Isaac } 127da43764aSToby Isaac } 128da43764aSToby Isaac ierr = PetscMalloc1((pEnd - pStart) + (pRefEnd - pRefStart),&permvals);CHKERRQ(ierr); 129da43764aSToby Isaac offset = 0; 130da43764aSToby Isaac /* stratify points in the union by topological dimension */ 131da43764aSToby Isaac for (d = 0; d <= dim; d++) { 132da43764aSToby Isaac PetscInt cStart, cEnd, c; 133da43764aSToby Isaac 134da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K, d, &cStart, &cEnd);CHKERRQ(ierr); 135da43764aSToby Isaac for (c = cStart; c < cEnd; c++) { 136da43764aSToby Isaac permvals[offset++] = c; 137da43764aSToby Isaac } 138da43764aSToby Isaac 139da43764aSToby Isaac ierr = DMPlexGetHeightStratum(Kref, d, &cStart, &cEnd);CHKERRQ(ierr); 140da43764aSToby Isaac for (c = cStart; c < cEnd; c++) { 141da43764aSToby Isaac permvals[offset++] = c + (pEnd - pStart); 142da43764aSToby Isaac } 143da43764aSToby Isaac } 144da43764aSToby Isaac ierr = ISCreateGeneral(comm, (pEnd - pStart) + (pRefEnd - pRefStart), permvals, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr); 145da43764aSToby Isaac ierr = PetscSectionSetPermutation(unionSection,perm);CHKERRQ(ierr); 146da43764aSToby Isaac ierr = PetscSectionSetUp(unionSection);CHKERRQ(ierr); 147da43764aSToby Isaac ierr = PetscSectionGetStorageSize(unionSection,&numUnionPoints);CHKERRQ(ierr); 148da43764aSToby Isaac ierr = PetscMalloc2(numUnionPoints,&coneSizes,dim+1,&numDimPoints);CHKERRQ(ierr); 149da43764aSToby Isaac /* count dimension points */ 150da43764aSToby Isaac for (d = 0; d <= dim; d++) { 151da43764aSToby Isaac PetscInt cStart, cOff, cOff2; 152da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K,d,&cStart,NULL);CHKERRQ(ierr); 153da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff);CHKERRQ(ierr); 154da43764aSToby Isaac if (d < dim) { 155da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K,d+1,&cStart,NULL);CHKERRQ(ierr); 156da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff2);CHKERRQ(ierr); 157da43764aSToby Isaac } 158da43764aSToby Isaac else { 159da43764aSToby Isaac cOff2 = numUnionPoints; 160da43764aSToby Isaac } 161da43764aSToby Isaac numDimPoints[dim - d] = cOff2 - cOff; 162da43764aSToby Isaac } 163da43764aSToby Isaac ierr = PetscSectionCreate(comm, &unionConeSection);CHKERRQ(ierr); 164da43764aSToby Isaac ierr = PetscSectionSetChart(unionConeSection, 0, numUnionPoints);CHKERRQ(ierr); 165da43764aSToby Isaac /* count the cones in the union */ 166da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 167da43764aSToby Isaac PetscInt dof, uOff; 168da43764aSToby Isaac 169da43764aSToby Isaac ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr); 170da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr); 171da43764aSToby Isaac ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr); 172da43764aSToby Isaac coneSizes[uOff] = dof; 173da43764aSToby Isaac } 174da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 175da43764aSToby Isaac PetscInt dof, uDof, uOff; 176da43764aSToby Isaac 177da43764aSToby Isaac ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr); 178da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 179da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 180da43764aSToby Isaac if (uDof) { 181da43764aSToby Isaac ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr); 182da43764aSToby Isaac coneSizes[uOff] = dof; 183da43764aSToby Isaac } 184da43764aSToby Isaac } 185da43764aSToby Isaac ierr = PetscSectionSetUp(unionConeSection);CHKERRQ(ierr); 186da43764aSToby Isaac ierr = PetscSectionGetStorageSize(unionConeSection,&numCones);CHKERRQ(ierr); 187da43764aSToby Isaac ierr = PetscMalloc2(numCones,&unionCones,numCones,&unionOrientations);CHKERRQ(ierr); 188da43764aSToby Isaac /* write the cones in the union */ 189da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 190da43764aSToby Isaac PetscInt dof, uOff, c, cOff; 191da43764aSToby Isaac const PetscInt *cone, *orientation; 192da43764aSToby Isaac 193da43764aSToby Isaac ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr); 194da43764aSToby Isaac ierr = DMPlexGetCone(K, p, &cone);CHKERRQ(ierr); 195da43764aSToby Isaac ierr = DMPlexGetConeOrientation(K, p, &orientation);CHKERRQ(ierr); 196da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr); 197da43764aSToby Isaac ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr); 198da43764aSToby Isaac for (c = 0; c < dof; c++) { 199da43764aSToby Isaac PetscInt e, eOff; 200da43764aSToby Isaac e = cone[c]; 201da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr); 202da43764aSToby Isaac unionCones[cOff + c] = eOff; 203da43764aSToby Isaac unionOrientations[cOff + c] = orientation[c]; 204da43764aSToby Isaac } 205da43764aSToby Isaac } 206da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 207da43764aSToby Isaac PetscInt dof, uDof, uOff, c, cOff; 208da43764aSToby Isaac const PetscInt *cone, *orientation; 209da43764aSToby Isaac 210da43764aSToby Isaac ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr); 211da43764aSToby Isaac ierr = DMPlexGetCone(Kref, p, &cone);CHKERRQ(ierr); 212da43764aSToby Isaac ierr = DMPlexGetConeOrientation(Kref, p, &orientation);CHKERRQ(ierr); 213da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 214da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 215da43764aSToby Isaac if (uDof) { 216da43764aSToby Isaac ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr); 217da43764aSToby Isaac for (c = 0; c < dof; c++) { 218da43764aSToby Isaac PetscInt e, eOff, eDof; 219da43764aSToby Isaac 220da43764aSToby Isaac e = cone[c]; 221da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, e - pRefStart + (pEnd - pStart),&eDof);CHKERRQ(ierr); 222da43764aSToby Isaac if (eDof) { 223da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pRefStart + (pEnd - pStart), &eOff);CHKERRQ(ierr); 224da43764aSToby Isaac } 225da43764aSToby Isaac else { 226da43764aSToby Isaac ierr = DMLabelGetValue(identityRef, e, &e);CHKERRQ(ierr); 227da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr); 228da43764aSToby Isaac } 229da43764aSToby Isaac unionCones[cOff + c] = eOff; 230da43764aSToby Isaac unionOrientations[cOff + c] = orientation[c]; 231da43764aSToby Isaac } 232da43764aSToby Isaac } 233da43764aSToby Isaac } 234da43764aSToby Isaac /* get the coordinates */ 235da43764aSToby Isaac { 236da43764aSToby Isaac PetscInt vStart, vEnd, vRefStart, vRefEnd, v, vDof, vOff; 237da43764aSToby Isaac PetscSection KcoordsSec, KrefCoordsSec; 238da43764aSToby Isaac Vec KcoordsVec, KrefCoordsVec; 239da43764aSToby Isaac PetscScalar *Kcoords; 240da43764aSToby Isaac 241da43764aSToby Isaac DMGetCoordinateSection(K, &KcoordsSec);CHKERRQ(ierr); 242da43764aSToby Isaac DMGetCoordinatesLocal(K, &KcoordsVec);CHKERRQ(ierr); 243da43764aSToby Isaac DMGetCoordinateSection(Kref, &KrefCoordsSec);CHKERRQ(ierr); 244da43764aSToby Isaac DMGetCoordinatesLocal(Kref, &KrefCoordsVec);CHKERRQ(ierr); 245da43764aSToby Isaac 246da43764aSToby Isaac numVerts = numDimPoints[0]; 247da43764aSToby Isaac ierr = PetscMalloc1(numVerts * dim,&unionCoords);CHKERRQ(ierr); 248da43764aSToby Isaac ierr = DMPlexGetDepthStratum(K,0,&vStart,&vEnd);CHKERRQ(ierr); 249da43764aSToby Isaac 250da43764aSToby Isaac offset = 0; 251da43764aSToby Isaac for (v = vStart; v < vEnd; v++) { 252da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,v - pStart,&vOff);CHKERRQ(ierr); 253da43764aSToby Isaac ierr = VecGetValuesSection(KcoordsVec, KcoordsSec, v, &Kcoords);CHKERRQ(ierr); 254da43764aSToby Isaac for (d = 0; d < dim; d++) { 255da43764aSToby Isaac unionCoords[offset * dim + d] = Kcoords[d]; 256da43764aSToby Isaac } 257da43764aSToby Isaac offset++; 258da43764aSToby Isaac } 259da43764aSToby Isaac ierr = DMPlexGetDepthStratum(Kref,0,&vRefStart,&vRefEnd);CHKERRQ(ierr); 260da43764aSToby Isaac for (v = vRefStart; v < vRefEnd; v++) { 261da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection,v - pRefStart + (pEnd - pStart),&vDof);CHKERRQ(ierr); 262da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,v - pRefStart + (pEnd - pStart),&vOff);CHKERRQ(ierr); 263da43764aSToby Isaac ierr = VecGetValuesSection(KrefCoordsVec, KrefCoordsSec, v, &Kcoords);CHKERRQ(ierr); 264da43764aSToby Isaac if (vDof) { 265da43764aSToby Isaac for (d = 0; d < dim; d++) { 266da43764aSToby Isaac unionCoords[offset * dim + d] = Kcoords[d]; 267da43764aSToby Isaac } 268da43764aSToby Isaac offset++; 269da43764aSToby Isaac } 270da43764aSToby Isaac } 271da43764aSToby Isaac } 272da43764aSToby Isaac ierr = DMCreate(comm,ref);CHKERRQ(ierr); 273da43764aSToby Isaac ierr = DMSetType(*ref,DMPLEX);CHKERRQ(ierr); 274da43764aSToby Isaac ierr = DMPlexSetDimension(*ref,dim);CHKERRQ(ierr); 275da43764aSToby Isaac ierr = DMPlexCreateFromDAG(*ref,dim,numDimPoints,coneSizes,unionCones,unionOrientations,unionCoords);CHKERRQ(ierr); 27610f7e118SToby Isaac /* set the tree */ 27710f7e118SToby Isaac ierr = PetscSectionCreate(comm,&parentSection);CHKERRQ(ierr); 27810f7e118SToby Isaac ierr = PetscSectionSetChart(parentSection,0,numUnionPoints);CHKERRQ(ierr); 27910f7e118SToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 28010f7e118SToby Isaac PetscInt uDof, uOff; 28110f7e118SToby Isaac 28210f7e118SToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 28310f7e118SToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 28410f7e118SToby Isaac if (uDof) { 28510f7e118SToby Isaac PetscSectionSetDof(parentSection,uOff,1);CHKERRQ(ierr); 28610f7e118SToby Isaac } 28710f7e118SToby Isaac } 28810f7e118SToby Isaac ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr); 28910f7e118SToby Isaac ierr = PetscSectionGetStorageSize(parentSection,&parentSize);CHKERRQ(ierr); 29010f7e118SToby Isaac ierr = PetscMalloc2(parentSize,&parents,parentSize,&childIDs);CHKERRQ(ierr); 29110f7e118SToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 29210f7e118SToby Isaac PetscInt uDof, uOff; 29310f7e118SToby Isaac 29410f7e118SToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 29510f7e118SToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 29610f7e118SToby Isaac if (uDof) { 29710f7e118SToby Isaac PetscInt pOff, parent, parentU; 29810f7e118SToby Isaac PetscSectionGetOffset(parentSection,uOff,&pOff);CHKERRQ(ierr); 29910f7e118SToby Isaac DMLabelGetValue(identityRef,p,&parent);CHKERRQ(ierr); 30010f7e118SToby Isaac ierr = PetscSectionGetOffset(unionSection, parent - pStart,&parentU);CHKERRQ(ierr); 30110f7e118SToby Isaac parents[pOff] = parentU; 30210f7e118SToby Isaac childIDs[pOff] = uOff; 30310f7e118SToby Isaac } 30410f7e118SToby Isaac } 30510f7e118SToby Isaac ierr = DMPlexSetTree(*ref,parentSection,parents,childIDs);CHKERRQ(ierr); 30610f7e118SToby Isaac ierr = PetscSectionDestroy(&parentSection);CHKERRQ(ierr); 30710f7e118SToby Isaac ierr = PetscFree2(parents,childIDs);CHKERRQ(ierr); 30810f7e118SToby Isaac 309da43764aSToby Isaac /* clean up */ 310da43764aSToby Isaac ierr = PetscSectionDestroy(&unionSection);CHKERRQ(ierr); 311da43764aSToby Isaac ierr = PetscSectionDestroy(&unionConeSection);CHKERRQ(ierr); 312da43764aSToby Isaac ierr = ISDestroy(&perm);CHKERRQ(ierr); 313da43764aSToby Isaac ierr = PetscFree(unionCoords);CHKERRQ(ierr); 314da43764aSToby Isaac ierr = PetscFree2(unionCones,unionOrientations);CHKERRQ(ierr); 315da43764aSToby Isaac ierr = PetscFree2(coneSizes,numDimPoints);CHKERRQ(ierr); 316da43764aSToby Isaac ierr = DMDestroy(&K);CHKERRQ(ierr); 317da43764aSToby Isaac ierr = DMDestroy(&Kref);CHKERRQ(ierr); 318da43764aSToby Isaac PetscFunctionReturn(0); 319da43764aSToby Isaac } 320da43764aSToby Isaac 321d961a43aSToby Isaac #undef __FUNCT__ 322878b19aaSToby Isaac #define __FUNCT__ "DMPlexTreeSymmetrize" 323878b19aaSToby Isaac static PetscErrorCode DMPlexTreeSymmetrize(DM dm) 324878b19aaSToby Isaac { 325878b19aaSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 326878b19aaSToby Isaac PetscSection childSec, pSec; 327878b19aaSToby Isaac PetscInt p, pSize, cSize, parMax = PETSC_MIN_INT, parMin = PETSC_MAX_INT; 328878b19aaSToby Isaac PetscInt *offsets, *children, pStart, pEnd; 329878b19aaSToby Isaac PetscErrorCode ierr; 330878b19aaSToby Isaac 331878b19aaSToby Isaac PetscFunctionBegin; 332878b19aaSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 333878b19aaSToby Isaac ierr = PetscSectionDestroy(&mesh->childSection);CHKERRQ(ierr); 334878b19aaSToby Isaac ierr = PetscFree(mesh->children);CHKERRQ(ierr); 335878b19aaSToby Isaac pSec = mesh->parentSection; 336878b19aaSToby Isaac if (!pSec) PetscFunctionReturn(0); 337878b19aaSToby Isaac ierr = PetscSectionGetStorageSize(pSec,&pSize);CHKERRQ(ierr); 338878b19aaSToby Isaac for (p = 0; p < pSize; p++) { 339878b19aaSToby Isaac PetscInt par = mesh->parents[p]; 340878b19aaSToby Isaac 341878b19aaSToby Isaac parMax = PetscMax(parMax,par+1); 342878b19aaSToby Isaac parMin = PetscMin(parMin,par); 343878b19aaSToby Isaac } 344878b19aaSToby Isaac if (parMin > parMax) { 345878b19aaSToby Isaac parMin = -1; 346878b19aaSToby Isaac parMax = -1; 347878b19aaSToby Isaac } 348878b19aaSToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)pSec),&childSec);CHKERRQ(ierr); 349878b19aaSToby Isaac ierr = PetscSectionSetChart(childSec,parMin,parMax);CHKERRQ(ierr); 350878b19aaSToby Isaac for (p = 0; p < pSize; p++) { 351878b19aaSToby Isaac PetscInt par = mesh->parents[p]; 352878b19aaSToby Isaac 353878b19aaSToby Isaac ierr = PetscSectionAddDof(childSec,par,1);CHKERRQ(ierr); 354878b19aaSToby Isaac } 355878b19aaSToby Isaac ierr = PetscSectionSetUp(childSec);CHKERRQ(ierr); 356878b19aaSToby Isaac ierr = PetscSectionGetStorageSize(childSec,&cSize);CHKERRQ(ierr); 357878b19aaSToby Isaac ierr = PetscMalloc1(cSize,&children);CHKERRQ(ierr); 358878b19aaSToby Isaac ierr = PetscCalloc1(parMax-parMin,&offsets);CHKERRQ(ierr); 359878b19aaSToby Isaac ierr = PetscSectionGetChart(pSec,&pStart,&pEnd);CHKERRQ(ierr); 360878b19aaSToby Isaac for (p = pStart; p < pEnd; p++) { 361878b19aaSToby Isaac PetscInt dof, off, i; 362878b19aaSToby Isaac 363878b19aaSToby Isaac ierr = PetscSectionGetDof(pSec,p,&dof);CHKERRQ(ierr); 364878b19aaSToby Isaac ierr = PetscSectionGetOffset(pSec,p,&off);CHKERRQ(ierr); 365878b19aaSToby Isaac for (i = 0; i < dof; i++) { 366878b19aaSToby Isaac PetscInt par = mesh->parents[off + i], cOff; 367878b19aaSToby Isaac 368878b19aaSToby Isaac ierr = PetscSectionGetOffset(childSec,par,&cOff);CHKERRQ(ierr); 369878b19aaSToby Isaac children[cOff + offsets[par-parMin]++] = p; 370878b19aaSToby Isaac } 371878b19aaSToby Isaac } 372878b19aaSToby Isaac mesh->childSection = childSec; 373878b19aaSToby Isaac mesh->children = children; 374878b19aaSToby Isaac ierr = PetscFree(offsets);CHKERRQ(ierr); 375878b19aaSToby Isaac PetscFunctionReturn(0); 376878b19aaSToby Isaac } 377878b19aaSToby Isaac 378878b19aaSToby Isaac #undef __FUNCT__ 3796dd5a8c8SToby Isaac #define __FUNCT__ "AnchorsFlatten" 3806dd5a8c8SToby Isaac static PetscErrorCode AnchorsFlatten (PetscSection section, IS is, PetscSection *sectionNew, IS *isNew) 3816dd5a8c8SToby Isaac { 3826dd5a8c8SToby Isaac PetscInt pStart, pEnd, size, sizeNew, i, p, *valsNew = NULL; 3836dd5a8c8SToby Isaac const PetscInt *vals; 3846dd5a8c8SToby Isaac PetscSection secNew; 3856dd5a8c8SToby Isaac PetscBool anyNew, globalAnyNew; 3866dd5a8c8SToby Isaac PetscBool compress; 3876dd5a8c8SToby Isaac PetscErrorCode ierr; 3886dd5a8c8SToby Isaac 3896dd5a8c8SToby Isaac PetscFunctionBegin; 3906dd5a8c8SToby Isaac ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr); 3916dd5a8c8SToby Isaac ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr); 3926dd5a8c8SToby Isaac ierr = ISGetIndices(is,&vals);CHKERRQ(ierr); 3936dd5a8c8SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secNew);CHKERRQ(ierr); 3946dd5a8c8SToby Isaac ierr = PetscSectionSetChart(secNew,pStart,pEnd);CHKERRQ(ierr); 3956dd5a8c8SToby Isaac for (i = 0; i < size; i++) { 3966dd5a8c8SToby Isaac PetscInt dof; 3976dd5a8c8SToby Isaac 3986dd5a8c8SToby Isaac p = vals[i]; 3996dd5a8c8SToby Isaac if (p < pStart || p >= pEnd) continue; 4006dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 4016dd5a8c8SToby Isaac if (dof) break; 4026dd5a8c8SToby Isaac } 4036dd5a8c8SToby Isaac if (i == size) { 4046dd5a8c8SToby Isaac ierr = PetscSectionSetUp(secNew);CHKERRQ(ierr); 4056dd5a8c8SToby Isaac anyNew = PETSC_FALSE; 4066dd5a8c8SToby Isaac compress = PETSC_FALSE; 4076dd5a8c8SToby Isaac sizeNew = 0; 4086dd5a8c8SToby Isaac } 4096dd5a8c8SToby Isaac else { 4106dd5a8c8SToby Isaac anyNew = PETSC_TRUE; 4116dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 4126dd5a8c8SToby Isaac PetscInt dof, off; 4136dd5a8c8SToby Isaac 4146dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 4156dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 4166dd5a8c8SToby Isaac for (i = 0; i < dof; i++) { 4176dd5a8c8SToby Isaac PetscInt q = vals[off + i], qDof = 0; 4186dd5a8c8SToby Isaac 4196dd5a8c8SToby Isaac if (q >= pStart && q < pEnd) { 4206dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr); 4216dd5a8c8SToby Isaac } 4226dd5a8c8SToby Isaac if (qDof) { 4236dd5a8c8SToby Isaac ierr = PetscSectionAddDof(secNew, p, qDof);CHKERRQ(ierr); 4246dd5a8c8SToby Isaac } 4256dd5a8c8SToby Isaac else { 4266dd5a8c8SToby Isaac ierr = PetscSectionAddDof(secNew, p, 1);CHKERRQ(ierr); 4276dd5a8c8SToby Isaac } 4286dd5a8c8SToby Isaac } 4296dd5a8c8SToby Isaac } 4306dd5a8c8SToby Isaac ierr = PetscSectionSetUp(secNew);CHKERRQ(ierr); 4316dd5a8c8SToby Isaac ierr = PetscSectionGetStorageSize(secNew,&sizeNew);CHKERRQ(ierr); 4326dd5a8c8SToby Isaac ierr = PetscMalloc1(sizeNew,&valsNew);CHKERRQ(ierr); 4336dd5a8c8SToby Isaac compress = PETSC_FALSE; 4346dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 4356dd5a8c8SToby Isaac PetscInt dof, off, count, offNew, dofNew; 4366dd5a8c8SToby Isaac 4376dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 4386dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 4396dd5a8c8SToby Isaac ierr = PetscSectionGetDof(secNew, p, &dofNew);CHKERRQ(ierr); 4406dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(secNew, p, &offNew);CHKERRQ(ierr); 4416dd5a8c8SToby Isaac count = 0; 4426dd5a8c8SToby Isaac for (i = 0; i < dof; i++) { 4436dd5a8c8SToby Isaac PetscInt q = vals[off + i], qDof = 0, qOff = 0, j; 4446dd5a8c8SToby Isaac 4456dd5a8c8SToby Isaac if (q >= pStart && q < pEnd) { 4466dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr); 4476dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(section, q, &qOff);CHKERRQ(ierr); 4486dd5a8c8SToby Isaac } 4496dd5a8c8SToby Isaac if (qDof) { 4506dd5a8c8SToby Isaac PetscInt oldCount = count; 4516dd5a8c8SToby Isaac 4526dd5a8c8SToby Isaac for (j = 0; j < qDof; j++) { 4536dd5a8c8SToby Isaac PetscInt k, r = vals[qOff + j]; 4546dd5a8c8SToby Isaac 4556dd5a8c8SToby Isaac for (k = 0; k < oldCount; k++) { 4566dd5a8c8SToby Isaac if (valsNew[offNew + k] == r) { 4576dd5a8c8SToby Isaac break; 4586dd5a8c8SToby Isaac } 4596dd5a8c8SToby Isaac } 4606dd5a8c8SToby Isaac if (k == oldCount) { 4616dd5a8c8SToby Isaac valsNew[offNew + count++] = r; 4626dd5a8c8SToby Isaac } 4636dd5a8c8SToby Isaac } 4646dd5a8c8SToby Isaac } 4656dd5a8c8SToby Isaac else { 4666dd5a8c8SToby Isaac PetscInt k, oldCount = count; 4676dd5a8c8SToby Isaac 4686dd5a8c8SToby Isaac for (k = 0; k < oldCount; k++) { 4696dd5a8c8SToby Isaac if (valsNew[offNew + k] == q) { 4706dd5a8c8SToby Isaac break; 4716dd5a8c8SToby Isaac } 4726dd5a8c8SToby Isaac } 4736dd5a8c8SToby Isaac if (k == oldCount) { 4746dd5a8c8SToby Isaac valsNew[offNew + count++] = q; 4756dd5a8c8SToby Isaac } 4766dd5a8c8SToby Isaac } 4776dd5a8c8SToby Isaac } 4786dd5a8c8SToby Isaac if (count < dofNew) { 4796dd5a8c8SToby Isaac ierr = PetscSectionSetDof(secNew, p, count);CHKERRQ(ierr); 4806dd5a8c8SToby Isaac compress = PETSC_TRUE; 4816dd5a8c8SToby Isaac } 4826dd5a8c8SToby Isaac } 4836dd5a8c8SToby Isaac } 4846dd5a8c8SToby Isaac ierr = ISRestoreIndices(is,&vals);CHKERRQ(ierr); 4856dd5a8c8SToby Isaac ierr = MPI_Allreduce(&anyNew,&globalAnyNew,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr); 4866dd5a8c8SToby Isaac if (!globalAnyNew) { 4876dd5a8c8SToby Isaac ierr = PetscSectionDestroy(&secNew);CHKERRQ(ierr); 4886dd5a8c8SToby Isaac *sectionNew = NULL; 4896dd5a8c8SToby Isaac *isNew = NULL; 4906dd5a8c8SToby Isaac } 4916dd5a8c8SToby Isaac else { 4926dd5a8c8SToby Isaac PetscBool globalCompress; 4936dd5a8c8SToby Isaac 4946dd5a8c8SToby Isaac ierr = MPI_Allreduce(&compress,&globalCompress,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr); 4956dd5a8c8SToby Isaac if (compress) { 4966dd5a8c8SToby Isaac PetscSection secComp; 4976dd5a8c8SToby Isaac PetscInt *valsComp = NULL; 4986dd5a8c8SToby Isaac 4996dd5a8c8SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secComp);CHKERRQ(ierr); 5006dd5a8c8SToby Isaac ierr = PetscSectionSetChart(secComp,pStart,pEnd);CHKERRQ(ierr); 5016dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 5026dd5a8c8SToby Isaac PetscInt dof; 5036dd5a8c8SToby Isaac 5046dd5a8c8SToby Isaac ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr); 5056dd5a8c8SToby Isaac ierr = PetscSectionSetDof(secComp, p, dof);CHKERRQ(ierr); 5066dd5a8c8SToby Isaac } 5076dd5a8c8SToby Isaac ierr = PetscSectionSetUp(secComp);CHKERRQ(ierr); 5086dd5a8c8SToby Isaac ierr = PetscSectionGetStorageSize(secComp,&sizeNew);CHKERRQ(ierr); 5096dd5a8c8SToby Isaac ierr = PetscMalloc1(sizeNew,&valsComp);CHKERRQ(ierr); 5106dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 5116dd5a8c8SToby Isaac PetscInt dof, off, offNew, j; 5126dd5a8c8SToby Isaac 5136dd5a8c8SToby Isaac ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr); 5146dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(secNew, p, &off);CHKERRQ(ierr); 5156dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(secComp, p, &offNew);CHKERRQ(ierr); 5166dd5a8c8SToby Isaac for (j = 0; j < dof; j++) { 5176dd5a8c8SToby Isaac valsComp[offNew + j] = valsNew[off + j]; 5186dd5a8c8SToby Isaac } 5196dd5a8c8SToby Isaac } 5206dd5a8c8SToby Isaac ierr = PetscSectionDestroy(&secNew);CHKERRQ(ierr); 5216dd5a8c8SToby Isaac secNew = secComp; 5226dd5a8c8SToby Isaac ierr = PetscFree(valsNew);CHKERRQ(ierr); 5236dd5a8c8SToby Isaac valsNew = valsComp; 5246dd5a8c8SToby Isaac } 5256dd5a8c8SToby Isaac ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),sizeNew,valsNew,PETSC_OWN_POINTER,isNew);CHKERRQ(ierr); 5266dd5a8c8SToby Isaac } 5276dd5a8c8SToby Isaac PetscFunctionReturn(0); 5286dd5a8c8SToby Isaac } 5296dd5a8c8SToby Isaac 5306dd5a8c8SToby Isaac #undef __FUNCT__ 53166af876cSToby Isaac #define __FUNCT__ "DMPlexComputeConstraints_Tree" 53266af876cSToby Isaac static PetscErrorCode DMPlexComputeConstraints_Tree(DM dm) 53366af876cSToby Isaac { 53466af876cSToby Isaac PetscInt p, pStart, pEnd, *anchors, size; 53566af876cSToby Isaac PetscInt aMin = PETSC_MAX_INT, aMax = PETSC_MIN_INT; 53666af876cSToby Isaac PetscSection aSec; 53766af876cSToby Isaac IS aIS; 53866af876cSToby Isaac PetscErrorCode ierr; 53966af876cSToby Isaac 54066af876cSToby Isaac PetscFunctionBegin; 54166af876cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 54266af876cSToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 54366af876cSToby Isaac for (p = pStart; p < pEnd; p++) { 54466af876cSToby Isaac PetscInt parent; 54566af876cSToby Isaac 54666af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 54766af876cSToby Isaac if (parent != p) { 54866af876cSToby Isaac aMin = PetscMin(aMin,p); 54966af876cSToby Isaac aMax = PetscMax(aMax,p+1); 55066af876cSToby Isaac } 55166af876cSToby Isaac } 55266af876cSToby Isaac if (aMin > aMax) { 55366af876cSToby Isaac aMin = -1; 55466af876cSToby Isaac aMax = -1; 55566af876cSToby Isaac } 55666af876cSToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&aSec);CHKERRQ(ierr); 55766af876cSToby Isaac ierr = PetscSectionSetChart(aSec,aMin,aMax);CHKERRQ(ierr); 55866af876cSToby Isaac for (p = aMin; p < aMax; p++) { 55966af876cSToby Isaac PetscInt parent, ancestor = p; 56066af876cSToby Isaac 56166af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 56266af876cSToby Isaac while (parent != ancestor) { 56366af876cSToby Isaac ancestor = parent; 56466af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr); 56566af876cSToby Isaac } 56666af876cSToby Isaac if (ancestor != p) { 56766af876cSToby Isaac PetscInt closureSize, *closure = NULL; 56866af876cSToby Isaac 56966af876cSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 57066af876cSToby Isaac ierr = PetscSectionSetDof(aSec,p,closureSize);CHKERRQ(ierr); 57166af876cSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 57266af876cSToby Isaac } 57366af876cSToby Isaac } 57466af876cSToby Isaac ierr = PetscSectionSetUp(aSec);CHKERRQ(ierr); 57566af876cSToby Isaac ierr = PetscSectionGetStorageSize(aSec,&size);CHKERRQ(ierr); 57666af876cSToby Isaac ierr = PetscMalloc1(size,&anchors);CHKERRQ(ierr); 57766af876cSToby Isaac for (p = aMin; p < aMax; p++) { 57866af876cSToby Isaac PetscInt parent, ancestor = p; 57966af876cSToby Isaac 58066af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 58166af876cSToby Isaac while (parent != ancestor) { 58266af876cSToby Isaac ancestor = parent; 58366af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr); 58466af876cSToby Isaac } 58566af876cSToby Isaac if (ancestor != p) { 58666af876cSToby Isaac PetscInt j, closureSize, *closure = NULL, aOff; 58766af876cSToby Isaac 58866af876cSToby Isaac ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr); 58966af876cSToby Isaac 59066af876cSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 59166af876cSToby Isaac for (j = 0; j < closureSize; j++) { 59266af876cSToby Isaac anchors[aOff + j] = closure[2*j]; 59366af876cSToby Isaac } 59466af876cSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 59566af876cSToby Isaac } 59666af876cSToby Isaac } 59766af876cSToby Isaac ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm),size,anchors,PETSC_OWN_POINTER,&aIS);CHKERRQ(ierr); 5986dd5a8c8SToby Isaac { 5996dd5a8c8SToby Isaac PetscSection aSecNew = aSec; 6006dd5a8c8SToby Isaac IS aISNew = aIS; 6016dd5a8c8SToby Isaac 6026dd5a8c8SToby Isaac ierr = PetscObjectReference((PetscObject)aSec);CHKERRQ(ierr); 6036dd5a8c8SToby Isaac ierr = PetscObjectReference((PetscObject)aIS);CHKERRQ(ierr); 6046dd5a8c8SToby Isaac while (aSecNew) { 6056dd5a8c8SToby Isaac ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr); 6066dd5a8c8SToby Isaac ierr = ISDestroy(&aIS);CHKERRQ(ierr); 6076dd5a8c8SToby Isaac aSec = aSecNew; 6086dd5a8c8SToby Isaac aIS = aISNew; 6096dd5a8c8SToby Isaac aSecNew = NULL; 6106dd5a8c8SToby Isaac aISNew = NULL; 6116dd5a8c8SToby Isaac ierr = AnchorsFlatten(aSec,aIS,&aSecNew,&aISNew);CHKERRQ(ierr); 6126dd5a8c8SToby Isaac } 6136dd5a8c8SToby Isaac } 61466af876cSToby Isaac ierr = DMPlexSetConstraints(dm,aSec,aIS);CHKERRQ(ierr); 61566af876cSToby Isaac ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr); 61666af876cSToby Isaac ierr = ISDestroy(&aIS);CHKERRQ(ierr); 61766af876cSToby Isaac PetscFunctionReturn(0); 61866af876cSToby Isaac } 61966af876cSToby Isaac 62066af876cSToby Isaac #undef __FUNCT__ 6210b7167a0SToby Isaac #define __FUNCT__ "DMPlexSetTree" 6220b7167a0SToby Isaac /*@ 6230b7167a0SToby Isaac DMPlexSetTree - set the tree that describes the hierarchy of non-conforming mesh points. This routine also creates 6240b7167a0SToby Isaac the point-to-point constraints determined by the tree: a point is constained to the points in the closure of its 6250b7167a0SToby Isaac tree root. 6260b7167a0SToby Isaac 6270b7167a0SToby Isaac Collective on dm 6280b7167a0SToby Isaac 6290b7167a0SToby Isaac Input Parameters: 6300b7167a0SToby Isaac + dm - the DMPlex object 6310b7167a0SToby Isaac . parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section 6320b7167a0SToby Isaac offset indexes the parent and childID list; the reference count of parentSection is incremented 6330b7167a0SToby Isaac . parents - a list of the point parents; copied, can be destroyed 6340b7167a0SToby Isaac - childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then 6350b7167a0SToby Isaac the child corresponds to the point in the reference tree with index childIDs; copied, can be destroyed 6360b7167a0SToby Isaac 6370b7167a0SToby Isaac Level: intermediate 6380b7167a0SToby Isaac 639b2f41788SToby Isaac .seealso: DMPlexGetTree(), DMPlexSetReferenceTree(), DMPlexSetConstraints(), DMPlexGetTreeParent(), DMPlexGetTreeChildren() 6400b7167a0SToby Isaac @*/ 641b2f41788SToby Isaac PetscErrorCode DMPlexSetTree(DM dm, PetscSection parentSection, PetscInt parents[], PetscInt childIDs[]) 6420b7167a0SToby Isaac { 6430b7167a0SToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 6440b7167a0SToby Isaac PetscInt size; 6450b7167a0SToby Isaac PetscErrorCode ierr; 6460b7167a0SToby Isaac 6470b7167a0SToby Isaac PetscFunctionBegin; 6480b7167a0SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6490b7167a0SToby Isaac PetscValidHeaderSpecific(parentSection, PETSC_SECTION_CLASSID, 2); 6500b7167a0SToby Isaac ierr = PetscObjectReference((PetscObject)parentSection);CHKERRQ(ierr); 6510b7167a0SToby Isaac ierr = PetscSectionDestroy(&mesh->parentSection);CHKERRQ(ierr); 6520b7167a0SToby Isaac mesh->parentSection = parentSection; 6530b7167a0SToby Isaac ierr = PetscSectionGetStorageSize(parentSection,&size);CHKERRQ(ierr); 6540b7167a0SToby Isaac if (parents != mesh->parents) { 6550b7167a0SToby Isaac ierr = PetscFree(mesh->parents);CHKERRQ(ierr); 6560b7167a0SToby Isaac ierr = PetscMalloc1(size,&mesh->parents);CHKERRQ(ierr); 6570b7167a0SToby Isaac ierr = PetscMemcpy(mesh->parents, parents, size * sizeof(*parents));CHKERRQ(ierr); 6580b7167a0SToby Isaac } 6590b7167a0SToby Isaac if (childIDs != mesh->childIDs) { 6600b7167a0SToby Isaac ierr = PetscFree(mesh->childIDs);CHKERRQ(ierr); 6610b7167a0SToby Isaac ierr = PetscMalloc1(size,&mesh->childIDs);CHKERRQ(ierr); 6620b7167a0SToby Isaac ierr = PetscMemcpy(mesh->childIDs, childIDs, size * sizeof(*childIDs));CHKERRQ(ierr); 6630b7167a0SToby Isaac } 664878b19aaSToby Isaac ierr = DMPlexTreeSymmetrize(dm);CHKERRQ(ierr); 66566af876cSToby Isaac ierr = DMPlexComputeConstraints_Tree(dm);CHKERRQ(ierr); 6660b7167a0SToby Isaac PetscFunctionReturn(0); 6670b7167a0SToby Isaac } 6680b7167a0SToby Isaac 6690b7167a0SToby Isaac #undef __FUNCT__ 670b2f41788SToby Isaac #define __FUNCT__ "DMPlexGetTree" 671b2f41788SToby Isaac /*@ 672b2f41788SToby Isaac DMPlexGetTree - get the tree that describes the hierarchy of non-conforming mesh points. 673b2f41788SToby Isaac Collective on dm 674b2f41788SToby Isaac 675b2f41788SToby Isaac Input Parameters: 676b2f41788SToby Isaac . dm - the DMPlex object 677b2f41788SToby Isaac 678b2f41788SToby Isaac Output Parameters: 679b2f41788SToby Isaac + parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section 680b2f41788SToby Isaac offset indexes the parent and childID list 681b2f41788SToby Isaac . parents - a list of the point parents 682b2f41788SToby Isaac . childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then 683b2f41788SToby Isaac the child corresponds to the point in the reference tree with index childID 684b2f41788SToby Isaac . childSection - the inverse of the parent section 685b2f41788SToby Isaac - children - a list of the point children 686b2f41788SToby Isaac 687b2f41788SToby Isaac Level: intermediate 688b2f41788SToby Isaac 689b2f41788SToby Isaac .seealso: DMPlexSetTree(), DMPlexSetReferenceTree(), DMPlexSetConstraints(), DMPlexGetTreeParent(), DMPlexGetTreeChildren() 690b2f41788SToby Isaac @*/ 691b2f41788SToby Isaac PetscErrorCode DMPlexGetTree(DM dm, PetscSection *parentSection, PetscInt *parents[], PetscInt *childIDs[], PetscSection *childSection, PetscInt *children[]) 692b2f41788SToby Isaac { 693b2f41788SToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 694b2f41788SToby Isaac 695b2f41788SToby Isaac PetscFunctionBegin; 696b2f41788SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 697b2f41788SToby Isaac if (parentSection) *parentSection = mesh->parentSection; 698b2f41788SToby Isaac if (parents) *parents = mesh->parents; 699b2f41788SToby Isaac if (childIDs) *childIDs = mesh->childIDs; 700b2f41788SToby Isaac if (childSection) *childSection = mesh->childSection; 701b2f41788SToby Isaac if (children) *children = mesh->children; 702b2f41788SToby Isaac PetscFunctionReturn(0); 703b2f41788SToby Isaac } 704b2f41788SToby Isaac 705b2f41788SToby Isaac #undef __FUNCT__ 706d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeParent" 707d961a43aSToby Isaac /*@ 708d961a43aSToby Isaac DMPlexGetTreeParent - get the parent of a point in the tree describing the point hierarchy (not the Sieve DAG) 709d961a43aSToby Isaac 710d961a43aSToby Isaac Input Parameters: 711d961a43aSToby Isaac + dm - the DMPlex object 712d961a43aSToby Isaac - point - the query point 713d961a43aSToby Isaac 714d961a43aSToby Isaac Output Parameters: 715d961a43aSToby Isaac + parent - if not NULL, set to the parent of the point, or the point itself if the point does not have a parent 716d961a43aSToby Isaac - childID - if not NULL, set to the child ID of the point with respect to its parent, or 0 if the point 717d961a43aSToby Isaac does not have a parent 718d961a43aSToby Isaac 719d961a43aSToby Isaac Level: intermediate 720d961a43aSToby Isaac 721d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeChildren() 722d961a43aSToby Isaac @*/ 723d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeParent(DM dm, PetscInt point, PetscInt *parent, PetscInt *childID) 724d961a43aSToby Isaac { 725d961a43aSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 726d961a43aSToby Isaac PetscSection pSec; 727d961a43aSToby Isaac PetscErrorCode ierr; 728d961a43aSToby Isaac 729d961a43aSToby Isaac PetscFunctionBegin; 730d961a43aSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 731d961a43aSToby Isaac pSec = mesh->parentSection; 732d961a43aSToby Isaac if (pSec && point >= pSec->pStart && point < pSec->pEnd) { 733d961a43aSToby Isaac PetscInt dof; 734d961a43aSToby Isaac 735d961a43aSToby Isaac ierr = PetscSectionGetDof (pSec, point, &dof);CHKERRQ(ierr); 736d961a43aSToby Isaac if (dof) { 737d961a43aSToby Isaac PetscInt off; 738d961a43aSToby Isaac 739d961a43aSToby Isaac ierr = PetscSectionGetOffset (pSec, point, &off);CHKERRQ(ierr); 740d961a43aSToby Isaac if (parent) *parent = mesh->parents[off]; 741d961a43aSToby Isaac if (childID) *childID = mesh->childIDs[off]; 742d961a43aSToby Isaac PetscFunctionReturn(0); 743d961a43aSToby Isaac } 744d961a43aSToby Isaac } 745d961a43aSToby Isaac if (parent) { 746d961a43aSToby Isaac *parent = point; 747d961a43aSToby Isaac } 748d961a43aSToby Isaac if (childID) { 749d961a43aSToby Isaac *childID = 0; 750d961a43aSToby Isaac } 751d961a43aSToby Isaac PetscFunctionReturn(0); 752d961a43aSToby Isaac } 753d961a43aSToby Isaac 754d961a43aSToby Isaac #undef __FUNCT__ 755d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeChildren" 756d961a43aSToby Isaac /*@C 757d961a43aSToby Isaac DMPlexGetTreeChildren - get the children of a point in the tree describing the point hierarchy (not the Sieve DAG) 758d961a43aSToby Isaac 759d961a43aSToby Isaac Input Parameters: 760d961a43aSToby Isaac + dm - the DMPlex object 761d961a43aSToby Isaac - point - the query point 762d961a43aSToby Isaac 763d961a43aSToby Isaac Output Parameters: 764d961a43aSToby Isaac + numChildren - if not NULL, set to the number of children 765d961a43aSToby Isaac - children - if not NULL, set to a list children, or set to NULL if the point has no children 766d961a43aSToby Isaac 767d961a43aSToby Isaac Level: intermediate 768d961a43aSToby Isaac 769d961a43aSToby Isaac Fortran Notes: 770d961a43aSToby Isaac Since it returns an array, this routine is only available in Fortran 90, and you must 771d961a43aSToby Isaac include petsc.h90 in your code. 772d961a43aSToby Isaac 773d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeParent() 774d961a43aSToby Isaac @*/ 775d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeChildren(DM dm, PetscInt point, PetscInt *numChildren, const PetscInt *children[]) 776d961a43aSToby Isaac { 777d961a43aSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 778d961a43aSToby Isaac PetscSection childSec; 779d961a43aSToby Isaac PetscInt dof = 0; 780d961a43aSToby Isaac PetscErrorCode ierr; 781d961a43aSToby Isaac 782d961a43aSToby Isaac PetscFunctionBegin; 783d961a43aSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 784d961a43aSToby Isaac childSec = mesh->childSection; 785d961a43aSToby Isaac if (childSec && point >= childSec->pStart && point < childSec->pEnd) { 786d961a43aSToby Isaac ierr = PetscSectionGetDof (childSec, point, &dof);CHKERRQ(ierr); 787d961a43aSToby Isaac } 788d961a43aSToby Isaac if (numChildren) *numChildren = dof; 789d961a43aSToby Isaac if (children) { 790d961a43aSToby Isaac if (dof) { 791d961a43aSToby Isaac PetscInt off; 792d961a43aSToby Isaac 793d961a43aSToby Isaac ierr = PetscSectionGetOffset (childSec, point, &off);CHKERRQ(ierr); 794d961a43aSToby Isaac *children = &mesh->children[off]; 795d961a43aSToby Isaac } 796d961a43aSToby Isaac else { 797d961a43aSToby Isaac *children = NULL; 798d961a43aSToby Isaac } 799d961a43aSToby Isaac } 800d961a43aSToby Isaac PetscFunctionReturn(0); 801d961a43aSToby Isaac } 8020c37af3bSToby Isaac 8030c37af3bSToby Isaac #undef __FUNCT__ 8040c37af3bSToby Isaac #define __FUNCT__ "DMPlexComputeConstraintMatrix_ReferenceTree" 8050c37af3bSToby Isaac PetscErrorCode DMPlexComputeConstraintMatrix_ReferenceTree(DM dm) 8060c37af3bSToby Isaac { 8070c37af3bSToby Isaac PetscDS ds; 8080c37af3bSToby Isaac PetscInt spdim; 8090c37af3bSToby Isaac PetscInt numFields, f, c, cStart, cEnd, pStart, pEnd, conStart, conEnd; 8100c37af3bSToby Isaac const PetscInt *anchors; 8110c37af3bSToby Isaac PetscSection section, cSec, aSec; 8120c37af3bSToby Isaac Mat cMat; 8130c37af3bSToby Isaac PetscReal *v0, *v0parent, *vtmp, *J, *Jparent, *invJparent, detJ, detJparent; 8140c37af3bSToby Isaac IS aIS; 8150c37af3bSToby Isaac PetscErrorCode ierr; 8160c37af3bSToby Isaac 8170c37af3bSToby Isaac PetscFunctionBegin; 8180c37af3bSToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 8190c37af3bSToby Isaac ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 8200c37af3bSToby Isaac ierr = PetscDSGetNumFields(ds,&numFields);CHKERRQ(ierr); 8210c37af3bSToby Isaac ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 8220c37af3bSToby Isaac ierr = DMGetDefaultSection(dm,§ion);CHKERRQ(ierr); 8230c37af3bSToby Isaac ierr = DMPlexGetConstraintSection(dm,&cSec);CHKERRQ(ierr); 8240c37af3bSToby Isaac ierr = DMPlexGetConstraints(dm,&aSec,&aIS);CHKERRQ(ierr); 8250c37af3bSToby Isaac ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr); 8260c37af3bSToby Isaac ierr = PetscSectionGetChart(cSec,&conStart,&conEnd);CHKERRQ(ierr); 8270c37af3bSToby Isaac ierr = DMPlexGetConstraintMatrix(dm,&cMat);CHKERRQ(ierr); 8280c37af3bSToby Isaac ierr = DMPlexGetDimension(dm,&spdim);CHKERRQ(ierr); 8290c37af3bSToby Isaac ierr = PetscMalloc6(spdim,&v0,spdim,&v0parent,spdim,&vtmp,spdim*spdim,&J,spdim*spdim,&Jparent,spdim*spdim,&invJparent);CHKERRQ(ierr); 8300c37af3bSToby Isaac 8310c37af3bSToby Isaac for (f = 0; f < numFields; f++) { 8320c37af3bSToby Isaac PetscFE fe; 8330c37af3bSToby Isaac PetscDualSpace space; 8340c37af3bSToby Isaac PetscInt i, j, k, nPoints, offset; 8350c37af3bSToby Isaac PetscInt fSize, fComp; 8360c37af3bSToby Isaac PetscScalar *B = NULL; 8370c37af3bSToby Isaac PetscReal *weights, *pointsRef, *pointsReal; 8380c37af3bSToby Isaac Mat Amat, Bmat, Xmat; 8390c37af3bSToby Isaac 8400c37af3bSToby Isaac ierr = PetscDSGetDiscretization(ds,f,(PetscObject *)(&fe));CHKERRQ(ierr); 8410c37af3bSToby Isaac ierr = PetscFEGetDualSpace(fe,&space);CHKERRQ(ierr); 8420c37af3bSToby Isaac ierr = PetscDualSpaceGetDimension(space,&fSize);CHKERRQ(ierr); 8430c37af3bSToby Isaac ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr); 8440c37af3bSToby Isaac ierr = MatCreate(PETSC_COMM_SELF,&Amat);CHKERRQ(ierr); 8450c37af3bSToby Isaac ierr = MatSetSizes(Amat,fSize,fSize,fSize,fSize);CHKERRQ(ierr); 8460c37af3bSToby Isaac ierr = MatSetType(Amat,MATSEQDENSE);CHKERRQ(ierr); 8470c37af3bSToby Isaac ierr = MatSetUp(Amat);CHKERRQ(ierr); 8480c37af3bSToby Isaac ierr = MatDuplicate(Amat,MAT_DO_NOT_COPY_VALUES,&Bmat);CHKERRQ(ierr); 8490c37af3bSToby Isaac ierr = MatDuplicate(Amat,MAT_DO_NOT_COPY_VALUES,&Xmat);CHKERRQ(ierr); 8500c37af3bSToby Isaac nPoints = 0; 8510c37af3bSToby Isaac for (i = 0; i < fSize; i++) { 8520c37af3bSToby Isaac PetscInt qPoints; 8530c37af3bSToby Isaac PetscQuadrature quad; 8540c37af3bSToby Isaac 8550c37af3bSToby Isaac ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 8560c37af3bSToby Isaac ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr); 8570c37af3bSToby Isaac nPoints += qPoints; 8580c37af3bSToby Isaac } 8590c37af3bSToby Isaac ierr = PetscMalloc3(nPoints,&weights,spdim*nPoints,&pointsRef,spdim*nPoints,&pointsReal);CHKERRQ(ierr); 8600c37af3bSToby Isaac offset = 0; 8610c37af3bSToby Isaac for (i = 0; i < fSize; i++) { 8620c37af3bSToby Isaac PetscInt qPoints; 8630c37af3bSToby Isaac const PetscReal *p, *w; 8640c37af3bSToby Isaac PetscQuadrature quad; 8650c37af3bSToby Isaac 8660c37af3bSToby Isaac ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 8670c37af3bSToby Isaac ierr = PetscQuadratureGetData(quad,NULL,&qPoints,&p,&w);CHKERRQ(ierr); 8680c37af3bSToby Isaac ierr = PetscMemcpy(weights+offset,w,qPoints*sizeof(*w));CHKERRQ(ierr); 8690c37af3bSToby Isaac ierr = PetscMemcpy(pointsRef+spdim*offset,p,spdim*qPoints*sizeof(*p));CHKERRQ(ierr); 8700c37af3bSToby Isaac offset += qPoints; 8710c37af3bSToby Isaac } 8720c37af3bSToby Isaac ierr = PetscFEGetTabulation(fe,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr); 8730c37af3bSToby Isaac offset = 0; 8740c37af3bSToby Isaac for (i = 0; i < fSize; i++) { 8750c37af3bSToby Isaac PetscInt qPoints; 8760c37af3bSToby Isaac PetscQuadrature quad; 8770c37af3bSToby Isaac 8780c37af3bSToby Isaac ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 8790c37af3bSToby Isaac ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr); 8800c37af3bSToby Isaac for (j = 0; j < fSize; j++) { 8810c37af3bSToby Isaac PetscScalar val = 0.; 8820c37af3bSToby Isaac 8830c37af3bSToby Isaac for (k = 0; k < qPoints; k++) { 8840c37af3bSToby Isaac val += B[((offset + k) * fSize + j) * fComp] * weights[k]; 8850c37af3bSToby Isaac } 8860c37af3bSToby Isaac ierr = MatSetValue(Amat,i,j,val,INSERT_VALUES);CHKERRQ(ierr); 8870c37af3bSToby Isaac } 8880c37af3bSToby Isaac offset += qPoints; 8890c37af3bSToby Isaac } 8900c37af3bSToby Isaac ierr = PetscFERestoreTabulation(fe,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr); 8910c37af3bSToby Isaac ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8920c37af3bSToby Isaac ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8930c37af3bSToby Isaac ierr = MatLUFactor(Amat,NULL,NULL,NULL);CHKERRQ(ierr); 8940c37af3bSToby Isaac for (c = cStart; c < cEnd; c++) { 8950c37af3bSToby Isaac PetscInt parent; 8960c37af3bSToby Isaac PetscInt closureSize, closureSizeP, *closure = NULL, *closureP = NULL; 8970c37af3bSToby Isaac PetscInt *childOffsets, *parentOffsets; 8980c37af3bSToby Isaac 8990c37af3bSToby Isaac ierr = DMPlexGetTreeParent(dm,c,&parent,NULL);CHKERRQ(ierr); 9000c37af3bSToby Isaac if (parent == c) continue; 9010c37af3bSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 9020c37af3bSToby Isaac for (i = 0; i < closureSize; i++) { 9030c37af3bSToby Isaac PetscInt p = closure[2*i]; 9040c37af3bSToby Isaac PetscInt conDof; 9050c37af3bSToby Isaac 9060c37af3bSToby Isaac if (p < conStart || p >= conEnd) continue; 9070c37af3bSToby Isaac if (numFields > 1) { 9080c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(cSec,p,f,&conDof);CHKERRQ(ierr); 9090c37af3bSToby Isaac } 9100c37af3bSToby Isaac else { 9110c37af3bSToby Isaac ierr = PetscSectionGetDof(cSec,p,&conDof);CHKERRQ(ierr); 9120c37af3bSToby Isaac } 9130c37af3bSToby Isaac if (conDof) break; 9140c37af3bSToby Isaac } 9150c37af3bSToby Isaac if (i == closureSize) { 9160c37af3bSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 9170c37af3bSToby Isaac continue; 9180c37af3bSToby Isaac } 9190c37af3bSToby Isaac 9200c37af3bSToby Isaac ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr); 9210c37af3bSToby Isaac ierr = DMPlexComputeCellGeometry(dm, parent, v0parent, Jparent, invJparent, &detJparent);CHKERRQ(ierr); 9220c37af3bSToby Isaac for (i = 0; i < nPoints; i++) { 9230c37af3bSToby Isaac CoordinatesRefToReal(spdim, spdim, v0, J, &pointsRef[i*spdim],vtmp);CHKERRQ(ierr); 9240c37af3bSToby Isaac CoordinatesRealToRef(spdim, spdim, v0parent, invJparent, vtmp, &pointsReal[i*spdim]);CHKERRQ(ierr); 9250c37af3bSToby Isaac } 9260c37af3bSToby Isaac ierr = PetscFEGetTabulation(fe,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr); 9270c37af3bSToby Isaac offset = 0; 9280c37af3bSToby Isaac for (i = 0; i < fSize; i++) { 9290c37af3bSToby Isaac PetscInt qPoints; 9300c37af3bSToby Isaac PetscQuadrature quad; 9310c37af3bSToby Isaac 9320c37af3bSToby Isaac ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 9330c37af3bSToby Isaac ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr); 9340c37af3bSToby Isaac for (j = 0; j < fSize; j++) { 9350c37af3bSToby Isaac PetscScalar val = 0.; 9360c37af3bSToby Isaac 9370c37af3bSToby Isaac for (k = 0; k < qPoints; k++) { 9380c37af3bSToby Isaac val += B[((offset + k) * fSize + j) * fComp] * weights[k]; 9390c37af3bSToby Isaac } 9400c37af3bSToby Isaac MatSetValue(Bmat,i,j,val,INSERT_VALUES);CHKERRQ(ierr); 9410c37af3bSToby Isaac } 9420c37af3bSToby Isaac offset += qPoints; 9430c37af3bSToby Isaac } 9440c37af3bSToby Isaac ierr = PetscFERestoreTabulation(fe,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr); 9450c37af3bSToby Isaac ierr = MatAssemblyBegin(Bmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9460c37af3bSToby Isaac ierr = MatAssemblyEnd(Bmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9470c37af3bSToby Isaac ierr = MatMatSolve(Amat,Bmat,Xmat);CHKERRQ(ierr); 9480c37af3bSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSizeP,&closureP);CHKERRQ(ierr); 9490c37af3bSToby Isaac ierr = PetscMalloc2(closureSize+1,&childOffsets,closureSizeP+1,&parentOffsets);CHKERRQ(ierr); 9500c37af3bSToby Isaac childOffsets[0] = 0; 9510c37af3bSToby Isaac for (i = 0; i < closureSize; i++) { 9520c37af3bSToby Isaac PetscInt p = closure[2*i]; 9530c37af3bSToby Isaac PetscInt dof; 9540c37af3bSToby Isaac 9550c37af3bSToby Isaac if (numFields > 1) { 9560c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(section,p,f,&dof);CHKERRQ(ierr); 9570c37af3bSToby Isaac } 9580c37af3bSToby Isaac else { 9590c37af3bSToby Isaac ierr = PetscSectionGetDof(section,p,&dof);CHKERRQ(ierr); 9600c37af3bSToby Isaac } 9610c37af3bSToby Isaac childOffsets[i+1]=childOffsets[i]+dof / fComp; 9620c37af3bSToby Isaac } 9630c37af3bSToby Isaac parentOffsets[0] = 0; 9640c37af3bSToby Isaac for (i = 0; i < closureSizeP; i++) { 9650c37af3bSToby Isaac PetscInt p = closureP[2*i]; 9660c37af3bSToby Isaac PetscInt dof; 9670c37af3bSToby Isaac 9680c37af3bSToby Isaac if (numFields > 1) { 9690c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(section,p,f,&dof);CHKERRQ(ierr); 9700c37af3bSToby Isaac } 9710c37af3bSToby Isaac else { 9720c37af3bSToby Isaac ierr = PetscSectionGetDof(section,p,&dof);CHKERRQ(ierr); 9730c37af3bSToby Isaac } 9740c37af3bSToby Isaac parentOffsets[i+1]=parentOffsets[i]+dof / fComp; 9750c37af3bSToby Isaac } 9760c37af3bSToby Isaac for (i = 0; i < closureSize; i++) { 9770c37af3bSToby Isaac PetscInt conDof, conOff, aDof, aOff; 9780c37af3bSToby Isaac PetscInt p = closure[2*i]; 9790c37af3bSToby Isaac PetscInt o = closure[2*i+1]; 9800c37af3bSToby Isaac 9810c37af3bSToby Isaac if (p < conStart || p >= conEnd) continue; 9820c37af3bSToby Isaac if (numFields > 1) { 9830c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(cSec,p,f,&conDof);CHKERRQ(ierr); 9840c37af3bSToby Isaac ierr = PetscSectionGetFieldOffset(cSec,p,f,&conOff);CHKERRQ(ierr); 9850c37af3bSToby Isaac } 9860c37af3bSToby Isaac else { 9870c37af3bSToby Isaac ierr = PetscSectionGetDof(cSec,p,&conDof);CHKERRQ(ierr); 9880c37af3bSToby Isaac ierr = PetscSectionGetOffset(cSec,p,&conOff);CHKERRQ(ierr); 9890c37af3bSToby Isaac } 9900c37af3bSToby Isaac if (!conDof) continue; 9910c37af3bSToby Isaac ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr); 9920c37af3bSToby Isaac ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr); 9930c37af3bSToby Isaac for (k = 0; k < aDof; k++) { 9940c37af3bSToby Isaac PetscInt a = anchors[aOff + k]; 9950c37af3bSToby Isaac PetscInt aSecDof, aSecOff; 9960c37af3bSToby Isaac 9970c37af3bSToby Isaac if (numFields > 1) { 9980c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(section,a,f,&aSecDof);CHKERRQ(ierr); 9990c37af3bSToby Isaac ierr = PetscSectionGetFieldOffset(section,a,f,&aSecOff);CHKERRQ(ierr); 10000c37af3bSToby Isaac } 10010c37af3bSToby Isaac else { 10020c37af3bSToby Isaac ierr = PetscSectionGetDof(section,a,&aSecDof);CHKERRQ(ierr); 10030c37af3bSToby Isaac ierr = PetscSectionGetOffset(section,a,&aSecOff);CHKERRQ(ierr); 10040c37af3bSToby Isaac } 10050c37af3bSToby Isaac if (!aSecDof) continue; 10060c37af3bSToby Isaac 10070c37af3bSToby Isaac for (j = 0; j < closureSizeP; j++) { 10080c37af3bSToby Isaac PetscInt q = closureP[2*j]; 10090c37af3bSToby Isaac PetscInt oq = closureP[2*j+1]; 10100c37af3bSToby Isaac 10110c37af3bSToby Isaac if (q == a) { 10120c37af3bSToby Isaac PetscInt r, s, t; 10130c37af3bSToby Isaac 10140c37af3bSToby Isaac for (r = childOffsets[i]; r < childOffsets[i+1]; r++) { 10150c37af3bSToby Isaac for (s = parentOffsets[j]; s < parentOffsets[j+1]; s++) { 10160c37af3bSToby Isaac PetscScalar val; 10170c37af3bSToby Isaac PetscInt insertCol, insertRow; 10180c37af3bSToby Isaac 10190c37af3bSToby Isaac ierr = MatGetValue(Xmat,r,s,&val);CHKERRQ(ierr); 10200c37af3bSToby Isaac if (o >= 0) { 10210c37af3bSToby Isaac insertRow = conOff + fComp * (r - childOffsets[i]); 10220c37af3bSToby Isaac } 10230c37af3bSToby Isaac else { 10240c37af3bSToby Isaac insertRow = conOff + fComp * (childOffsets[i + 1] - 1 - r); 10250c37af3bSToby Isaac } 10260c37af3bSToby Isaac if (oq >= 0) { 10270c37af3bSToby Isaac insertCol = aSecOff + fComp * (s - parentOffsets[j]); 10280c37af3bSToby Isaac } 10290c37af3bSToby Isaac else { 10300c37af3bSToby Isaac insertCol = aSecOff + fComp * (parentOffsets[j + 1] - 1 - s); 10310c37af3bSToby Isaac } 10320c37af3bSToby Isaac for (t = 0; t < fComp; t++) { 10330c37af3bSToby Isaac ierr = MatSetValue(cMat,insertRow + t,insertCol + t,val,INSERT_VALUES);CHKERRQ(ierr); 10340c37af3bSToby Isaac } 10350c37af3bSToby Isaac } 10360c37af3bSToby Isaac } 10370c37af3bSToby Isaac } 10380c37af3bSToby Isaac } 10390c37af3bSToby Isaac } 10400c37af3bSToby Isaac } 10410c37af3bSToby Isaac ierr = PetscFree2(childOffsets,parentOffsets);CHKERRQ(ierr); 10420c37af3bSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 10430c37af3bSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSizeP,&closureP);CHKERRQ(ierr); 10440c37af3bSToby Isaac } 10450c37af3bSToby Isaac ierr = MatDestroy(&Amat);CHKERRQ(ierr); 10460c37af3bSToby Isaac ierr = MatDestroy(&Bmat);CHKERRQ(ierr); 10470c37af3bSToby Isaac ierr = MatDestroy(&Xmat);CHKERRQ(ierr); 10480c37af3bSToby Isaac ierr = PetscFree3(weights,pointsRef,pointsReal);CHKERRQ(ierr); 10490c37af3bSToby Isaac } 10500c37af3bSToby Isaac ierr = MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10510c37af3bSToby Isaac ierr = MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10520c37af3bSToby Isaac ierr = PetscFree6(v0,v0parent,vtmp,J,Jparent,invJparent);CHKERRQ(ierr); 10530c37af3bSToby Isaac ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr); 10540c37af3bSToby Isaac 10550c37af3bSToby Isaac PetscFunctionReturn(0); 10560c37af3bSToby Isaac } 1057*95a0b26dSToby Isaac 1058*95a0b26dSToby Isaac #undef __FUNCT__ 1059*95a0b26dSToby Isaac #define __FUNCT__ "DMPlexComputeConstraintMatrix_Tree" 1060*95a0b26dSToby Isaac PetscErrorCode DMPlexComputeConstraintMatrix_Tree(DM dm) 1061*95a0b26dSToby Isaac { 1062*95a0b26dSToby Isaac DM refTree; 1063*95a0b26dSToby Isaac PetscDS ds; 1064*95a0b26dSToby Isaac Mat refCmat, cMat; 1065*95a0b26dSToby Isaac PetscInt numFields, f, pRefStart, pRefEnd, p, *rows, *cols, maxDof, maxAnDof, *perm, *iperm, pStart, pEnd, conStart, conEnd, **refPointFieldN; 1066*95a0b26dSToby Isaac PetscScalar ***refPointFieldMats, *pointWork; 1067*95a0b26dSToby Isaac PetscSection refConSec, conSec, refAnSec, anSec, section, refSection; 1068*95a0b26dSToby Isaac IS refAnIS, anIS; 1069*95a0b26dSToby Isaac const PetscInt *refAnchors, *anchors; 1070*95a0b26dSToby Isaac PetscErrorCode ierr; 1071*95a0b26dSToby Isaac 1072*95a0b26dSToby Isaac PetscFunctionBegin; 1073*95a0b26dSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1074*95a0b26dSToby Isaac ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 1075*95a0b26dSToby Isaac ierr = PetscDSGetNumFields(ds,&numFields);CHKERRQ(ierr); 1076*95a0b26dSToby Isaac ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr); 1077*95a0b26dSToby Isaac ierr = DMSetDS(refTree,ds);CHKERRQ(ierr); 1078*95a0b26dSToby Isaac ierr = DMPlexComputeConstraintMatrix_ReferenceTree(refTree);CHKERRQ(ierr); 1079*95a0b26dSToby Isaac ierr = DMPlexGetConstraintMatrix(refTree,&refCmat);CHKERRQ(ierr); 1080*95a0b26dSToby Isaac ierr = DMPlexGetChart(refTree,&pRefStart,&pRefEnd);CHKERRQ(ierr); 1081*95a0b26dSToby Isaac ierr = DMPlexGetConstraints(refTree,&refAnSec,&refAnIS);CHKERRQ(ierr); 1082*95a0b26dSToby Isaac ierr = DMPlexGetConstraints(dm,&anSec,&anIS);CHKERRQ(ierr); 1083*95a0b26dSToby Isaac ierr = ISGetIndices(refAnIS,&refAnchors);CHKERRQ(ierr); 1084*95a0b26dSToby Isaac ierr = ISGetIndices(anIS,&anchors);CHKERRQ(ierr); 1085*95a0b26dSToby Isaac ierr = DMPlexGetConstraintSection(refTree,&refConSec);CHKERRQ(ierr); 1086*95a0b26dSToby Isaac ierr = DMPlexGetConstraintSection(dm,&conSec);CHKERRQ(ierr); 1087*95a0b26dSToby Isaac ierr = DMGetDefaultSection(refTree,&refSection);CHKERRQ(ierr); 1088*95a0b26dSToby Isaac ierr = DMGetDefaultSection(dm,§ion);CHKERRQ(ierr); 1089*95a0b26dSToby Isaac ierr = PetscSectionGetChart(refConSec,&pRefStart,&pRefEnd);CHKERRQ(ierr); 1090*95a0b26dSToby Isaac ierr = PetscSectionGetChart(conSec,&conStart,&conEnd);CHKERRQ(ierr); 1091*95a0b26dSToby Isaac ierr = PetscMalloc1(pRefEnd-pRefStart,&refPointFieldMats);CHKERRQ(ierr); 1092*95a0b26dSToby Isaac ierr = PetscMalloc1(pRefEnd-pRefStart,&refPointFieldN);CHKERRQ(ierr); 1093*95a0b26dSToby Isaac ierr = PetscSectionGetMaxDof(refConSec,&maxDof);CHKERRQ(ierr); 1094*95a0b26dSToby Isaac ierr = PetscSectionGetMaxDof(refAnSec,&maxAnDof);CHKERRQ(ierr); 1095*95a0b26dSToby Isaac ierr = PetscMalloc1(maxDof,&rows);CHKERRQ(ierr); 1096*95a0b26dSToby Isaac ierr = PetscMalloc1(maxDof*maxAnDof,&cols);CHKERRQ(ierr); 1097*95a0b26dSToby Isaac ierr = PetscMalloc1(maxDof*maxDof*maxAnDof,&pointWork);CHKERRQ(ierr); 1098*95a0b26dSToby Isaac 1099*95a0b26dSToby Isaac /* step 1: get submats for every constrained point in the reference tree */ 1100*95a0b26dSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 1101*95a0b26dSToby Isaac PetscInt parent, closureSize, *closure = NULL, pDof; 1102*95a0b26dSToby Isaac 1103*95a0b26dSToby Isaac ierr = DMPlexGetTreeParent(refTree,p,&parent,NULL);CHKERRQ(ierr); 1104*95a0b26dSToby Isaac ierr = PetscSectionGetDof(refConSec,p,&pDof);CHKERRQ(ierr); 1105*95a0b26dSToby Isaac if (!pDof || parent == p) continue; 1106*95a0b26dSToby Isaac 1107*95a0b26dSToby Isaac ierr = PetscMalloc1(numFields,&refPointFieldMats[p-pRefStart]);CHKERRQ(ierr); 1108*95a0b26dSToby Isaac ierr = PetscCalloc1(numFields,&refPointFieldN[p-pRefStart]);CHKERRQ(ierr); 1109*95a0b26dSToby Isaac ierr = DMPlexGetTransitiveClosure(refTree,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1110*95a0b26dSToby Isaac for (f = 0; f < numFields; f++) { 1111*95a0b26dSToby Isaac PetscInt cDof, cOff, numCols, r, i, fComp; 1112*95a0b26dSToby Isaac PetscFE fe; 1113*95a0b26dSToby Isaac 1114*95a0b26dSToby Isaac ierr = PetscDSGetDiscretization(ds,f,(PetscObject *) &fe);CHKERRQ(ierr); 1115*95a0b26dSToby Isaac ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr); 1116*95a0b26dSToby Isaac 1117*95a0b26dSToby Isaac if (numFields > 1) { 1118*95a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(refConSec,p,f,&cDof);CHKERRQ(ierr); 1119*95a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(refConSec,p,f,&cOff);CHKERRQ(ierr); 1120*95a0b26dSToby Isaac } 1121*95a0b26dSToby Isaac else { 1122*95a0b26dSToby Isaac ierr = PetscSectionGetDof(refConSec,p,&cDof);CHKERRQ(ierr); 1123*95a0b26dSToby Isaac ierr = PetscSectionGetOffset(refConSec,p,&cOff);CHKERRQ(ierr); 1124*95a0b26dSToby Isaac } 1125*95a0b26dSToby Isaac 1126*95a0b26dSToby Isaac if (!cDof) continue; 1127*95a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 1128*95a0b26dSToby Isaac rows[r] = cOff + r; 1129*95a0b26dSToby Isaac } 1130*95a0b26dSToby Isaac numCols = 0; 1131*95a0b26dSToby Isaac for (i = 0; i < closureSize; i++) { 1132*95a0b26dSToby Isaac PetscInt q = closure[2*i]; 1133*95a0b26dSToby Isaac PetscInt o = closure[2*i+1]; 1134*95a0b26dSToby Isaac PetscInt aDof, aOff, j; 1135*95a0b26dSToby Isaac 1136*95a0b26dSToby Isaac if (numFields > 1) { 1137*95a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(refSection,q,f,&aDof);CHKERRQ(ierr); 1138*95a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(refSection,q,f,&aOff);CHKERRQ(ierr); 1139*95a0b26dSToby Isaac } 1140*95a0b26dSToby Isaac else { 1141*95a0b26dSToby Isaac ierr = PetscSectionGetDof(refSection,q,&aDof);CHKERRQ(ierr); 1142*95a0b26dSToby Isaac ierr = PetscSectionGetOffset(refSection,q,&aOff);CHKERRQ(ierr); 1143*95a0b26dSToby Isaac } 1144*95a0b26dSToby Isaac 1145*95a0b26dSToby Isaac for (j = 0; j < aDof; j++) { 1146*95a0b26dSToby Isaac PetscInt node = (o >= 0) ? (j / fComp) : ((aDof - 1 - j) / fComp); 1147*95a0b26dSToby Isaac PetscInt comp = (j % fComp); 1148*95a0b26dSToby Isaac 1149*95a0b26dSToby Isaac cols[numCols++] = aOff + node * fComp + comp; 1150*95a0b26dSToby Isaac } 1151*95a0b26dSToby Isaac } 1152*95a0b26dSToby Isaac refPointFieldN[p-pRefStart][f] = numCols; 1153*95a0b26dSToby Isaac ierr = PetscMalloc1(cDof*numCols,&refPointFieldMats[p-pRefStart][f]);CHKERRQ(ierr); 1154*95a0b26dSToby Isaac ierr = MatGetValues(refCmat,cDof,rows,numCols,cols,refPointFieldMats[p-pRefStart][f]);CHKERRQ(ierr); 1155*95a0b26dSToby Isaac } 1156*95a0b26dSToby Isaac ierr = DMPlexRestoreTransitiveClosure(refTree,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1157*95a0b26dSToby Isaac } 1158*95a0b26dSToby Isaac 1159*95a0b26dSToby Isaac /* step 2: compute the preorder */ 1160*95a0b26dSToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 1161*95a0b26dSToby Isaac ierr = PetscMalloc2(pEnd-pStart,&perm,pEnd-pStart,&iperm);CHKERRQ(ierr); 1162*95a0b26dSToby Isaac for (p = pStart; p < pEnd; p++) { 1163*95a0b26dSToby Isaac perm[p - pStart] = p; 1164*95a0b26dSToby Isaac iperm[p - pStart] = p-pStart; 1165*95a0b26dSToby Isaac } 1166*95a0b26dSToby Isaac for (p = 0; p < pEnd - pStart;) { 1167*95a0b26dSToby Isaac PetscInt point = perm[p]; 1168*95a0b26dSToby Isaac PetscInt parent; 1169*95a0b26dSToby Isaac 1170*95a0b26dSToby Isaac ierr = DMPlexGetTreeParent(dm,point,&parent,NULL);CHKERRQ(ierr); 1171*95a0b26dSToby Isaac if (parent == point) { 1172*95a0b26dSToby Isaac p++; 1173*95a0b26dSToby Isaac } 1174*95a0b26dSToby Isaac else { 1175*95a0b26dSToby Isaac PetscInt size, closureSize, *closure = NULL, i; 1176*95a0b26dSToby Isaac 1177*95a0b26dSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1178*95a0b26dSToby Isaac for (i = 0; i < closureSize; i++) { 1179*95a0b26dSToby Isaac PetscInt q = closure[2*i]; 1180*95a0b26dSToby Isaac if (iperm[q-pStart] > iperm[point-pStart]) { 1181*95a0b26dSToby Isaac /* swap */ 1182*95a0b26dSToby Isaac perm[p] = q; 1183*95a0b26dSToby Isaac perm[iperm[q-pStart]] = point; 1184*95a0b26dSToby Isaac iperm[point-pStart] = iperm[q-pStart]; 1185*95a0b26dSToby Isaac iperm[q-pStart] = p; 1186*95a0b26dSToby Isaac break; 1187*95a0b26dSToby Isaac } 1188*95a0b26dSToby Isaac } 1189*95a0b26dSToby Isaac size = closureSize; 1190*95a0b26dSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1191*95a0b26dSToby Isaac if (i == size) { 1192*95a0b26dSToby Isaac p++; 1193*95a0b26dSToby Isaac } 1194*95a0b26dSToby Isaac } 1195*95a0b26dSToby Isaac } 1196*95a0b26dSToby Isaac 1197*95a0b26dSToby Isaac /* step 3: fill the constraint matrix */ 1198*95a0b26dSToby Isaac /* we are going to use a preorder progressive fill strategy. Mat doesn't 1199*95a0b26dSToby Isaac * allow progressive fill without assembly, so we are going to set up the 1200*95a0b26dSToby Isaac * values outside of the Mat first. 1201*95a0b26dSToby Isaac */ 1202*95a0b26dSToby Isaac { 1203*95a0b26dSToby Isaac PetscInt nRows, row, nnz; 1204*95a0b26dSToby Isaac PetscBool done; 1205*95a0b26dSToby Isaac const PetscInt *ia, *ja; 1206*95a0b26dSToby Isaac PetscScalar *vals; 1207*95a0b26dSToby Isaac 1208*95a0b26dSToby Isaac ierr = DMPlexGetConstraintMatrix(dm,&cMat);CHKERRQ(ierr); 1209*95a0b26dSToby Isaac ierr = MatGetRowIJ(cMat,0,PETSC_FALSE,PETSC_FALSE,&nRows,&ia,&ja,&done);CHKERRQ(ierr); 1210*95a0b26dSToby Isaac if (!done) SETERRQ(PetscObjectComm((PetscObject)cMat),PETSC_ERR_PLIB,"Could not get RowIJ of constraint matrix"); 1211*95a0b26dSToby Isaac nnz = ia[nRows]; 1212*95a0b26dSToby Isaac /* malloc and then zero rows right before we fill them: this way valgrind 1213*95a0b26dSToby Isaac * can tell if we are doing progressive fill in the wrong order */ 1214*95a0b26dSToby Isaac ierr = PetscMalloc1(nnz,&vals);CHKERRQ(ierr); 1215*95a0b26dSToby Isaac for (p = 0; p < pEnd - pStart; p++) { 1216*95a0b26dSToby Isaac PetscInt parent, childid, closureSize, *closure = NULL; 1217*95a0b26dSToby Isaac PetscInt point = perm[p], pointDof; 1218*95a0b26dSToby Isaac 1219*95a0b26dSToby Isaac ierr = DMPlexGetTreeParent(dm,point,&parent,&childid);CHKERRQ(ierr); 1220*95a0b26dSToby Isaac if ((point < conStart) || (point >= conEnd) || (parent == point)) continue; 1221*95a0b26dSToby Isaac ierr = PetscSectionGetDof(conSec,point,&pointDof);CHKERRQ(ierr); 1222*95a0b26dSToby Isaac if (!pointDof) continue; 1223*95a0b26dSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1224*95a0b26dSToby Isaac for (f = 0; f < numFields; f++) { 1225*95a0b26dSToby Isaac PetscInt cDof, cOff, numCols, numFillCols, i, r, fComp, matOffset, offset; 1226*95a0b26dSToby Isaac PetscScalar *pointMat; 1227*95a0b26dSToby Isaac PetscFE fe; 1228*95a0b26dSToby Isaac 1229*95a0b26dSToby Isaac ierr = PetscDSGetDiscretization(ds,f,(PetscObject *) &fe);CHKERRQ(ierr); 1230*95a0b26dSToby Isaac ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr); 1231*95a0b26dSToby Isaac 1232*95a0b26dSToby Isaac if (numFields > 1) { 1233*95a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(conSec,point,f,&cDof);CHKERRQ(ierr); 1234*95a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(conSec,point,f,&cOff);CHKERRQ(ierr); 1235*95a0b26dSToby Isaac } 1236*95a0b26dSToby Isaac else { 1237*95a0b26dSToby Isaac ierr = PetscSectionGetDof(conSec,point,&cDof);CHKERRQ(ierr); 1238*95a0b26dSToby Isaac ierr = PetscSectionGetOffset(conSec,point,&cOff);CHKERRQ(ierr); 1239*95a0b26dSToby Isaac } 1240*95a0b26dSToby Isaac if (!cDof) continue; 1241*95a0b26dSToby Isaac 1242*95a0b26dSToby Isaac /* make sure that every row for this point is the same size */ 1243*95a0b26dSToby Isaac #if defined(PETSC_USE_DEBUG) 1244*95a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 1245*95a0b26dSToby Isaac if (cDof > 1 && r) { 1246*95a0b26dSToby Isaac if ((ia[cOff+r+1]-ia[cOff+r]) != (ia[cOff+r]-ia[cOff+r-1])) { 1247*95a0b26dSToby 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])); 1248*95a0b26dSToby Isaac } 1249*95a0b26dSToby Isaac } 1250*95a0b26dSToby Isaac } 1251*95a0b26dSToby Isaac #endif 1252*95a0b26dSToby Isaac /* zero rows */ 1253*95a0b26dSToby Isaac for (i = ia[cOff] ; i< ia[cOff+cDof];i++) { 1254*95a0b26dSToby Isaac vals[i] = 0.; 1255*95a0b26dSToby Isaac } 1256*95a0b26dSToby Isaac matOffset = ia[cOff]; 1257*95a0b26dSToby Isaac numFillCols = ia[cOff+1] - matOffset; 1258*95a0b26dSToby Isaac pointMat = refPointFieldMats[childid-pRefStart][f]; 1259*95a0b26dSToby Isaac numCols = refPointFieldN[childid-pRefStart][f]; 1260*95a0b26dSToby Isaac offset = 0; 1261*95a0b26dSToby Isaac for (i = 0; i < closureSize; i++) { 1262*95a0b26dSToby Isaac PetscInt q = closure[2*i]; 1263*95a0b26dSToby Isaac PetscInt o = closure[2*i+1]; 1264*95a0b26dSToby Isaac PetscInt aDof, aOff, j, k, qConDof, qConOff; 1265*95a0b26dSToby Isaac 1266*95a0b26dSToby Isaac qConDof = qConOff = 0; 1267*95a0b26dSToby Isaac if (numFields > 1) { 1268*95a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(section,q,f,&aDof);CHKERRQ(ierr); 1269*95a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(section,q,f,&aOff);CHKERRQ(ierr); 1270*95a0b26dSToby Isaac if (q >= conStart && q < conEnd) { 1271*95a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(conSec,q,f,&qConDof);CHKERRQ(ierr); 1272*95a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(conSec,q,f,&qConOff);CHKERRQ(ierr); 1273*95a0b26dSToby Isaac } 1274*95a0b26dSToby Isaac } 1275*95a0b26dSToby Isaac else { 1276*95a0b26dSToby Isaac ierr = PetscSectionGetDof(section,q,&aDof);CHKERRQ(ierr); 1277*95a0b26dSToby Isaac ierr = PetscSectionGetOffset(section,q,&aOff);CHKERRQ(ierr); 1278*95a0b26dSToby Isaac if (q >= conStart && q < conEnd) { 1279*95a0b26dSToby Isaac ierr = PetscSectionGetDof(conSec,q,&qConDof);CHKERRQ(ierr); 1280*95a0b26dSToby Isaac ierr = PetscSectionGetOffset(conSec,q,&qConOff);CHKERRQ(ierr); 1281*95a0b26dSToby Isaac } 1282*95a0b26dSToby Isaac } 1283*95a0b26dSToby Isaac if (!aDof) continue; 1284*95a0b26dSToby Isaac if (qConDof) { 1285*95a0b26dSToby Isaac /* this point has anchors: its rows of the matrix should already 1286*95a0b26dSToby Isaac * be filled, thanks to preordering */ 1287*95a0b26dSToby Isaac /* first multiply into pointWork, then set in matrix */ 1288*95a0b26dSToby Isaac PetscInt aMatOffset = ia[qConOff]; 1289*95a0b26dSToby Isaac PetscInt aNumFillCols = ia[qConOff + 1] - aMatOffset; 1290*95a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 1291*95a0b26dSToby Isaac for (j = 0; j < aNumFillCols; j++) { 1292*95a0b26dSToby Isaac PetscScalar inVal = 0; 1293*95a0b26dSToby Isaac for (k = 0; k < aDof; k++) { 1294*95a0b26dSToby Isaac PetscInt node = (o >= 0) ? (k / fComp) : ((aDof - 1 - k) / fComp); 1295*95a0b26dSToby Isaac PetscInt comp = (k % fComp); 1296*95a0b26dSToby Isaac PetscInt col = node * fComp + comp; 1297*95a0b26dSToby Isaac 1298*95a0b26dSToby Isaac inVal += pointMat[r * numCols + offset + col] * vals[aMatOffset + aNumFillCols * k + j]; 1299*95a0b26dSToby Isaac } 1300*95a0b26dSToby Isaac pointWork[r * aNumFillCols + j] = inVal; 1301*95a0b26dSToby Isaac } 1302*95a0b26dSToby Isaac } 1303*95a0b26dSToby Isaac /* assume that the columns are sorted, spend less time searching */ 1304*95a0b26dSToby Isaac for (j = 0, k = 0; j < aNumFillCols; j++) { 1305*95a0b26dSToby Isaac PetscInt col = ja[aMatOffset + j]; 1306*95a0b26dSToby Isaac for (;k < numFillCols; k++) { 1307*95a0b26dSToby Isaac if (ja[matOffset + k] == col) { 1308*95a0b26dSToby Isaac break; 1309*95a0b26dSToby Isaac } 1310*95a0b26dSToby Isaac } 1311*95a0b26dSToby Isaac if (k == numFillCols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"No nonzero space for (%d, %d)", cOff, col); 1312*95a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 1313*95a0b26dSToby Isaac vals[matOffset + numFillCols * r + k] = pointWork[r * aNumFillCols + j]; 1314*95a0b26dSToby Isaac } 1315*95a0b26dSToby Isaac } 1316*95a0b26dSToby Isaac } 1317*95a0b26dSToby Isaac else { 1318*95a0b26dSToby Isaac /* find where to put this portion of pointMat into the matrix */ 1319*95a0b26dSToby Isaac for (k = 0; k < numFillCols; k++) { 1320*95a0b26dSToby Isaac if (ja[matOffset + k] == aOff) { 1321*95a0b26dSToby Isaac break; 1322*95a0b26dSToby Isaac } 1323*95a0b26dSToby Isaac } 1324*95a0b26dSToby Isaac if (k == numFillCols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"No nonzero space for (%d, %d)", cOff, aOff); 1325*95a0b26dSToby Isaac for (j = 0; j < aDof; j++) { 1326*95a0b26dSToby Isaac PetscInt node = (o >= 0) ? (j / fComp) : ((aDof - 1 - j) / fComp); 1327*95a0b26dSToby Isaac PetscInt comp = (j % fComp); 1328*95a0b26dSToby Isaac PetscInt col = node * fComp + comp; 1329*95a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 1330*95a0b26dSToby Isaac vals[matOffset + numFillCols * r + k + col] += pointMat[r * numCols + offset + col]; 1331*95a0b26dSToby Isaac } 1332*95a0b26dSToby Isaac } 1333*95a0b26dSToby Isaac } 1334*95a0b26dSToby Isaac offset += aDof; 1335*95a0b26dSToby Isaac } 1336*95a0b26dSToby Isaac } 1337*95a0b26dSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1338*95a0b26dSToby Isaac } 1339*95a0b26dSToby Isaac for (row = 0; row < nRows; row++) { 1340*95a0b26dSToby Isaac ierr = MatSetValues(cMat,1,&row,ia[row+1]-ia[row],&ja[ia[row]],&vals[ia[row]],INSERT_VALUES);CHKERRQ(ierr); 1341*95a0b26dSToby Isaac } 1342*95a0b26dSToby Isaac ierr = MatRestoreRowIJ(cMat,0,PETSC_FALSE,PETSC_FALSE,&nRows,&ia,&ja,&done);CHKERRQ(ierr); 1343*95a0b26dSToby Isaac if (!done) SETERRQ(PetscObjectComm((PetscObject)cMat),PETSC_ERR_PLIB,"Could not restore RowIJ of constraint matrix"); 1344*95a0b26dSToby Isaac ierr = MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1345*95a0b26dSToby Isaac ierr = MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1346*95a0b26dSToby Isaac ierr = PetscFree(vals);CHKERRQ(ierr); 1347*95a0b26dSToby Isaac } 1348*95a0b26dSToby Isaac 1349*95a0b26dSToby Isaac /* clean up */ 1350*95a0b26dSToby Isaac ierr = ISRestoreIndices(refAnIS,&refAnchors);CHKERRQ(ierr); 1351*95a0b26dSToby Isaac ierr = ISRestoreIndices(anIS,&anchors);CHKERRQ(ierr); 1352*95a0b26dSToby Isaac ierr = PetscFree2(perm,iperm);CHKERRQ(ierr); 1353*95a0b26dSToby Isaac ierr = PetscFree(rows);CHKERRQ(ierr); 1354*95a0b26dSToby Isaac ierr = PetscFree(cols);CHKERRQ(ierr); 1355*95a0b26dSToby Isaac ierr = PetscFree(pointWork);CHKERRQ(ierr); 1356*95a0b26dSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 1357*95a0b26dSToby Isaac PetscInt parent, pDof; 1358*95a0b26dSToby Isaac 1359*95a0b26dSToby Isaac ierr = DMPlexGetTreeParent(refTree,p,&parent,NULL);CHKERRQ(ierr); 1360*95a0b26dSToby Isaac ierr = PetscSectionGetDof(refConSec,p,&pDof);CHKERRQ(ierr); 1361*95a0b26dSToby Isaac if (!pDof || parent == p) continue; 1362*95a0b26dSToby Isaac 1363*95a0b26dSToby Isaac for (f = 0; f < numFields; f++) { 1364*95a0b26dSToby Isaac PetscInt cDof; 1365*95a0b26dSToby Isaac 1366*95a0b26dSToby Isaac if (numFields > 1) { 1367*95a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(refConSec,p,f,&cDof);CHKERRQ(ierr); 1368*95a0b26dSToby Isaac } 1369*95a0b26dSToby Isaac else { 1370*95a0b26dSToby Isaac ierr = PetscSectionGetDof(refConSec,p,&cDof);CHKERRQ(ierr); 1371*95a0b26dSToby Isaac } 1372*95a0b26dSToby Isaac 1373*95a0b26dSToby Isaac if (!cDof) continue; 1374*95a0b26dSToby Isaac ierr = PetscFree(refPointFieldMats[p - pRefStart][f]);CHKERRQ(ierr); 1375*95a0b26dSToby Isaac } 1376*95a0b26dSToby Isaac ierr = PetscFree(refPointFieldMats[p - pRefStart]);CHKERRQ(ierr); 1377*95a0b26dSToby Isaac ierr = PetscFree(refPointFieldN[p - pRefStart]);CHKERRQ(ierr); 1378*95a0b26dSToby Isaac } 1379*95a0b26dSToby Isaac ierr = PetscFree(refPointFieldMats);CHKERRQ(ierr); 1380*95a0b26dSToby Isaac ierr = PetscFree(refPointFieldN);CHKERRQ(ierr); 1381*95a0b26dSToby Isaac PetscFunctionReturn(0); 1382*95a0b26dSToby Isaac } 1383*95a0b26dSToby Isaac 1384