PetscErrorCode PETSCDM_DLLEXPORT MeshRefineSingularity_Fichera(ALE::Obj<ALE::Mesh> mesh, MPI_Comm comm, double * singularity, double factor, ALE::Obj<ALE::Mesh> * refinedMesh, PetscTruth interpolate = PETSC_FALSE) { ALE::Obj<ALE::Mesh> oldMesh = mesh; double oldLimit; PetscErrorCode ierr; PetscFunctionBegin; //ierr = MeshGetMesh(mesh, oldMesh);CHKERRQ(ierr); //ierr = MeshCreate(comm, refinedMesh);CHKERRQ(ierr); int dim = oldMesh->getDimension(); oldLimit = oldMesh->getMaxVolume(); //double oldLimInv = 1./oldLimit; double curLimit, tmpLimit; double minLimit = oldLimit/16384.; //arbitrary; const ALE::Obj<ALE::Mesh::real_section_type>& coordinates = oldMesh->getRealSection("coordinates"); const ALE::Obj<ALE::Mesh::real_section_type>& volume_limits = oldMesh->getRealSection("volume_limits"); volume_limits->setFiberDimension(oldMesh->heightStratum(0), 1); oldMesh->allocate(volume_limits); const ALE::Obj<ALE::Mesh::label_sequence>& cells = oldMesh->heightStratum(0); ALE::Mesh::label_sequence::iterator c_iter = cells->begin(); ALE::Mesh::label_sequence::iterator c_iter_end = cells->end(); double centerCoords[dim]; while (c_iter != c_iter_end) { const double * coords = oldMesh->restrictClosure(coordinates, *c_iter); for (int i = 0; i < dim; i++) { centerCoords[i] = 0; for (int j = 0; j < dim+1; j++) { centerCoords[i] += coords[j*dim+i]; } centerCoords[i] = centerCoords[i]/(dim+1); //PetscPrintf(oldMesh->comm(), "%f, ", centerCoords[i]); } //PetscPrintf(oldMesh->comm(), "\n"); double dist = 0.; double cornerdist = 0.; //HERE'S THE DIFFERENCE: if centercoords is less than the singularity coordinate for each direction, include that direction in the distance /* for (int i = 0; i < dim; i++) { if (centerCoords[i] <= singularity[i]) { dist += (centerCoords[i] - singularity[i])*(centerCoords[i] - singularity[i]); } } */ //determine: the per-dimension distance: cases if (dim > 2) { for (int i = 0; i < dim; i++) { cornerdist = 0.; if (centerCoords[i] > singularity[i]) { for (int j = 0; j < dim; j++) { if (j != i) cornerdist += (centerCoords[j] - singularity[j])*(centerCoords[j] - singularity[j]); } if (cornerdist < dist || dist == 0.) dist = cornerdist; } } } //patch up AROUND the corner by minimizing between the distance from the relevant axis and the singular vertex double singdist = 0.; for (int i = 0; i < dim; i++) { singdist += (centerCoords[i] - singularity[i])*(centerCoords[i] - singularity[i]); } if (singdist < dist || dist == 0.) dist = singdist; if (dist > 0.) { dist = sqrt(dist); double mu = pow(dist, factor); //PetscPrintf(oldMesh->comm(), "%f, %f\n", mu, dist); tmpLimit = oldLimit*pow(mu, dim); if (tmpLimit > minLimit) { curLimit = tmpLimit; } else curLimit = minLimit; } else curLimit = minLimit; //PetscPrintf(oldMesh->comm(), "%f, %f\n", dist, tmpLimit); volume_limits->updatePoint(*c_iter, &curLimit); c_iter++; } ALE::Obj<ALE::Mesh> newMesh = ALE::Generator::refineMesh(oldMesh, volume_limits, interpolate); //ierr = MeshSetMesh(*refinedMesh, newMesh);CHKERRQ(ierr); *refinedMesh = newMesh; const ALE::Obj<ALE::Mesh::real_section_type>& s = newMesh->getRealSection("default"); const Obj<std::set<std::string> >& discs = oldMesh->getDiscretizations(); for(std::set<std::string>::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter) { newMesh->setDiscretization(*f_iter, oldMesh->getDiscretization(*f_iter)); } newMesh->setupField(s); // PetscPrintf(newMesh->comm(), "refined\n"); PetscFunctionReturn(0); }
EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "DMConvert_DA_Mesh" PetscErrorCode DMConvert_DA_Mesh(DM dm, const DMType newtype, DM *dmNew) { PetscSection section; DM cda; DMDALocalInfo info; Vec coordinates; PetscInt *cone, *coneO; PetscInt dim, M, N, P, numCells, numGlobalCells, numCorners, numVertices, c = 0, v = 0; PetscInt ye, ze; PetscInt debug = 0; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscOptionsGetInt(PETSC_NULL, "-dm_mesh_debug", &debug, PETSC_NULL);CHKERRQ(ierr); ierr = DMDAGetInfo(dm, &dim, &M, &N, &P, 0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr); if (info.sw > 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Currently, only DMDAs with unti stencil width can be converted to DMMeshes."); /* In order to get a partition of cells, rather than vertices, we give each process the cells between vertices it owns and also higher numbered ghost vertices (vertices to the right and up) */ numCorners = 1 << dim; numCells = ((info.gxm+info.gxs - info.xs) - 1); if (dim > 1) {numCells *= ((info.gym+info.gys - info.ys) - 1);} if (dim > 2) {numCells *= ((info.gzm+info.gzs - info.zs) - 1);} numVertices = (info.gxm+info.gxs - info.xs); if (dim > 1) {numVertices *= (info.gym+info.gys - info.ys);} if (dim > 2) {numVertices *= (info.gzm+info.gzs - info.zs);} numGlobalCells = M-1; if (dim > 1) {numGlobalCells *= N-1;} if (dim > 2) {numGlobalCells *= P-1;} ALE::Obj<PETSC_MESH_TYPE> mesh = new PETSC_MESH_TYPE(((PetscObject) dm)->comm, info.dim, debug); ALE::Obj<PETSC_MESH_TYPE::sieve_type> sieve = new PETSC_MESH_TYPE::sieve_type(((PetscObject) dm)->comm, 0, numCells+numVertices, debug); PETSC_MESH_TYPE::renumbering_type renumbering; mesh->setSieve(sieve); /* Number each cell for the vertex in the lower left corner */ if (dim < 3) {ze = 1; P = 1;} else {ze = info.gzs+info.gzm-1;} if (dim < 2) {ye = 1; N = 1;} else {ye = info.gys+info.gym-1;} for(PetscInt k = info.zs; k < ze; ++k) { for(PetscInt j = info.ys; j < ye; ++j) { for(PetscInt i = info.xs; i < info.gxs+info.gxm-1; ++i, ++c) { PetscInt globalC = (k*(N-1) + j)*(M-1) + i; renumbering[globalC] = c; sieve->setConeSize(c, numCorners); } } } if (c != numCells) {SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "Error in generated cell numbering, %d should be %d", c, numCells);} /* Get vertex renumbering */ for(PetscInt k = info.zs; k < info.gzs+info.gzm; ++k) { for(PetscInt j = info.ys; j < info.gys+info.gym; ++j) { for(PetscInt i = info.xs; i < info.gxs+info.gxm; ++i, ++v) { PetscInt globalV = (k*N + j)*M + i + numGlobalCells; renumbering[globalV] = v+numCells; } } } if (v != numVertices) {SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "Error in generated vertex numbering, %d should be %d", v, numVertices);} /* Calculate support sizes */ for(PetscInt k = info.zs; k < ze; ++k, ++c) { for(PetscInt j = info.ys; j < ye; ++j) { for(PetscInt i = info.xs; i < info.gxs+info.gxm-1; ++i) { for(PetscInt kp = k; kp <= k+(dim>2); ++kp) { for(PetscInt jp = j; jp <= j+(dim>1); ++jp) { for(PetscInt ip = i; ip <= i+1; ++ip) { PetscInt globalV = (kp*N + jp)*M + ip + numGlobalCells; sieve->addSupportSize(renumbering[globalV], 1); } } } } } } sieve->allocate(); ierr = PetscMalloc2(numCorners,PetscInt,&cone,numCorners,PetscInt,&coneO);CHKERRQ(ierr); for(PetscInt v = 0; v < numCorners; ++v) { coneO[v] = 1; } for(PetscInt k = info.zs; k < ze; ++k) { for(PetscInt j = info.ys; j < ye; ++j) { for(PetscInt i = info.xs; i < info.gxs+info.gxm-1; ++i) { PetscInt globalC = (k*(N-1) + j)*(M-1) + i; PetscInt v = 0; cone[v++] = renumbering[(k*N + j)*M + i+0 + numGlobalCells]; cone[v++] = renumbering[(k*N + j)*M + i+1 + numGlobalCells]; if (dim > 1) { cone[v++] = renumbering[(k*N + j+1)*M + i+0 + numGlobalCells]; cone[v++] = renumbering[(k*N + j+1)*M + i+1 + numGlobalCells]; } if (dim > 2) { cone[v++] = renumbering[((k+1)*N + j+0)*M + i+0 + numGlobalCells]; cone[v++] = renumbering[((k+1)*N + j+0)*M + i+1 + numGlobalCells]; cone[v++] = renumbering[((k+1)*N + j+1)*M + i+0 + numGlobalCells]; cone[v++] = renumbering[((k+1)*N + j+1)*M + i+1 + numGlobalCells]; } sieve->setCone(cone, renumbering[globalC]); sieve->setConeOrientation(coneO, renumbering[globalC]); } } } ierr = PetscFree2(cone,coneO);CHKERRQ(ierr); sieve->symmetrize(); mesh->stratify(); /* Create boundary marker */ { const Obj<PETSC_MESH_TYPE::label_type>& boundary = mesh->createLabel("marker"); for(PetscInt k = info.zs; k < info.gzs+info.gzm; ++k) { for(PetscInt j = info.ys; j < info.gys+info.gym; ++j) { if (info.xs == 0) { PetscInt globalV = (k*N + j)*M + info.xs + numGlobalCells; mesh->setValue(boundary, renumbering[globalV], 1); } if (info.gxs+info.gxm-1 == M-1) { PetscInt globalV = (k*N + j)*M + info.gxs+info.gxm-1 + numGlobalCells; mesh->setValue(boundary, renumbering[globalV], 1); } } } if (dim > 1) { for(PetscInt k = info.zs; k < info.gzs+info.gzm; ++k) { for(PetscInt i = info.xs; i < info.gxs+info.gxm; ++i) { if (info.ys == 0) { PetscInt globalV = (k*N + info.ys)*M + i + numGlobalCells; mesh->setValue(boundary, renumbering[globalV], 1); } if (info.gys+info.gym-1 == N-1) { PetscInt globalV = (k*N + info.gys+info.gym-1)*M + i + numGlobalCells; mesh->setValue(boundary, renumbering[globalV], 1); } } } } if (dim > 2) { for(PetscInt j = info.ys; j < info.gys+info.gym; ++j) { for(PetscInt i = info.xs; i < info.gxs+info.gxm; ++i) { if (info.zs == 0) { PetscInt globalV = (info.zs*N + j)*M + i + numGlobalCells; mesh->setValue(boundary, renumbering[globalV], 1); } if (info.gzs+info.gzm-1 == P-1) { PetscInt globalV = ((info.gzs+info.gzm-1)*N + j)*M + i + numGlobalCells; mesh->setValue(boundary, renumbering[globalV], 1); } } } } } /* Create new DM */ ierr = DMMeshCreate(((PetscObject) dm)->comm, dmNew);CHKERRQ(ierr); ierr = DMMeshSetMesh(*dmNew, mesh);CHKERRQ(ierr); /* Set coordinates */ ierr = PetscSectionCreate(((PetscObject) dm)->comm, §ion);CHKERRQ(ierr); ierr = PetscSectionSetChart(section, numCells, numCells+numVertices);CHKERRQ(ierr); for(PetscInt v = numCells; v < numCells+numVertices; ++v) { ierr = PetscSectionSetDof(section, v, dim);CHKERRQ(ierr); } ierr = PetscSectionSetUp(section);CHKERRQ(ierr); ierr = DMMeshSetCoordinateSection(*dmNew, section);CHKERRQ(ierr); ierr = DMDAGetCoordinateDA(dm, &cda);CHKERRQ(ierr); ierr = DMDAGetGhostedCoordinates(dm, &coordinates);CHKERRQ(ierr); { Obj<PETSC_MESH_TYPE::real_section_type> coordSection = mesh->getRealSection("coordinates"); switch(dim) { case 1: { PetscScalar **coords; ierr = DMDAVecGetArrayDOF(cda, coordinates, &coords);CHKERRQ(ierr); for(PetscInt i = info.xs; i < info.gxs+info.gxm; ++i) { PetscInt globalV = i + numGlobalCells; coordSection->updatePoint(renumbering[globalV], coords[i]); } ierr = DMDAVecRestoreArrayDOF(cda, coordinates, &coords);CHKERRQ(ierr); break; } case 2: { PetscScalar ***coords; ierr = DMDAVecGetArrayDOF(cda, coordinates, &coords);CHKERRQ(ierr); for(PetscInt j = info.ys; j < info.gys+info.gym; ++j) { for(PetscInt i = info.xs; i < info.gxs+info.gxm; ++i) { PetscInt globalV = j*M + i + numGlobalCells; coordSection->updatePoint(renumbering[globalV], coords[j][i]); } } ierr = DMDAVecRestoreArrayDOF(cda, coordinates, &coords);CHKERRQ(ierr); break; } case 3: { PetscScalar ****coords; ierr = DMDAVecGetArrayDOF(cda, coordinates, &coords);CHKERRQ(ierr); for(PetscInt k = info.zs; k < info.gzs+info.gzm; ++k, ++v) { for(PetscInt j = info.ys; j < info.gys+info.gym; ++j) { for(PetscInt i = info.xs; i < info.gxs+info.gxm; ++i) { PetscInt globalV = (k*N + j)*M + i + numGlobalCells; coordSection->updatePoint(renumbering[globalV], coords[k][j][i]); } } } ierr = DMDAVecRestoreArrayDOF(cda, coordinates, &coords);CHKERRQ(ierr); break; } default: SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid DMDA dimension %d", dim); } } /* Get overlap for interdomain communication */ { typedef PETSC_MESH_TYPE::point_type point_type; PETSc::Log::Event("CreateOverlap").begin(); ALE::Obj<PETSC_MESH_TYPE::send_overlap_type> sendParallelMeshOverlap = mesh->getSendOverlap(); ALE::Obj<PETSC_MESH_TYPE::recv_overlap_type> recvParallelMeshOverlap = mesh->getRecvOverlap(); // Can I figure this out in a nicer way? ALE::SetFromMap<std::map<point_type,point_type> > globalPoints(renumbering); ALE::OverlapBuilder<>::constructOverlap(globalPoints, renumbering, sendParallelMeshOverlap, recvParallelMeshOverlap); if (debug) { sendParallelMeshOverlap->view("Send Overlap"); recvParallelMeshOverlap->view("Recieve Overlap"); } mesh->setCalculatedOverlap(true); PETSc::Log::Event("CreateOverlap").end(); } PetscFunctionReturn(0); }