Experiment separating out Ptex from TopologyRefiner.

Remove the ptex-specific code from the Far::TopologyRefiner and instead provide it in a separate class Far::PtexIndices.  Clients who need to use the Ptex API need to first build a Far::PtexIndices object by providing it with a refiner.

This has the advantage of keeping the API on the TopologyRefiner a little cleaner.  The ptex methods were const but were mutating state with const_casts.  The new mechanism still achieves the same lazy initialization behavior by forcing clients to instantiate them exactly when needed.

A disadvantage of this approach is that the PatchTablesFactory creates its own PtexIndices and throws it out after the patch tables are created.  This is great if you're never going to need the ptex indices again, but not so great if you will need them again.
This commit is contained in:
George ElKoura 2015-04-19 23:25:43 -07:00
parent 7e3bd65432
commit 420473b45b
12 changed files with 355 additions and 215 deletions

View File

@ -24,12 +24,16 @@
#include "particles.h"
#include <far/ptexIndices.h>
#include <cassert>
STParticles::STParticles(Refiner const & refiner, int nparticles, bool centered) :
_speed(1.0f) {
int nptexfaces = refiner.GetNumPtexFaces(),
OpenSubdiv::Far::PtexIndices ptexIndices(refiner);
int nptexfaces = ptexIndices.GetNumPtexFaces(),
nsamples = nptexfaces * nparticles;
srand(static_cast<int>(2147483647));
@ -76,12 +80,12 @@ STParticles::STParticles(Refiner const & refiner, int nparticles, bool centered)
refiner.GetFaceVertices(0, face);
if (fverts.size()==4) {
refiner.GetPtexAdjacency(face, 0, adjfaces, adjedges);
ptexIndices.GetPtexAdjacency(refiner, face, 0, adjfaces, adjedges);
_adjacency[ptexface] = FaceInfo(adjfaces, adjedges, false);
++ptexface;
} else {
for (int vert=0; vert<fverts.size(); ++vert) {
refiner.GetPtexAdjacency(face, vert, adjfaces, adjedges);
ptexIndices.GetPtexAdjacency(refiner, face, vert, adjfaces, adjedges);
_adjacency[ptexface+vert] =
FaceInfo(adjfaces, adjedges, true);
}

View File

@ -46,6 +46,7 @@ GLFWmonitor* g_primary=0;
#include <osd/glDrawContext.h>
#include <osd/glDrawRegistry.h>
#include <far/error.h>
#include <far/ptexIndices.h>
#include <osd/cpuGLVertexBuffer.h>
#include <osd/cpuComputeContext.h>
@ -229,7 +230,8 @@ createOsdMesh() {
OpenSubdiv::Far::TopologyRefinerFactory<Shape>::Options(sdctype, sdcoptions));
// count ptex face id
int numPtexFaces = refiner->GetNumPtexFaces();
OpenSubdiv::Far::PtexIndices ptexIndices(*refiner);
int numPtexFaces = ptexIndices.GetNumPtexFaces();
delete g_mesh;
g_mesh = NULL;

View File

@ -47,6 +47,7 @@ GLFWmonitor* g_primary=0;
#include <osd/glDrawRegistry.h>
#include <osd/glMesh.h>
#include <far/error.h>
#include <far/ptexIndices.h>
#include <osd/cpuGLVertexBuffer.h>
#include <osd/cpuComputeContext.h>
@ -588,7 +589,8 @@ createOsdMesh( const std::string &shapeStr, int level, Scheme scheme=kCatmark )
}
// create ptex index to coarse face index mapping
int numPtexFaces = refiner->GetNumPtexFaces();
Far::PtexIndices ptexIndices(*refiner);
int numPtexFaces = ptexIndices.GetNumPtexFaces();
// XXX: duped logic to simpleHbr
std::vector<int> ptexIndexToFaceMapping(numPtexFaces);

View File

@ -49,6 +49,7 @@ GLFWmonitor* g_primary=0;
#include "../common/gl_hud.h"
#include <far/patchTablesFactory.h>
#include <far/ptexIndices.h>
#include <far/stencilTablesFactory.h>
#include <osd/cpuGLVertexBuffer.h>
@ -327,7 +328,8 @@ createMesh(ShapeDesc const & shapeDesc, int level) {
refiner->RefineAdaptive(options);
}
int nfaces = refiner->GetNumPtexFaces();
Far::PtexIndices ptexIndices(*refiner);
int nfaces = ptexIndices.GetNumPtexFaces();
float * u = new float[g_nsamples*nfaces], * uPtr = u,
* v = new float[g_nsamples*nfaces], * vPtr = v;

View File

@ -35,6 +35,7 @@ set(SOURCE_FILES
patchMap.cpp
patchTables.cpp
patchTablesFactory.cpp
ptexIndices.cpp
stencilTablesFactory.cpp
topologyRefiner.cpp
topologyRefinerFactory.cpp
@ -56,6 +57,7 @@ set(PUBLIC_HEADER_FILES
patchMap.h
patchTables.h
patchTablesFactory.h
ptexIndices.h
stencilTables.h
stencilTablesFactory.h
topologyRefiner.h

View File

@ -23,6 +23,7 @@
//
#include "../far/patchTablesFactory.h"
#include "../far/error.h"
#include "../far/ptexIndices.h"
#include "../far/topologyRefiner.h"
#include "../vtr/level.h"
#include "../vtr/refinement.h"
@ -639,8 +640,9 @@ PatchTablesFactoryBase::gatherFVarData(AdaptiveContext & context, int level,
//
PatchParam *
PatchTablesFactoryBase::computePatchParam(
TopologyRefiner const & refiner, int depth, Vtr::Index faceIndex,
int rotation, int boundaryMask, int transitionMask, PatchParam *coord) {
TopologyRefiner const & refiner, PtexIndices const &ptexIndices,
int depth, Vtr::Index faceIndex, int rotation, int boundaryMask,
int transitionMask, PatchParam *coord) {
if (coord == NULL) return NULL;
@ -682,7 +684,7 @@ PatchTablesFactoryBase::computePatchParam(
faceIndex = parentFaceIndex;
}
Vtr::Index ptexIndex = refiner.GetPtexIndex(faceIndex);
Vtr::Index ptexIndex = ptexIndices.GetPtexIndex(faceIndex);
assert(ptexIndex!=-1);
if (nonquad) {
@ -743,6 +745,8 @@ PatchTablesFactoryBase::createUniform(TopologyRefiner const & refiner, Options o
firstlevel = options.generateAllLevels ? 0 : maxlevel,
nlevels = maxlevel-firstlevel+1;
PtexIndices ptexIndices(refiner);
PatchDescriptor::Type ptype = PatchDescriptor::NON_PATCH;
if (options.triangulateQuads) {
ptype = PatchDescriptor::TRIANGLES;
@ -760,7 +764,7 @@ PatchTablesFactoryBase::createUniform(TopologyRefiner const & refiner, Options o
//
PatchTables * tables = new PatchTables(maxvalence);
tables->_numPtexFaces = refiner.GetNumPtexFaces();
tables->_numPtexFaces = ptexIndices.GetNumPtexFaces();
tables->reservePatchArrays(nlevels);
@ -832,7 +836,7 @@ PatchTablesFactoryBase::createUniform(TopologyRefiner const & refiner, Options o
*iptr++ = levelVertOffset + fverts[vert];
}
pptr = computePatchParam(refiner, level, face, /*rot*/0, /*boundary*/0, /*transition*/0, pptr);
pptr = computePatchParam(refiner, ptexIndices, level, face, /*rot*/0, /*boundary*/0, /*transition*/0, pptr);
if (generateFVarPatches) {
for (fvc=fvc.begin(); fvc!=fvc.end(); ++fvc) {
@ -901,6 +905,8 @@ PatchTablesFactoryT<ENDCAP_FACTORY>::createAdaptive(TopologyRefiner const & refi
assert(not refiner.IsUniform());
PtexIndices ptexIndices(refiner);
AdaptiveContext context(refiner, options);
//
@ -932,7 +938,7 @@ PatchTablesFactoryT<ENDCAP_FACTORY>::createAdaptive(TopologyRefiner const & refi
context.patchInventory.getValue(desc), &voffset, &poffset, &qoffset );
}
context.tables->_numPtexFaces = refiner.GetNumPtexFaces();
context.tables->_numPtexFaces = ptexIndices.GetNumPtexFaces();
// Allocate various tables
bool hasSharpness = context.options.useSingleCreasePatch;
@ -954,7 +960,7 @@ PatchTablesFactoryT<ENDCAP_FACTORY>::createAdaptive(TopologyRefiner const & refi
//
// Now populate the patches:
//
populateAdaptivePatches(context, endCapFactory);
populateAdaptivePatches(context, ptexIndices, endCapFactory);
return context.tables;
}
@ -1228,8 +1234,9 @@ PatchTablesFactoryT<ENDCAP_FACTORY>::identifyAdaptivePatches(
//
template<typename ENDCAP_FACTORY>
void
PatchTablesFactoryT<ENDCAP_FACTORY>::populateAdaptivePatches(AdaptiveContext & context,
ENDCAP_FACTORY *endCapFactory) {
PatchTablesFactoryT<ENDCAP_FACTORY>::populateAdaptivePatches(
AdaptiveContext & context, PtexIndices const & ptexIndices,
ENDCAP_FACTORY *endCapFactory) {
TopologyRefiner const & refiner = context.refiner;
@ -1326,7 +1333,7 @@ PatchTablesFactoryT<ENDCAP_FACTORY>::populateAdaptivePatches(AdaptiveContext & c
offsetAndPermuteIndices(patchVerts, 16, levelVertOffset, permuteInterior, iptrs.R);
iptrs.R += 16;
pptrs.R = computePatchParam(refiner, i, faceIndex, /*rotation*/0, /*boundary*/0, transitionMask, pptrs.R);
pptrs.R = computePatchParam(refiner, ptexIndices, i, faceIndex, /*rotation*/0, /*boundary*/0, transitionMask, pptrs.R);
// XXX: sharpness will be integrated into patch param soon.
if (sptrs.R) *sptrs.R++ = assignSharpnessIndex(0, tables->_sharpnessValues);
@ -1362,7 +1369,7 @@ PatchTablesFactoryT<ENDCAP_FACTORY>::populateAdaptivePatches(AdaptiveContext & c
sharpness = std::min(sharpness, (float)(context.options.maxIsolationLevel-i));
iptrs.R += 16;
pptrs.R = computePatchParam(refiner, i, faceIndex, bIndex, /*boundary*/0, transitionMask, pptrs.R);
pptrs.R = computePatchParam(refiner, ptexIndices, i, faceIndex, bIndex, /*boundary*/0, transitionMask, pptrs.R);
if (sptrs.R) *sptrs.R++ = assignSharpnessIndex(sharpness, tables->_sharpnessValues);
fofss.R += gatherFVarData(context,
@ -1376,7 +1383,7 @@ PatchTablesFactoryT<ENDCAP_FACTORY>::populateAdaptivePatches(AdaptiveContext & c
bIndex = (bIndex+1)%4;
iptrs.R += 16;
pptrs.R = computePatchParam(refiner, i, faceIndex, bIndex, boundaryMask, transitionMask, pptrs.R);
pptrs.R = computePatchParam(refiner, ptexIndices, i, faceIndex, bIndex, boundaryMask, transitionMask, pptrs.R);
if (sptrs.R) *sptrs.R++ = assignSharpnessIndex(0, tables->_sharpnessValues);
@ -1389,7 +1396,7 @@ PatchTablesFactoryT<ENDCAP_FACTORY>::populateAdaptivePatches(AdaptiveContext & c
offsetAndPermuteIndices(patchVerts, 16, levelVertOffset, permuteCorner, iptrs.R);
iptrs.R += 16;
pptrs.R = computePatchParam(refiner, i, faceIndex, bIndex, boundaryMask, transitionMask, pptrs.R);
pptrs.R = computePatchParam(refiner, ptexIndices, i, faceIndex, bIndex, boundaryMask, transitionMask, pptrs.R);
if (sptrs.R) *sptrs.R++ = assignSharpnessIndex(0, tables->_sharpnessValues);
@ -1417,7 +1424,7 @@ PatchTablesFactoryT<ENDCAP_FACTORY>::populateAdaptivePatches(AdaptiveContext & c
for (int j = 0; j < cvs.size(); ++j) iptrs.R[j] = cvs[j];
iptrs.R += cvs.size();
pptrs.R = computePatchParam(
refiner, i, faceIndex, 0, /*boundary*/0, /*transition*/0, pptrs.R);
refiner, ptexIndices, i, faceIndex, 0, /*boundary*/0, /*transition*/0, pptrs.R);
if (sptrs.R) *sptrs.R++ = assignSharpnessIndex(0, tables->_sharpnessValues);
fofss.R += gatherFVarData(context,
i, faceIndex, levelFaceOffset,
@ -1429,7 +1436,7 @@ PatchTablesFactoryT<ENDCAP_FACTORY>::populateAdaptivePatches(AdaptiveContext & c
for (int j = 0; j < cvs.size(); ++j) iptrs.GP[j] = cvs[j];
iptrs.GP += cvs.size();
pptrs.GP = computePatchParam(
refiner, i, faceIndex, 0, /*boundary*/0, /*transition*/0, pptrs.GP);
refiner, ptexIndices, i, faceIndex, 0, /*boundary*/0, /*transition*/0, pptrs.GP);
if (sptrs.R) *sptrs.R++ = assignSharpnessIndex(0, tables->_sharpnessValues);
fofss.GP += gatherFVarData(context,
i, faceIndex, levelFaceOffset,
@ -1441,7 +1448,7 @@ PatchTablesFactoryT<ENDCAP_FACTORY>::populateAdaptivePatches(AdaptiveContext & c
for (int j = 0; j < cvs.size(); ++j) iptrs.G[j] = cvs[j];
iptrs.G += cvs.size();
pptrs.G = computePatchParam(
refiner, i, faceIndex, 0, /*boundary*/0, /*transition*/0, pptrs.G);
refiner, ptexIndices, i, faceIndex, 0, /*boundary*/0, /*transition*/0, pptrs.G);
if (sptrs.R) *sptrs.R++ = assignSharpnessIndex(0, tables->_sharpnessValues);
fofss.G += gatherFVarData(context,
i, faceIndex, levelFaceOffset,
@ -1453,7 +1460,7 @@ PatchTablesFactoryT<ENDCAP_FACTORY>::populateAdaptivePatches(AdaptiveContext & c
for (int j = 0; j < cvs.size(); ++j) iptrs.GB[j] = cvs[j];
iptrs.GB += cvs.size();
pptrs.GB = computePatchParam(
refiner, i, faceIndex, 0, /*boundary*/0, /*transition*/0, pptrs.GB);
refiner, ptexIndices, i, faceIndex, 0, /*boundary*/0, /*transition*/0, pptrs.GB);
if (sptrs.R) *sptrs.R++ = assignSharpnessIndex(0, tables->_sharpnessValues);
fofss.GB += gatherFVarData(context,
i, faceIndex, levelFaceOffset,

View File

@ -37,6 +37,7 @@ namespace Vtr { class Level; }
namespace Far {
class PtexIndices;
class TopologyRefiner;
class PatchTablesFactoryBase {
@ -130,6 +131,7 @@ protected:
Options options, int npatches, PatchTables * tables);
static PatchParam * computePatchParam(TopologyRefiner const & refiner,
PtexIndices const & ptexIndices,
int level, int face, int rotation,
int boundaryMask, int transitionMask, PatchParam * coord);
@ -172,7 +174,8 @@ private:
static void identifyAdaptivePatches(AdaptiveContext & state, ENDCAP_FACTORY *);
static void populateAdaptivePatches(AdaptiveContext & state, ENDCAP_FACTORY *);
static void populateAdaptivePatches(AdaptiveContext & state,
PtexIndices const & ptexIndices, ENDCAP_FACTORY *);
private:
};

View File

@ -0,0 +1,202 @@
//
// Copyright 2014 DreamWorks Animation LLC.
//
// Licensed under the Apache License, Version 2.0 (the "Apache License")
// with the following modification; you may not use this file except in
// compliance with the Apache License and the following modification to it:
// Section 6. Trademarks. is deleted and replaced with:
//
// 6. Trademarks. This License does not grant permission to use the trade
// names, trademarks, service marks, or product names of the Licensor
// and its affiliates, except as required to comply with Section 4(c) of
// the License and to reproduce the content of the NOTICE file.
//
// You may obtain a copy of the Apache License at
//
// http://www.apache.org/licenses/LICENSE-2.0
//
// Unless required by applicable law or agreed to in writing, software
// distributed under the Apache License with the above modification is
// distributed on an "AS IS" BASIS, WITHOUT WARRANTIES OR CONDITIONS OF ANY
// KIND, either express or implied. See the Apache License for the specific
// language governing permissions and limitations under the Apache License.
//
#include "../far/ptexIndices.h"
#include "../vtr/level.h"
#include <cassert>
namespace OpenSubdiv {
namespace OPENSUBDIV_VERSION {
namespace Far {
PtexIndices::PtexIndices(TopologyRefiner const &refiner) {
initializePtexIndices(refiner);
}
PtexIndices::~PtexIndices() {
}
void
PtexIndices::initializePtexIndices(TopologyRefiner const &refiner) {
Vtr::Level const & coarseLevel = refiner.getLevel(0);
int nfaces = coarseLevel.getNumFaces();
_ptexIndices.resize(nfaces+1);
int ptexID=0;
int regFaceSize = Sdc::SchemeTypeTraits::GetRegularFaceSize(
refiner.GetSchemeType());
for (int i = 0; i < nfaces; ++i) {
_ptexIndices[i] = ptexID;
Vtr::ConstIndexArray fverts = coarseLevel.getFaceVertices(i);
ptexID += fverts.size()==regFaceSize ? 1 : fverts.size();
}
// last entry contains the number of ptex texture faces
_ptexIndices[nfaces]=ptexID;
}
int
PtexIndices::GetNumPtexFaces() const {
return _ptexIndices.back();
}
int
PtexIndices::GetPtexIndex(Index f) const {
assert(f<(int)_ptexIndices.size());
return _ptexIndices[f];
}
namespace {
// Returns the face adjacent to 'face' along edge 'edge'
inline Index
getAdjacentFace(Vtr::Level const & level, Index edge, Index face) {
Far::ConstIndexArray adjFaces = level.getEdgeFaces(edge);
if (adjFaces.size()!=2) {
return -1;
}
return (adjFaces[0]==face) ? adjFaces[1] : adjFaces[0];
}
}
void
PtexIndices::GetPtexAdjacency(
TopologyRefiner const &refiner,
int face, int quadrant,
int adjFaces[4], int adjEdges[4]) const {
assert(refiner.GetSchemeType()==Sdc::SCHEME_CATMARK);
Vtr::Level const & level = refiner.getLevel(0);
ConstIndexArray fedges = level.getFaceEdges(face);
if (fedges.size()==4) {
// Regular ptex quad face
for (int i=0; i<4; ++i) {
int edge = fedges[i];
Index adjface = getAdjacentFace(level, edge, face);
if (adjface==-1) {
adjFaces[i] = -1; // boundary or non-manifold
adjEdges[i] = 0;
} else {
ConstIndexArray aedges = level.getFaceEdges(adjface);
if (aedges.size()==4) {
adjFaces[i] = _ptexIndices[adjface];
adjEdges[i] = aedges.FindIndexIn4Tuple(edge);
assert(adjEdges[i]!=-1);
} else {
// neighbor is a sub-face
adjFaces[i] = _ptexIndices[adjface] +
(aedges.FindIndex(edge)+1)%aedges.size();
adjEdges[i] = 3;
}
assert(adjFaces[i]!=-1);
}
}
} else {
// Ptex sub-face 'quadrant' (non-quad)
//
// Ptex adjacency pattern for non-quads:
//
// v2
/* o
// / \
// / \
// /0 3\
// / \
// o_ 1 2 _o
// / -_ _- \
// / 2 -o- 1 \
// /3 | 0\
// / 1|2 \
// / 0 | 3 \
// o----------o----------o
// v0 v1
*/
assert(quadrant>=0 and quadrant<fedges.size());
int nextQuadrant = (quadrant+1) % fedges.size(),
prevQuadrant = (quadrant+fedges.size()-1) % fedges.size();
{ // resolve neighbors within the sub-face (edges 1 & 2)
adjFaces[1] = _ptexIndices[face] + nextQuadrant;
adjEdges[1] = 2;
adjFaces[2] = _ptexIndices[face] + prevQuadrant;
adjEdges[2] = 1;
}
{ // resolve neighbor outisde the sub-face (edge 0)
int edge0 = fedges[quadrant];
Index adjface0 = getAdjacentFace(level, edge0, face);
if (adjface0==-1) {
adjFaces[0] = -1; // boundary or non-manifold
adjEdges[0] = 0;
} else {
ConstIndexArray afedges = level.getFaceEdges(adjface0);
if (afedges.size()==4) {
adjFaces[0] = _ptexIndices[adjface0];
adjEdges[0] = afedges.FindIndexIn4Tuple(edge0);
} else {
int subedge = (afedges.FindIndex(edge0)+1)%afedges.size();
adjFaces[0] = _ptexIndices[adjface0] + subedge;
adjEdges[0] = 3;
}
assert(adjFaces[0]!=-1);
}
// resolve neighbor outisde the sub-face (edge 3)
int edge3 = fedges[prevQuadrant];
Index adjface3 = getAdjacentFace(level, edge3, face);
if (adjface3==-1) {
adjFaces[3]=-1; // boundary or non-manifold
adjEdges[3]=0;
} else {
ConstIndexArray afedges = level.getFaceEdges(adjface3);
if (afedges.size()==4) {
adjFaces[3] = _ptexIndices[adjface3];
adjEdges[3] = afedges.FindIndexIn4Tuple(edge3);
} else {
int subedge = afedges.FindIndex(edge3);
adjFaces[3] = _ptexIndices[adjface3] + subedge;
adjEdges[3] = 0;
}
assert(adjFaces[3]!=-1);
}
}
}
}
} // end namespace Far
} // end namespace OPENSUBDIV_VERSION
} // end namespace OpenSubdiv

View File

@ -0,0 +1,102 @@
//
// Copyright 2014 DreamWorks Animation LLC.
//
// Licensed under the Apache License, Version 2.0 (the "Apache License")
// with the following modification; you may not use this file except in
// compliance with the Apache License and the following modification to it:
// Section 6. Trademarks. is deleted and replaced with:
//
// 6. Trademarks. This License does not grant permission to use the trade
// names, trademarks, service marks, or product names of the Licensor
// and its affiliates, except as required to comply with Section 4(c) of
// the License and to reproduce the content of the NOTICE file.
//
// You may obtain a copy of the Apache License at
//
// http://www.apache.org/licenses/LICENSE-2.0
//
// Unless required by applicable law or agreed to in writing, software
// distributed under the Apache License with the above modification is
// distributed on an "AS IS" BASIS, WITHOUT WARRANTIES OR CONDITIONS OF ANY
// KIND, either express or implied. See the Apache License for the specific
// language governing permissions and limitations under the Apache License.
//
#ifndef FAR_PTEX_INDICES_H
#define FAR_PTEX_INDICES_H
#include "../version.h"
#include "../far/topologyRefiner.h"
#include "../far/types.h"
#include <vector>
namespace OpenSubdiv {
namespace OPENSUBDIV_VERSION {
namespace Far {
///
/// \brief
///
class PtexIndices {
public:
/// \brief Constructor
PtexIndices(TopologyRefiner const &refiner);
/// \brief Destructor
~PtexIndices();
//@{
///
/// Ptex
///
/// \brief Returns the number of ptex faces in the mesh
///
int GetNumPtexFaces() const;
/// \brief Returns the ptex face index given a coarse face 'f' or -1
///
int GetPtexIndex(Index f) const;
/// \brief Returns ptex face adjacency information for a given coarse face
///
/// @param refiner refiner used to build this PtexIndices object.
///
/// @param face coarse face index
///
/// @param quadrant quadrant index if 'face' is not a quad (the local ptex
/// sub-face index). Must be less than the number of face
/// vertices.
///
/// @param adjFaces ptex face indices of adjacent faces
///
/// @param adjEdges ptex edge indices of adjacent faces
///
void GetPtexAdjacency(
TopologyRefiner const &refiner,
int face, int quadrant,
int adjFaces[4], int adjEdges[4]) const;
//@}
private:
void initializePtexIndices(TopologyRefiner const &refiner);
private:
std::vector<Index> _ptexIndices;
};
} // end namespace Far
} // end namespace OPENSUBDIV_VERSION
using namespace OPENSUBDIV_VERSION;
} // end namespace OpenSubdiv
#endif /* FAR_PTEX_INDICES_H */

View File

@ -175,167 +175,6 @@ TopologyRefiner::GetNumHoles(int level) const {
return sum;
}
//
// Ptex information accessors
//
void
TopologyRefiner::initializePtexIndices() const {
Vtr::Level const & coarseLevel = getLevel(0);
std::vector<int> & ptexIndices = const_cast<std::vector<int> &>(_ptexIndices);
int nfaces = coarseLevel.getNumFaces();
ptexIndices.resize(nfaces+1);
int ptexID=0;
int regFaceSize = Sdc::SchemeTypeTraits::GetRegularFaceSize(GetSchemeType());
for (int i = 0; i < nfaces; ++i) {
ptexIndices[i] = ptexID;
Vtr::ConstIndexArray fverts = coarseLevel.getFaceVertices(i);
ptexID += fverts.size()==regFaceSize ? 1 : fverts.size();
}
// last entry contains the number of ptex texture faces
ptexIndices[nfaces]=ptexID;
}
int
TopologyRefiner::GetNumPtexFaces() const {
if (_ptexIndices.empty()) {
initializePtexIndices();
}
return _ptexIndices.back();
}
int
TopologyRefiner::GetPtexIndex(Index f) const {
if (_ptexIndices.empty()) {
initializePtexIndices();
}
assert(f<(int)_ptexIndices.size());
return _ptexIndices[f];
}
namespace {
// Returns the face adjacent to 'face' along edge 'edge'
inline Index
getAdjacentFace(Vtr::Level const & level, Index edge, Index face) {
Far::ConstIndexArray adjFaces = level.getEdgeFaces(edge);
if (adjFaces.size()!=2) {
return -1;
}
return (adjFaces[0]==face) ? adjFaces[1] : adjFaces[0];
}
}
void
TopologyRefiner::GetPtexAdjacency(int face, int quadrant,
int adjFaces[4], int adjEdges[4]) const {
assert(GetSchemeType()==Sdc::SCHEME_CATMARK);
if (_ptexIndices.empty()) {
initializePtexIndices();
}
Vtr::Level const & level = getLevel(0);
ConstIndexArray fedges = level.getFaceEdges(face);
if (fedges.size()==4) {
// Regular ptex quad face
for (int i=0; i<4; ++i) {
int edge = fedges[i];
Index adjface = getAdjacentFace(level, edge, face);
if (adjface==-1) {
adjFaces[i] = -1; // boundary or non-manifold
adjEdges[i] = 0;
} else {
ConstIndexArray aedges = level.getFaceEdges(adjface);
if (aedges.size()==4) {
adjFaces[i] = _ptexIndices[adjface];
adjEdges[i] = aedges.FindIndexIn4Tuple(edge);
assert(adjEdges[i]!=-1);
} else {
// neighbor is a sub-face
adjFaces[i] = _ptexIndices[adjface] +
(aedges.FindIndex(edge)+1)%aedges.size();
adjEdges[i] = 3;
}
assert(adjFaces[i]!=-1);
}
}
} else {
// Ptex sub-face 'quadrant' (non-quad)
//
// Ptex adjacency pattern for non-quads:
//
// v2
/* o
// / \
// / \
// /0 3\
// / \
// o_ 1 2 _o
// / -_ _- \
// / 2 -o- 1 \
// /3 | 0\
// / 1|2 \
// / 0 | 3 \
// o----------o----------o
// v0 v1
*/
assert(quadrant>=0 and quadrant<fedges.size());
int nextQuadrant = (quadrant+1) % fedges.size(),
prevQuadrant = (quadrant+fedges.size()-1) % fedges.size();
{ // resolve neighbors within the sub-face (edges 1 & 2)
adjFaces[1] = _ptexIndices[face] + nextQuadrant;
adjEdges[1] = 2;
adjFaces[2] = _ptexIndices[face] + prevQuadrant;
adjEdges[2] = 1;
}
{ // resolve neighbor outisde the sub-face (edge 0)
int edge0 = fedges[quadrant];
Index adjface0 = getAdjacentFace(level, edge0, face);
if (adjface0==-1) {
adjFaces[0] = -1; // boundary or non-manifold
adjEdges[0] = 0;
} else {
ConstIndexArray afedges = level.getFaceEdges(adjface0);
if (afedges.size()==4) {
adjFaces[0] = _ptexIndices[adjface0];
adjEdges[0] = afedges.FindIndexIn4Tuple(edge0);
} else {
int subedge = (afedges.FindIndex(edge0)+1)%afedges.size();
adjFaces[0] = _ptexIndices[adjface0] + subedge;
adjEdges[0] = 3;
}
assert(adjFaces[0]!=-1);
}
// resolve neighbor outisde the sub-face (edge 3)
int edge3 = fedges[prevQuadrant];
Index adjface3 = getAdjacentFace(level, edge3, face);
if (adjface3==-1) {
adjFaces[3]=-1; // boundary or non-manifold
adjEdges[3]=0;
} else {
ConstIndexArray afedges = level.getFaceEdges(adjface3);
if (afedges.size()==4) {
adjFaces[3] = _ptexIndices[adjface3];
adjEdges[3] = afedges.FindIndexIn4Tuple(edge3);
} else {
int subedge = afedges.FindIndex(edge3);
adjFaces[3] = _ptexIndices[adjface3] + subedge;
adjEdges[3] = 0;
}
assert(adjFaces[3]!=-1);
}
}
}
}
//
// Main refinement method -- allocating and initializing levels and refinements:

View File

@ -488,32 +488,6 @@ public:
//@}
//@{
/// Ptex
///
/// \brief Returns the number of ptex faces in the mesh
int GetNumPtexFaces() const;
/// \brief Returns the ptex face index given a coarse face 'f' or -1
int GetPtexIndex(Index f) const;
/// \brief Returns ptex face adjacency information for a given coarse face
///
/// @param face coarse face index
///
/// @param quadrant quadrant index if 'face' is not a quad (the local ptex
// sub-face index). Must be less than the number of face
// vertices.
///
/// @param adjFaces ptex face indices of adjacent faces
///
/// @param adjEdges ptex edge indices of adjacent faces
///
void GetPtexAdjacency(int face, int quadrant,
int adjFaces[4], int adjEdges[4]) const;
//@}
//@{
/// @name Debugging aides
@ -592,6 +566,7 @@ protected:
friend class EndCapLegacyGregoryPatchFactory;
friend class EndCapGregoryBasisPatchFactory;
friend class EndCapRegularPatchFactory;
friend class PtexIndices;
Vtr::Level & getLevel(int l) { return *_levels[l]; }
Vtr::Level const & getLevel(int l) const { return *_levels[l]; }
@ -618,7 +593,6 @@ private:
template <Sdc::SchemeType SCHEME, class T, class U> void faceVaryingLimit(T const & src, U * dst, int channel) const;
void initializePtexIndices() const;
void initializeInventory();
void updateInventory(Vtr::Level const & newLevel);
@ -648,8 +622,6 @@ private:
std::vector<Vtr::Level *> _levels;
std::vector<Vtr::Refinement *> _refinements;
std::vector<Index> _ptexIndices;
};

View File

@ -43,6 +43,7 @@
#include <opensubdiv/far/topologyRefinerFactory.h>
#include <opensubdiv/far/patchTablesFactory.h>
#include <opensubdiv/far/patchMap.h>
#include <opensubdiv/far/ptexIndices.h>
#include <cassert>
#include <cstdio>
@ -160,10 +161,12 @@ int main(int, char **) {
// Create a Far::PatchMap to help locating patches in the table
Far::PatchMap patchmap(*patchTables);
// Create a Far::PtexIndices to help find indices of ptex faces.
Far::PtexIndices ptexIndices(*refiner);
// Generate random samples on each ptex face
int nsamples = 200,
nfaces = refiner->GetNumPtexFaces();
nfaces = ptexIndices.GetNumPtexFaces();
std::vector<LimitFrame> samples(nsamples * nfaces);