forked from FreyrS/dMaSIF
-
Notifications
You must be signed in to change notification settings - Fork 0
/
geometry_processing.py
822 lines (663 loc) · 32.8 KB
/
geometry_processing.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
import numpy as np
from math import pi
import torch
from pykeops.torch import LazyTensor
from plyfile import PlyData, PlyElement
from helper import *
import torch.nn as nn
import torch.nn.functional as F
# from matplotlib import pyplot as plt
from pykeops.torch.cluster import grid_cluster, cluster_ranges_centroids, from_matrix
from math import pi, sqrt
# Input-Output for tests =======================================================
import os
from pyvtk import PolyData, PointData, CellData, Scalars, Vectors, VtkData, PointData
def save_vtk(
fname, xyz, triangles=None, values=None, vectors=None, triangle_values=None
):
"""Saves a point cloud or triangle mesh as a .vtk file.
Files can be opened with Paraview or displayed using the PyVista library.
Args:
fname (string): filename.
xyz (Tensor): (N,3) point cloud or vertices.
triangles (integer Tensor, optional): (T,3) mesh connectivity. Defaults to None.
values (Tensor, optional): (N,D) values, supported by the vertices. Defaults to None.
vectors (Tensor, optional): (N,3) vectors, supported by the vertices. Defaults to None.
triangle_values (Tensor, optional): (T,D) values, supported by the triangles. Defaults to None.
"""
# Encode the points/vertices as a VTK structure:
if triangles is None: # Point cloud
structure = PolyData(points=numpy(xyz), vertices=np.arange(len(xyz)))
else: # Surface mesh
structure = PolyData(points=numpy(xyz), polygons=numpy(triangles))
data = [structure]
pointdata, celldata = [], []
# Point values - one channel per column of the `values` array:
if values is not None:
values = numpy(values)
if len(values.shape) == 1:
values = values[:, None]
features = values.T
pointdata += [
Scalars(f, name=f"features_{i:02d}") for i, f in enumerate(features)
]
# Point vectors - one vector per point:
if vectors is not None:
pointdata += [Vectors(numpy(vectors), name="vectors")]
# Store in the VTK object:
if pointdata != []:
pointdata = PointData(*pointdata)
data.append(pointdata)
# Triangle values - one channel per column of the `triangle_values` array:
if triangle_values is not None:
triangle_values = numpy(triangle_values)
if len(triangle_values.shape) == 1:
triangle_values = triangle_values[:, None]
features = triangle_values.T
celldata += [
Scalars(f, name=f"features_{i:02d}") for i, f in enumerate(features)
]
celldata = CellData(*celldata)
data.append(celldata)
# Write to hard drive:
vtk = VtkData(*data)
os.makedirs(os.path.dirname(fname), exist_ok=True)
vtk.tofile(fname)
# On-the-fly generation of the surfaces ========================================
def subsample(x, batch=None, scale=1.0):
"""Subsamples the point cloud using a grid (cubic) clustering scheme.
The function returns one average sample per cell, as described in Fig. 3.e)
of the paper.
Args:
x (Tensor): (N,3) point cloud.
batch (integer Tensor, optional): (N,) batch vector, as in PyTorch_geometric.
Defaults to None.
scale (float, optional): side length of the cubic grid cells. Defaults to 1 (Angstrom).
Returns:
(M,3): sub-sampled point cloud, with M <= N.
"""
if batch is None: # Single protein case:
if True: # Use a fast scatter_add_ implementation
labels = grid_cluster(x, scale).long()
C = labels.max() + 1
# We append a "1" to the input vectors, in order to
# compute both the numerator and denominator of the "average"
# fraction in one pass through the data.
x_1 = torch.cat((x, torch.ones_like(x[:, :1])), dim=1)
D = x_1.shape[1]
points = torch.zeros_like(x_1[:C])
points.scatter_add_(0, labels[:, None].repeat(1, D), x_1)
return (points[:, :-1] / points[:, -1:]).contiguous()
else: # Older implementation;
points = scatter(points * weights[:, None], labels, dim=0)
weights = scatter(weights, labels, dim=0)
points = points / weights[:, None]
else: # We process proteins using a for loop.
# This is probably sub-optimal, but I don't really know
# how to do more elegantly (this type of computation is
# not super well supported by PyTorch).
batch_size = torch.max(batch).item() + 1 # Typically, =32
points, batches = [], []
for b in range(batch_size):
p = subsample(x[batch == b], scale=scale)
points.append(p)
batches.append(b * torch.ones_like(batch[: len(p)]))
return torch.cat(points, dim=0), torch.cat(batches, dim=0)
def soft_distances(x, y, batch_x, batch_y, smoothness=0.01, atomtypes=None):
"""Computes a soft distance function to the atom centers of a protein.
Implements Eq. (1) of the paper in a fast and numerically stable way.
Args:
x (Tensor): (N,3) atom centers.
y (Tensor): (M,3) sampling locations.
batch_x (integer Tensor): (N,) batch vector for x, as in PyTorch_geometric.
batch_y (integer Tensor): (M,) batch vector for y, as in PyTorch_geometric.
smoothness (float, optional): atom radii if atom types are not provided. Defaults to .01.
atomtypes (integer Tensor, optional): (N,6) one-hot encoding of the atom chemical types. Defaults to None.
Returns:
Tensor: (M,) values of the soft distance function on the points `y`.
"""
# Build the (N, M, 1) symbolic matrix of squared distances:
x_i = LazyTensor(x[:, None, :]) # (N, 1, 3) atoms
y_j = LazyTensor(y[None, :, :]) # (1, M, 3) sampling points
D_ij = ((x_i - y_j) ** 2).sum(-1) # (N, M, 1) squared distances
# Use a block-diagonal sparsity mask to support heterogeneous batch processing:
D_ij.ranges = diagonal_ranges(batch_x, batch_y)
if atomtypes is not None:
# Turn the one-hot encoding "atomtypes" into a vector of diameters "smoothness_i":
# (N, 6) -> (N, 1, 1) (There are 6 atom types)
atomic_radii = torch.cuda.FloatTensor(
[170, 110, 152, 155, 180, 190], device=x.device
)
atomic_radii = atomic_radii / atomic_radii.min()
atomtype_radii = atomtypes * atomic_radii[None, :] # n_atoms, n_atomtypes
# smoothness = atomtypes @ atomic_radii # (N, 6) @ (6,) = (N,)
smoothness = torch.sum(
smoothness * atomtype_radii, dim=1, keepdim=False
) # n_atoms, 1
smoothness_i = LazyTensor(smoothness[:, None, None])
# Compute an estimation of the mean smoothness in a neighborhood
# of each sampling point:
# density = (-D_ij.sqrt()).exp().sum(0).view(-1) # (M,) local density of atoms
# smooth = (smoothness_i * (-D_ij.sqrt()).exp()).sum(0).view(-1) # (M,)
# mean_smoothness = smooth / density # (M,)
# soft_dists = -mean_smoothness * (
# (-D_ij.sqrt() / smoothness_i).logsumexp(dim=0)
# ).view(-1)
mean_smoothness = (-D_ij.sqrt()).exp().sum(0)
mean_smoothness_j = LazyTensor(mean_smoothness[None, :, :])
mean_smoothness = (
smoothness_i * (-D_ij.sqrt()).exp() / mean_smoothness_j
) # n_atoms, n_points, 1
mean_smoothness = mean_smoothness.sum(0).view(-1)
soft_dists = -mean_smoothness * (
(-D_ij.sqrt() / smoothness_i).logsumexp(dim=0)
).view(-1)
else:
soft_dists = -smoothness * ((-D_ij.sqrt() / smoothness).logsumexp(dim=0)).view(
-1
)
return soft_dists
def atoms_to_points_normals(
atoms,
batch,
distance=1.05,
smoothness=0.5,
resolution=1.0,
nits=4,
atomtypes=None,
sup_sampling=20,
variance=0.1,
):
"""Turns a collection of atoms into an oriented point cloud.
Sampling algorithm for protein surfaces, described in Fig. 3 of the paper.
Args:
atoms (Tensor): (N,3) coordinates of the atom centers `a_k`.
batch (integer Tensor): (N,) batch vector, as in PyTorch_geometric.
distance (float, optional): value of the level set to sample from
the smooth distance function. Defaults to 1.05.
smoothness (float, optional): radii of the atoms, if atom types are
not provided. Defaults to 0.5.
resolution (float, optional): side length of the cubic cells in
the final sub-sampling pass. Defaults to 1.0.
nits (int, optional): number of iterations . Defaults to 4.
atomtypes (Tensor, optional): (N,6) one-hot encoding of the atom
chemical types. Defaults to None.
Returns:
(Tensor): (M,3) coordinates for the surface points `x_i`.
(Tensor): (M,3) unit normals `n_i`.
(integer Tensor): (M,) batch vector, as in PyTorch_geometric.
"""
# a) Parameters for the soft distance function and its level set:
T = distance
N, D = atoms.shape
B = sup_sampling # Sup-sampling ratio
# Batch vectors:
batch_atoms = batch
batch_z = batch[:, None].repeat(1, B).view(N * B)
# b) Draw N*B points at random in the neighborhood of our atoms
z = atoms[:, None, :] + 10 * T * torch.randn(N, B, D).type_as(atoms)
z = z.view(-1, D) # (N*B, D)
# We don't want to backprop through a full network here!
atoms = atoms.detach().contiguous()
z = z.detach().contiguous()
# N.B.: Test mode disables the autograd engine: we must switch it on explicitely.
with torch.enable_grad():
if z.is_leaf:
z.requires_grad = True
# c) Iterative loop: gradient descent along the potential
# ".5 * (dist - T)^2" with respect to the positions z of our samples
for it in range(nits):
dists = soft_distances(
atoms,
z,
batch_atoms,
batch_z,
smoothness=smoothness,
atomtypes=atomtypes,
)
Loss = ((dists - T) ** 2).sum()
g = torch.autograd.grad(Loss, z)[0]
z.data -= 0.5 * g
# d) Only keep the points which are reasonably close to the level set:
dists = soft_distances(
atoms, z, batch_atoms, batch_z, smoothness=smoothness, atomtypes=atomtypes
)
margin = (dists - T).abs()
mask = margin < variance * T
# d') And remove the points that are trapped *inside* the protein:
zz = z.detach()
zz.requires_grad = True
for it in range(nits):
dists = soft_distances(
atoms,
zz,
batch_atoms,
batch_z,
smoothness=smoothness,
atomtypes=atomtypes,
)
Loss = (1.0 * dists).sum()
g = torch.autograd.grad(Loss, zz)[0]
normals = F.normalize(g, p=2, dim=-1) # (N, 3)
zz = zz + 1.0 * T * normals
dists = soft_distances(
atoms, zz, batch_atoms, batch_z, smoothness=smoothness, atomtypes=atomtypes
)
mask = mask & (dists > 1.5 * T)
z = z[mask].contiguous().detach()
batch_z = batch_z[mask].contiguous().detach()
# e) Subsample the point cloud:
points, batch_points = subsample(z, batch_z, scale=resolution)
# f) Compute the normals on this smaller point cloud:
p = points.detach()
p.requires_grad = True
dists = soft_distances(
atoms,
p,
batch_atoms,
batch_points,
smoothness=smoothness,
atomtypes=atomtypes,
)
Loss = (1.0 * dists).sum()
g = torch.autograd.grad(Loss, p)[0]
normals = F.normalize(g, p=2, dim=-1) # (N, 3)
points = points - 0.5 * normals
return points.detach(), normals.detach(), batch_points.detach()
# Surface mesh -> Normals ======================================================
def mesh_normals_areas(vertices, triangles=None, scale=[1.0], batch=None, normals=None):
"""Returns a smooth field of normals, possibly at different scales.
points, triangles or normals, scale(s) -> normals
(N, 3), (3, T) or (N,3), (S,) -> (N, 3) or (N, S, 3)
Simply put - if `triangles` are provided:
1. Normals are first computed for every triangle using simple 3D geometry
and are weighted according to surface area.
2. The normal at any given vertex is then computed as the weighted average
of the normals of all triangles in a neighborhood specified
by Gaussian windows whose radii are given in the list of "scales".
If `normals` are provided instead, we simply smooth the discrete vector
field using Gaussian windows whose radii are given in the list of "scales".
If more than one scale is provided, normal fields are computed in parallel
and returned in a single 3D tensor.
Args:
vertices (Tensor): (N,3) coordinates of mesh vertices or 3D points.
triangles (integer Tensor, optional): (3,T) mesh connectivity. Defaults to None.
scale (list of floats, optional): (S,) radii of the Gaussian smoothing windows. Defaults to [1.].
batch (integer Tensor, optional): batch vector, as in PyTorch_geometric. Defaults to None.
normals (Tensor, optional): (N,3) raw normals vectors on the vertices. Defaults to None.
Returns:
(Tensor): (N,3) or (N,S,3) point normals.
(Tensor): (N,) point areas, if triangles were provided.
"""
# Single- or Multi-scale mode:
if hasattr(scale, "__len__"):
scales, single_scale = scale, False
else:
scales, single_scale = [scale], True
scales = torch.Tensor(scales).type_as(vertices) # (S,)
# Compute the "raw" field of normals:
if triangles is not None:
# Vertices of all triangles in the mesh:
A = vertices[triangles[0, :]] # (N, 3)
B = vertices[triangles[1, :]] # (N, 3)
C = vertices[triangles[2, :]] # (N, 3)
# Triangle centers and normals (length = surface area):
centers = (A + B + C) / 3 # (N, 3)
V = (B - A).cross(C - A) # (N, 3)
# Vertice areas:
S = (V ** 2).sum(-1).sqrt() / 6 # (N,) 1/3 of a triangle area
areas = torch.zeros(len(vertices)).type_as(vertices) # (N,)
areas.scatter_add_(0, triangles[0, :], S) # Aggregate from "A's"
areas.scatter_add_(0, triangles[1, :], S) # Aggregate from "B's"
areas.scatter_add_(0, triangles[2, :], S) # Aggregate from "C's"
else: # Use "normals" instead
areas = None
V = normals
centers = vertices
# Normal of a vertex = average of all normals in a ball of size "scale":
x_i = LazyTensor(vertices[:, None, :]) # (N, 1, 3)
y_j = LazyTensor(centers[None, :, :]) # (1, M, 3)
v_j = LazyTensor(V[None, :, :]) # (1, M, 3)
s = LazyTensor(scales[None, None, :]) # (1, 1, S)
D_ij = ((x_i - y_j) ** 2).sum(-1) # (N, M, 1)
K_ij = (-D_ij / (2 * s ** 2)).exp() # (N, M, S)
# Support for heterogeneous batch processing:
if batch is not None:
batch_vertices = batch
batch_centers = batch[triangles[0, :]] if triangles is not None else batch
K_ij.ranges = diagonal_ranges(batch_vertices, batch_centers)
if single_scale:
U = (K_ij * v_j).sum(dim=1) # (N, 3)
else:
U = (K_ij.tensorprod(v_j)).sum(dim=1) # (N, S*3)
U = U.view(-1, len(scales), 3) # (N, S, 3)
normals = F.normalize(U, p=2, dim=-1) # (N, 3) or (N, S, 3)
return normals, areas
# Compute tangent planes and curvatures ========================================
def tangent_vectors(normals):
"""Returns a pair of vector fields u and v to complete the orthonormal basis [n,u,v].
normals -> uv
(N, 3) or (N, S, 3) -> (N, 2, 3) or (N, S, 2, 3)
This routine assumes that the 3D "normal" vectors are normalized.
It is based on the 2017 paper from Pixar, "Building an orthonormal basis, revisited".
Args:
normals (Tensor): (N,3) or (N,S,3) normals `n_i`, i.e. unit-norm 3D vectors.
Returns:
(Tensor): (N,2,3) or (N,S,2,3) unit vectors `u_i` and `v_i` to complete
the tangent coordinate systems `[n_i,u_i,v_i].
"""
x, y, z = normals[..., 0], normals[..., 1], normals[..., 2]
s = (2 * (z >= 0)) - 1.0 # = z.sign(), but =1. if z=0.
a = -1 / (s + z)
b = x * y * a
uv = torch.stack((1 + s * x * x * a, s * b, -s * x, b, s + y * y * a, -y), dim=-1)
uv = uv.view(uv.shape[:-1] + (2, 3))
return uv
def curvatures(
vertices, triangles=None, scales=[1.0], batch=None, normals=None, reg=0.01
):
"""Returns a collection of mean (H) and Gauss (K) curvatures at different scales.
points, faces, scales -> (H_1, K_1, ..., H_S, K_S)
(N, 3), (3, N), (S,) -> (N, S*2)
We rely on a very simple linear regression method, for all vertices:
1. Estimate normals and surface areas.
2. Compute a local tangent frame.
3. In a pseudo-geodesic Gaussian neighborhood at scale s,
compute the two (2, 2) covariance matrices PPt and PQt
between the displacement vectors "P = x_i - x_j" and
the normals "Q = n_i - n_j", projected on the local tangent plane.
4. Up to the sign, the shape operator S at scale s is then approximated
as "S = (reg**2 * I_2 + PPt)^-1 @ PQt".
5. The mean and Gauss curvatures are the trace and determinant of
this (2, 2) matrix.
As of today, this implementation does not weigh points by surface areas:
this could make a sizeable difference if protein surfaces were not
sub-sampled to ensure uniform sampling density.
For convergence analysis, see for instance
"Efficient curvature estimation for oriented point clouds",
Cao, Li, Sun, Assadi, Zhang, 2019.
Args:
vertices (Tensor): (N,3) coordinates of the points or mesh vertices.
triangles (integer Tensor, optional): (3,T) mesh connectivity. Defaults to None.
scales (list of floats, optional): list of (S,) smoothing scales. Defaults to [1.].
batch (integer Tensor, optional): batch vector, as in PyTorch_geometric. Defaults to None.
normals (Tensor, optional): (N,3) field of "raw" unit normals. Defaults to None.
reg (float, optional): small amount of Tikhonov/ridge regularization
in the estimation of the shape operator. Defaults to .01.
Returns:
(Tensor): (N, S*2) tensor of mean and Gauss curvatures computed for
every point at the required scales.
"""
# Number of points, number of scales:
N, S = vertices.shape[0], len(scales)
ranges = diagonal_ranges(batch)
# Compute the normals at different scales + vertice areas:
normals_s, _ = mesh_normals_areas(
vertices, triangles=triangles, normals=normals, scale=scales, batch=batch
) # (N, S, 3), (N,)
# Local tangent bases:
uv_s = tangent_vectors(normals_s) # (N, S, 2, 3)
features = []
for s, scale in enumerate(scales):
# Extract the relevant descriptors at the current scale:
normals = normals_s[:, s, :].contiguous() # (N, 3)
uv = uv_s[:, s, :, :].contiguous() # (N, 2, 3)
# Encode as symbolic tensors:
# Points:
x_i = LazyTensor(vertices.view(N, 1, 3))
x_j = LazyTensor(vertices.view(1, N, 3))
# Normals:
n_i = LazyTensor(normals.view(N, 1, 3))
n_j = LazyTensor(normals.view(1, N, 3))
# Tangent bases:
uv_i = LazyTensor(uv.view(N, 1, 6))
# Pseudo-geodesic squared distance:
d2_ij = ((x_j - x_i) ** 2).sum(-1) * ((2 - (n_i | n_j)) ** 2) # (N, N, 1)
# Gaussian window:
window_ij = (-d2_ij / (2 * (scale ** 2))).exp() # (N, N, 1)
# Project on the tangent plane:
P_ij = uv_i.matvecmult(x_j - x_i) # (N, N, 2)
Q_ij = uv_i.matvecmult(n_j - n_i) # (N, N, 2)
# Concatenate:
PQ_ij = P_ij.concat(Q_ij) # (N, N, 2+2)
# Covariances, with a scale-dependent weight:
PPt_PQt_ij = P_ij.tensorprod(PQ_ij) # (N, N, 2*(2+2))
PPt_PQt_ij = window_ij * PPt_PQt_ij # (N, N, 2*(2+2))
# Reduction - with batch support:
PPt_PQt_ij.ranges = ranges
PPt_PQt = PPt_PQt_ij.sum(1) # (N, 2*(2+2))
# Reshape to get the two covariance matrices:
PPt_PQt = PPt_PQt.view(N, 2, 2, 2)
PPt, PQt = PPt_PQt[:, :, 0, :], PPt_PQt[:, :, 1, :] # (N, 2, 2), (N, 2, 2)
# Add a small ridge regression:
PPt[:, 0, 0] += reg
PPt[:, 1, 1] += reg
# (minus) Shape operator, i.e. the differential of the Gauss map:
# = (PPt^-1 @ PQt) : simple estimation through linear regression
S = torch.solve(PQt, PPt).solution
a, b, c, d = S[:, 0, 0], S[:, 0, 1], S[:, 1, 0], S[:, 1, 1] # (N,)
# Normalization
mean_curvature = a + d
gauss_curvature = a * d - b * c
features += [mean_curvature.clamp(-1, 1), gauss_curvature.clamp(-1, 1)]
features = torch.stack(features, dim=-1)
return features
# Fast tangent convolution layer ===============================================
class ContiguousBackward(torch.autograd.Function):
"""
Function to ensure contiguous gradient in backward pass. To be applied after PyKeOps reduction.
N.B.: This workaround fixes a bug that will be fixed in ulterior KeOp releases.
"""
@staticmethod
def forward(ctx, input):
return input
@staticmethod
def backward(ctx, grad_output):
return grad_output.contiguous()
class dMaSIFConv(nn.Module):
def __init__(
self, in_channels=1, out_channels=1, radius=1.0, hidden_units=None, cheap=False
):
"""Creates the KeOps convolution layer.
I = in_channels is the dimension of the input features
O = out_channels is the dimension of the output features
H = hidden_units is the dimension of the intermediate representation
radius is the size of the pseudo-geodesic Gaussian window w_ij = W(d_ij)
This affordable layer implements an elementary "convolution" operator
on a cloud of N points (x_i) in dimension 3 that we decompose in three steps:
1. Apply the MLP "net_in" on the input features "f_i". (N, I) -> (N, H)
2. Compute H interaction terms in parallel with:
f_i = sum_j [ w_ij * conv(P_ij) * f_j ]
In the equation above:
- w_ij is a pseudo-geodesic window with a set radius.
- P_ij is a vector of dimension 3, equal to "x_j-x_i"
in the local oriented basis at x_i.
- "conv" is an MLP from R^3 to R^H:
- with 1 linear layer if "cheap" is True;
- with 2 linear layers and C=8 intermediate "cuts" otherwise.
- "*" is coordinate-wise product.
- f_j is the vector of transformed features.
3. Apply the MLP "net_out" on the output features. (N, H) -> (N, O)
A more general layer would have implemented conv(P_ij) as a full
(H, H) matrix instead of a mere (H,) vector... At a much higher
computational cost. The reasoning behind the code below is that
a given time budget is better spent on using a larger architecture
and more channels than on a very complex convolution operator.
Interactions between channels happen at steps 1. and 3.,
whereas the (costly) point-to-point interaction step 2.
lets the network aggregate information in spatial neighborhoods.
Args:
in_channels (int, optional): numper of input features per point. Defaults to 1.
out_channels (int, optional): number of output features per point. Defaults to 1.
radius (float, optional): deviation of the Gaussian window on the
quasi-geodesic distance `d_ij`. Defaults to 1..
hidden_units (int, optional): number of hidden features per point.
Defaults to out_channels.
cheap (bool, optional): shall we use a 1-layer deep Filter,
instead of a 2-layer deep MLP? Defaults to False.
"""
super(dMaSIFConv, self).__init__()
self.Input = in_channels
self.Output = out_channels
self.Radius = radius
self.Hidden = self.Output if hidden_units is None else hidden_units
self.Cuts = 8 # Number of hidden units for the 3D MLP Filter.
self.cheap = cheap
# For performance reasons, we cut our "hidden" vectors
# in n_heads "independent heads" of dimension 8.
self.heads_dim = 8 # 4 is probably too small; 16 is certainly too big
# We accept "Hidden" dimensions of size 1, 2, 3, 4, 5, 6, 7, 8, 16, 32, 64, ...
if self.Hidden < self.heads_dim:
self.heads_dim = self.Hidden
if self.Hidden % self.heads_dim != 0:
raise ValueError(f"The dimension of the hidden units ({self.Hidden})"\
+ f"should be a multiple of the heads dimension ({self.heads_dim}).")
else:
self.n_heads = self.Hidden // self.heads_dim
# Transformation of the input features:
self.net_in = nn.Sequential(
nn.Linear(self.Input, self.Hidden), # (H, I) + (H,)
nn.LeakyReLU(negative_slope=0.2),
nn.Linear(self.Hidden, self.Hidden), # (H, H) + (H,)
# nn.LayerNorm(self.Hidden),#nn.BatchNorm1d(self.Hidden),
nn.LeakyReLU(negative_slope=0.2),
) # (H,)
self.norm_in = nn.GroupNorm(4, self.Hidden)
# self.norm_in = nn.LayerNorm(self.Hidden)
# self.norm_in = nn.Identity()
# 3D convolution filters, encoded as an MLP:
if cheap:
self.conv = nn.Sequential(
nn.Linear(3, self.Hidden), nn.ReLU() # (H, 3) + (H,)
) # KeOps does not support well LeakyReLu
else:
self.conv = nn.Sequential(
nn.Linear(3, self.Cuts), # (C, 3) + (C,)
nn.ReLU(), # KeOps does not support well LeakyReLu
nn.Linear(self.Cuts, self.Hidden),
) # (H, C) + (H,)
# Transformation of the output features:
self.net_out = nn.Sequential(
nn.Linear(self.Hidden, self.Output), # (O, H) + (O,)
nn.LeakyReLU(negative_slope=0.2),
nn.Linear(self.Output, self.Output), # (O, O) + (O,)
# nn.LayerNorm(self.Output),#nn.BatchNorm1d(self.Output),
nn.LeakyReLU(negative_slope=0.2),
) # (O,)
self.norm_out = nn.GroupNorm(4, self.Output)
# self.norm_out = nn.LayerNorm(self.Output)
# self.norm_out = nn.Identity()
# Custom initialization for the MLP convolution filters:
# we get interesting piecewise affine cuts on a normalized neighborhood.
with torch.no_grad():
nn.init.normal_(self.conv[0].weight)
nn.init.uniform_(self.conv[0].bias)
self.conv[0].bias *= 0.8 * (self.conv[0].weight ** 2).sum(-1).sqrt()
if not cheap:
nn.init.uniform_(
self.conv[2].weight,
a=-1 / np.sqrt(self.Cuts),
b=1 / np.sqrt(self.Cuts),
)
nn.init.normal_(self.conv[2].bias)
self.conv[2].bias *= 0.5 * (self.conv[2].weight ** 2).sum(-1).sqrt()
def forward(self, points, nuv, features, ranges=None):
"""Performs a quasi-geodesic interaction step.
points, local basis, in features -> out features
(N, 3), (N, 3, 3), (N, I) -> (N, O)
This layer computes the interaction step of Eq. (7) in the paper,
in-between the application of two MLP networks independently on all
feature vectors.
Args:
points (Tensor): (N,3) point coordinates `x_i`.
nuv (Tensor): (N,3,3) local coordinate systems `[n_i,u_i,v_i]`.
features (Tensor): (N,I) input feature vectors `f_i`.
ranges (6-uple of integer Tensors, optional): low-level format
to support batch processing, as described in the KeOps documentation.
In practice, this will be built by a higher-level object
to encode the relevant "batch vectors" in a way that is convenient
for the KeOps CUDA engine. Defaults to None.
Returns:
(Tensor): (N,O) output feature vectors `f'_i`.
"""
# 1. Transform the input features: -------------------------------------
features = self.net_in(features) # (N, I) -> (N, H)
features = features.transpose(1, 0)[None, :, :] # (1,H,N)
features = self.norm_in(features)
features = features[0].transpose(1, 0).contiguous() # (1, H, N) -> (N, H)
# 2. Compute the local "shape contexts": -------------------------------
# 2.a Normalize the kernel radius:
points = points / (sqrt(2.0) * self.Radius) # (N, 3)
# 2.b Encode the variables as KeOps LazyTensors
# Vertices:
x_i = LazyTensor(points[:, None, :]) # (N, 1, 3)
x_j = LazyTensor(points[None, :, :]) # (1, N, 3)
# WARNING - Here, we assume that the normals are fixed:
normals = (
nuv[:, 0, :].contiguous().detach()
) # (N, 3) - remove the .detach() if needed
# Local bases:
nuv_i = LazyTensor(nuv.view(-1, 1, 9)) # (N, 1, 9)
# Normals:
n_i = nuv_i[:3] # (N, 1, 3)
n_j = LazyTensor(normals[None, :, :]) # (1, N, 3)
# To avoid register spilling when using large embeddings, we perform our KeOps reduction
# over the vector of length "self.Hidden = self.n_heads * self.heads_dim"
# as self.n_heads reduction over vectors of length self.heads_dim (= "Hd" in the comments).
head_out_features = []
for head in range(self.n_heads):
# Extract a slice of width Hd from the feature array
head_start = head * self.heads_dim
head_end = head_start + self.heads_dim
head_features = features[:, head_start:head_end].contiguous() # (N, H) -> (N, Hd)
# Features:
f_j = LazyTensor(head_features[None, :, :]) # (1, N, Hd)
# Convolution parameters:
if self.cheap:
# Extract a slice of Hd lines: (H, 3) -> (Hd, 3)
A = self.conv[0].weight[head_start:head_end, :].contiguous()
# Extract a slice of Hd coefficients: (H,) -> (Hd,)
B = self.conv[0].bias[head_start:head_end].contiguous()
AB = torch.cat((A, B[:, None]), dim=1) # (Hd, 4)
ab = LazyTensor(AB.view(1, 1, -1)) # (1, 1, Hd*4)
else:
A_1, B_1 = self.conv[0].weight, self.conv[0].bias # (C, 3), (C,)
# Extract a slice of Hd lines: (H, C) -> (Hd, C)
A_2 = self.conv[2].weight[head_start:head_end, :].contiguous()
# Extract a slice of Hd coefficients: (H,) -> (Hd,)
B_2 = self.conv[2].bias[head_start:head_end].contiguous()
a_1 = LazyTensor(A_1.view(1, 1, -1)) # (1, 1, C*3)
b_1 = LazyTensor(B_1.view(1, 1, -1)) # (1, 1, C)
a_2 = LazyTensor(A_2.view(1, 1, -1)) # (1, 1, Hd*C)
b_2 = LazyTensor(B_2.view(1, 1, -1)) # (1, 1, Hd)
# 2.c Pseudo-geodesic window:
# Pseudo-geodesic squared distance:
d2_ij = ((x_j - x_i) ** 2).sum(-1) * ((2 - (n_i | n_j)) ** 2) # (N, N, 1)
# Gaussian window:
window_ij = (-d2_ij).exp() # (N, N, 1)
# 2.d Local MLP:
# Local coordinates:
X_ij = nuv_i.matvecmult(x_j - x_i) # (N, N, 9) "@" (N, N, 3) = (N, N, 3)
# MLP:
if self.cheap:
X_ij = ab.matvecmult(
X_ij.concat(LazyTensor(1))
) # (N, N, Hd*4) @ (N, N, 3+1) = (N, N, Hd)
X_ij = X_ij.relu() # (N, N, Hd)
else:
X_ij = a_1.matvecmult(X_ij) + b_1 # (N, N, C)
X_ij = X_ij.relu() # (N, N, C)
X_ij = a_2.matvecmult(X_ij) + b_2 # (N, N, Hd)
X_ij = X_ij.relu()
# 2.e Actual computation:
F_ij = window_ij * X_ij * f_j # (N, N, Hd)
F_ij.ranges = ranges # Support for batches and/or block-sparsity
head_out_features.append(ContiguousBackward().apply(F_ij.sum(dim=1))) # (N, Hd)
# Concatenate the result of our n_heads "attention heads":
features = torch.cat(head_out_features, dim=1) # n_heads * (N, Hd) -> (N, H)
# 3. Transform the output features: ------------------------------------
features = self.net_out(features) # (N, H) -> (N, O)
features = features.transpose(1, 0)[None, :, :] # (1,O,N)
features = self.norm_out(features)
features = features[0].transpose(1, 0).contiguous()
return features