raytraverse.lightpoint¶
LightPointKD¶
- class raytraverse.lightpoint.LightPointKD(scene, vec=None, lum=None, vm=None, pt=(0, 0, 0), posidx=0, src='sky', srcn=1, srcdir=(0, 0, 1), calcomega=True, write=True, omega=None, filterviews=True, srcviews=None, parent=None, srcviewidxs=None, features=1)[source]¶
Bases:
object
light distribution from a point with KDtree structure for directional query
- Parameters:
scene (raytraverse.scene.BaseScene) –
vec (np.array, optional) – shape (N, >=3) where last three columns are normalized direction vectors of samples. If not given, tries to load from scene.outdir
lum (np.array, optional) – reshapeable to (N, srcn). sample values for each source corresponding to vec. If not given, tries to load from scene.outdir
vm (raytraverse.mapper.ViewMapper, optional) – a default viewmapper for image and metric calculations, should match viewmapper of sampler.run() if possible.
pt (tuple list np.array) – 3 item point location of light distribution
posidx (int, optional) – index position of point, will govern file naming so must be set to avoid clobbering writes. also used by spacemapper for planar sampling
src (str, optional) – name of source group. will govern file naming so must be set to avoid clobbering writes.
srcn (int, optional) – must match lum, does not need to be set if reloading from scene.outdir
calcomega (bool, optional) – if True (default) calculate solid angle of rays. This is unnecessary if point will be combined before calculating any metrics. setting to False will save some computation time.
write (bool, optional) – whether to save ray data to disk.
omega (np.array, optional) – provide precomputed omega values, if given, overrides calcomega
- vm¶
raytraverse.mapper.ViewMapper
- scene¶
raytraverse.scene.Scene
- posidx¶
index for point
- Type:
int
- pt¶
point location
- Type:
np.array
- src¶
source key
- Type:
str
- file¶
relative path to disk storage
- Type:
str
- srcdir¶
direction to source(s)
- property vec¶
direction vector (N,3)
- property lum¶
luminance (N,srcn)
- property d_kd¶
kd tree for spatial query
- Getter:
Returns kd tree structure
- Type:
scipy.spatial.cKDTree
- property omega¶
solid angle (N)
- Getter:
Returns array of solid angles
- Setter:
sets soolid angles with viewmapper
- Type:
np.array
- calc_omega(write=True)[source]¶
calculate solid angle
- Parameters:
write (bool, optional) – update/write kdtree data to file
- apply_coef(coefs)[source]¶
apply coefficient vector to self.lum
- Parameters:
coefs (np.array int float list) – shape (N, self.srcn) or broadcastable
- Returns:
alum – shape (N, self.vec.shape[0])
- Return type:
np.array
- add_to_img(img, vecs, mask=None, skyvec=1, interp=False, idx=None, interpweights=None, order=False, omega=False, vm=None, rnd=False, engine=None, srcrnd=False, **kwargs)[source]¶
add luminance contributions to image array (updates in place)
- Parameters:
img (np.array) – 2D image array to add to (either zeros or with other source)
vecs (np.array) – vectors corresponding to img pixels shape (N, 3)
mask (np.array, optional) – indices to img that correspond to vec (in case where whole image is not being updated, such as corners of fisheye)
skyvec (int float np.array, optional) – source coefficients, shape is (1,) or (srcn,)
interp (Union[bool, str], optional) –
if “precomp”, use index and interpweights
if True and engine is None, linearinterpolation
if “fastc” and engine: uses content_interp (best after sampling w/o detail)
if “highc” and engine: uses content_interp_wedge (best after sampling w/o detail)
if “fast”: use interp_fast (pair with sampling w/ detail)
if “high”: use interp_wedge (pair with sampling w/ detail)
idx (np.array, optional) – precomputed query/interpolation result
interpweights (np.array, optional) – precomputted interpolation weights
omega (bool) – if true, add value of ray solid angle instead of luminance
vm (raytraverse.mapper.ViewMapper, optional) –
rnd (bool, optional) – use random values as contribution (for visualizing data shape)
engine (raytraverse.renderer.Rtrace, optional) – engine for content aware interpolation
kwargs (dict, optional) – passed to interpolationn functions
- evaluate(skyvec, vm=None, idx=None, srconly=False, blursun=False, includeviews=True)[source]¶
return rays within view with skyvec applied. this is the analog to add_to_img for metric calculations
- Parameters:
skyvec (int float np.array, optional) – source coefficients, shape is (1,) or (srcn,)
vm (raytraverse.mapper.ViewMapper, optional) –
idx (np.array, optional) – precomputed query_ball result
srconly (bool, optional) – only evaluate direct sources (stored in self.srcviews)
includeviews (bool, optional) – include src views in returned results
- Returns:
rays (np.array) – shape (N, 3) rays falling within view
omega (np.array) – shape (N,) associated solid angles
lum (np.array) – shape (N,) associated luminances
- query_ray(vecs)[source]¶
return the index and distance of the nearest ray to each of vecs
- Parameters:
vecs (np.array) – shape (N, 3) normalized vectors to query, could represent image pixels for example.
- Returns:
i (np.array) – integer indices of closest ray to each query
d (np.array) – distance (corresponds to chord length on unit sphere) from query to ray in lightpoint. use translate.chord2theta to convert to angle.
- query_ball(vecs, viewangle=180)[source]¶
return set of rays within a view cone
- Parameters:
vecs (np.array) – shape (N, 3) vectors to query.
viewangle (int float) – opening angle of view cone
- Returns:
i – if vecs is a single point, a list of vector indices of rays within view cone. if vecs is a set of point an array of lists, one for each vec is returned.
- Return type:
list np.array
- direct_view(res=512, showsample=False, showweight=True, rnd=False, order=False, srcidx=None, interp=False, omega=False, scalefactor=1, vm=None, fisheye=True, grow=1, srcrnd=False)[source]¶
create an unweighted summary image of lightpoint
- add(lf2, src=None, calcomega=True, write=False, sumsrc=False)[source]¶
add light points of distinct sources together results in a new lightpoint with srcn=self.srcn+srcn2 and vector size=self.vecsize+vecsize2
- Parameters:
src (str, optional) – if None (default), src is “{lf1.src}_{lf2.src}”
calcomega (bool, optional) – passed to LightPointKD constructor
write (bool, optional) – passed to LightPointKD constructor
sumsrc (bool, optional) – if True adds matching source indices together (must be same shape) this assumes that the two lightpoints represent the same source but different components (such as direct/indirect)
- Returns:
will be subtyped according to self, unless lf2 is needed to preserve data
- Return type:
- update(vec, lum, omega=None, calcomega=True, write=True, filterviews=False)[source]¶
add additional rays to lightpoint in place
- Parameters:
vec (np.array, optional) – shape (N, >=3) where last three columns are normalized direction vectors of samples.
lum (np.array, optional) – reshapeable to (N, srcn). sample values for each source corresponding to vec.
omega (np.array, optional) – provide precomputed omega values, if given, overrides calcomega
calcomega (bool, optional) – if True (default) calculate solid angle of rays. This is unnecessary if point will be combined before calculating any metrics. setting to False will save some computation time. If False, resets omega to None!
write (bool, optional) – whether to save updated ray data to disk.
filterviews (bool, optional) – delete rays near sourceviews
SrcViewPoint¶
- class raytraverse.lightpoint.SrcViewPoint(scene, vecs, lum, pt=(0, 0, 0), posidx=0, src='sunview', res=64, srcomega=6.796702357283834e-05)[source]¶
Bases:
object
interface for sun view data
- scene¶
raytraverse.scene.Scene
- posidx¶
index for point
- Type:
int
- pt¶
point location
- Type:
np.array
- src¶
source key
- Type:
str
- radius¶
source radius
- Type:
float
- lum¶
source luminance (average)
- Type:
float
- property vm¶
CompressedPointKD¶
- class raytraverse.lightpoint.CompressedPointKD(scene, vec=None, lum=None, write=True, src=None, dist=0.0981, lerr=0.01, plotc=False, **kwargs)[source]¶
Bases:
LightPointKD
compressed data needs special methods for making images.
can be initialized either like LightPointKD (but with required omega argument), or if ‘scene’ is a LightPointKD then a compressed output is calculated from the input
- Parameters:
scene (BaseScene LightpointKD) –
src (str, optional) – new name for src passed to LightPointKD constructor
dist (float, optional) – translate.theta2chord(np.pi/32), primary clustering distance using the birch algorithm, for lossy compression of lf. this is the maximum radius of a cluster, preserving important directional information. clustering acts on ray direction and luminance, with weight of luminance dimension controlled by the lweight parameter.
lerr (float, optional) – min-max normalized error in luminance grouping.
plotc (bool, optional) – make directview plot of compressed output showing source vectors
- add_to_img(img, vecs, mask=None, skyvec=1, vm=None, fisheye=True, order=False, omega=False, rnd=False, **kwargs)[source]¶
add luminance contributions to image array (updates in place)
- Parameters:
img (np.array) – 2D image array to add to (either zeros or with other source)
vecs (np.array) – vectors corresponding to img pixels shape (N, 3)
mask (np.array, optional) – indices to img that correspond to vec (in case where whole image is not being updated, such as corners of fisheye)
skyvec (int float np.array, optional) – source coefficients, shape is (1,) or (srcn,)
vm (raytraverse.mapper.ViewMapper, optional) –
- compress(lp, src=None, dist=0.0981, lerr=0.01)[source]¶
A lossy compression based on clustering. Rays are clustered using the birch algoritm on a 4D vector (x,y,z,lum) where lum is the sum of contributions from all sources in the LightPoint. In the optional second stage (activated with secondary=True) sources are further grouped through agglomerative cluster using an average linkage. this is to help with source indentification/matching between LightPoints, but can introduce significant errors to computing non energy conserving metrics in cases where the applied sky vectors have large relative differences between adjacent patches (> 1.5:1) or if the variance in peak luminance above the lthreshold parameter is significant. These include cases where nearby transmitting materials is varied (example: a trans upper above a clear lower), or lthreshold is set too low. For this reason, it is better to use single stage compression for metric computation and only do glare source grouping for interpolation between LightPoints.
- Parameters:
lp (LightPointKD) –
src (str, optional) – new name for src passed to LightPointKD constructor
dist (float, optional) – translate.theta2chord(np.pi/32), primary clustering distance using the birch algorithm, for lossy compression of lf. this is the maximum radius of a cluster, preserving important directional information. clustering acts on ray direction and luminance, with weight of luminance dimension controlled by the lweight parameter.
lerr (float, optional) – min-max normalized error in luminance grouping.
plotc (bool, optional) – make directview plot of compressed output showing source vectors
- Return type:
arguments for initializing a CompressedPointKD