Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Efficient interior-only embedding #179

Merged
merged 8 commits into from
Feb 3, 2023
Merged
Show file tree
Hide file tree
Changes from 7 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 0 additions & 1 deletion pygeo/mphys/mphys_dvgeo.py
Original file line number Diff line number Diff line change
@@ -1,6 +1,5 @@
import openmdao.api as om
from .. import DVGeometry, DVConstraints
from ..constraints.baseConstraint import LinearConstraint

try:
from .. import DVGeometryVSP
Expand Down
1 change: 0 additions & 1 deletion pygeo/parameterization/DVGeoSketch.py
Original file line number Diff line number Diff line change
Expand Up @@ -12,7 +12,6 @@
from abc import abstractmethod
from collections import OrderedDict
import numpy as np
import copy
from mpi4py import MPI
from .BaseDVGeo import BaseDVGeometry
from .designVars import geoDVComposite
Expand Down
1 change: 0 additions & 1 deletion pygeo/parameterization/DVGeoVSP.py
Original file line number Diff line number Diff line change
Expand Up @@ -9,7 +9,6 @@
from .DVGeoSketch import DVGeoSketch
from pyspline.utils import searchQuads
from .designVars import vspDV
import copy

# openvsp python interface
try:
Expand Down
43 changes: 32 additions & 11 deletions pygeo/pyBlock.py
Original file line number Diff line number Diff line change
Expand Up @@ -6,6 +6,7 @@
import numpy as np
from scipy import sparse
from scipy.sparse import linalg
from scipy.spatial import ConvexHull
from pyspline import Volume
from pyspline.utils import openTecplot, writeTecplot3D, closeTecplot
from .geo_utils import readNValues, blendKnotVectors
Expand Down Expand Up @@ -784,7 +785,7 @@ def attachPoints(self, coordinates, ptSetName, interiorOnly=False, embTol=1e-10,
ptSetName : str
The name given to this set of coordinates.
interiorOnly : bool
Project only points that lie fully inside the volume
Only embed points that lie fully inside the volume
embTol : float
Tolerance on the distance between projected and closest point.
Determines if a point is embedded or not in the FFD volume if interiorOnly is True.
Expand All @@ -798,9 +799,8 @@ def attachPoints(self, coordinates, ptSetName, interiorOnly=False, embTol=1e-10,

# Project Points, if some were actually passed in:
if coordinates is not None:
checkErrors = not interiorOnly
marcomangano marked this conversation as resolved.
Show resolved Hide resolved
mask = None
volID, u, v, w, D = self.projectPoints(coordinates, checkErrors, embTol, eps, nIter)
volID, u, v, w, D = self.projectPoints(coordinates, interiorOnly, embTol, eps, nIter)

if interiorOnly:
# Create the mask before creating the embedded volume
Expand All @@ -817,7 +817,7 @@ def attachPoints(self, coordinates, ptSetName, interiorOnly=False, embTol=1e-10,
# Geometric Functions
# ----------------------------------------------------------------------

def projectPoints(self, x0, checkErrors, embTol, eps, nIter):
def projectPoints(self, x0, interiorOnly, embTol, eps, nIter):
"""Project a set of points x0, into any one of the volumes. It
returns the the volume ID, u, v, w, D of the point in volID or
closest to it.
Expand All @@ -834,9 +834,6 @@ def projectPoints(self, x0, checkErrors, embTol, eps, nIter):
----------
x0 : array of points (Nx3 array)
The list or array of points to use
checkErrors : bool
Flag to print out the error is points have not been projected
to the tolerance defined by embTol.

See Also
--------
Expand All @@ -858,7 +855,32 @@ def projectPoints(self, x0, checkErrors, embTol, eps, nIter):
v0 = 0.0
w0 = 0.0

# If we are only interested in interior points, we skip projecting exterior points to save time.
# We identify exterior points by checking if they are outside the convex hull of the control points.
# A point can be inside the convex hull but still outside the FFD volume(s).
# In this case, we rely on the more costly projection approach to identify the exterior points.
if interiorOnly:
# Compute the convex hull of all control points
hull = ConvexHull(self.coef)

# ConvexHull.equations describes the planes that define the faces of the convex hull
# Extract the normal vector and offset for each plane
hullNormals = hull.equations[:, :-1]
hullOffsets = hull.equations[:, -1]

# The normals point outside the convex hull, so a point is inside the convex hull if the distance in the
# normal direction from the point to every plane defining the convex hull is negative.
# This is computed in a vectorized manner below.
# The offset is negative in the normal direction, so we add the offset instead of subtracting.
distanceToPlanes = np.dot(x0, hullNormals.T) + hullOffsets
sseraj marked this conversation as resolved.
Show resolved Hide resolved
isInsideHull = np.all(distanceToPlanes <= eps, axis=1)

for i in range(N):

# Do not project this point if it is outside the convex hull and we are only interested in interior points
if interiorOnly and not isInsideHull[i]:
continue

for j in range(self.nVol):
iVol = volList[j]
u0, v0, w0, D0 = self.vols[iVol].projectPoint(x0[i], eps=eps, nIter=nIter)
Expand Down Expand Up @@ -887,10 +909,9 @@ def projectPoints(self, x0, checkErrors, embTol, eps, nIter):
volList = np.hstack([iVol, volList[:j], volList[j + 1 :]])
# end for (length of x0)

# If desired check the errors and print warnings:
if checkErrors:
# Loop back through the points and determine which ones are
# bad (> 50*eps) and print them to the screen:
# If we are interested in all points, we need to check whether they were all projected properly
if not interiorOnly:
# Loop back through the points and determine which ones are bad and print them to the screen
counter = 0
DMax = 0.0
DRms = 0.0
Expand Down