-
Notifications
You must be signed in to change notification settings - Fork 18
/
g2.py
349 lines (293 loc) · 14.4 KB
/
g2.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
import numpy as np
from itertools import chain, product
from splipy import Curve, Surface, Volume, SplineObject, BSplineBasis, TrimmedSurface
from splipy.utils import flip_and_move_plane_geometry, rotate_local_x_axis
from .master import MasterIO
import splipy.surface_factory as SurfaceFactory
import splipy.curve_factory as CurveFactory
import splipy.state as state
from numpy import sqrt, pi
class G2(MasterIO):
def read_next_non_whitespace(self):
line = next(self.fstream).strip()
while not line:
line = next(self.fstream).strip()
return line
def circle(self):
dim = int( self.read_next_non_whitespace().strip())
r = float( next(self.fstream).strip())
center= np.array(next(self.fstream).split(), dtype=float)
normal= np.array(next(self.fstream).split(), dtype=float)
xaxis = np.array(next(self.fstream).split(), dtype=float)
param = np.array(next(self.fstream).split(), dtype=float)
reverse = next(self.fstream).strip() != '0'
result = CurveFactory.circle(r=r, center=center, normal=normal, xaxis=xaxis)
result.reparam(param)
if reverse:
result.reverse()
return result
def ellipse(self):
dim = int( self.read_next_non_whitespace().strip())
r1 = float( next(self.fstream).strip())
r2 = float( next(self.fstream).strip())
center= np.array(next(self.fstream).split(), dtype=float)
normal= np.array(next(self.fstream).split(), dtype=float)
xaxis = np.array(next(self.fstream).split(), dtype=float)
param = np.array(next(self.fstream).split(), dtype=float)
reverse = next(self.fstream).strip() != '0'
result = CurveFactory.ellipse(r1=r1, r2=r2, center=center, normal=normal, xaxis=xaxis)
result.reparam(param)
if reverse:
result.reverse()
return result
def line(self):
dim = int( self.read_next_non_whitespace().strip())
start = np.array(next(self.fstream).split(), dtype=float)
direction= np.array(next(self.fstream).split(), dtype=float)
finite = next(self.fstream).strip() != '0'
param = np.array(next(self.fstream).split(), dtype=float)
reverse = next(self.fstream).strip() != '0'
d = np.array(direction)
s = np.array(start)
# d /= np.linalg.norm(d)
if not finite:
param = [-state.unlimited, +state.unlimited]
result = CurveFactory.line(s+d*param[0], s+d*param[1])
if reverse:
result.reverse()
return result
# def cone(self):
# dim = int( self.read_next_non_whitespace().strip())
# r = float( next(self.fstream).strip())
# center = np.array(next(self.fstream).split(' '), dtype=float)
# z_axis = np.array(next(self.fstream).split(' '), dtype=float)
# x_axis = np.array(next(self.fstream).split(' '), dtype=float)
# angle = float( next(self.fstream).strip())
# finite = next(self.fstream).strip() != '0'
# param_u = np.array(next(self.fstream).split(' '), dtype=float)
# if finite:
# param_v=np.array(next(self.fstream).split(' '), dtype=float)
def cylinder(self):
dim = int( self.read_next_non_whitespace().strip())
r = float( next(self.fstream).strip())
center = np.array(next(self.fstream).split(), dtype=float)
z_axis = np.array(next(self.fstream).split(), dtype=float)
x_axis = np.array(next(self.fstream).split(), dtype=float)
finite = next(self.fstream).strip() != '0'
param_u = np.array(next(self.fstream).split(), dtype=float)
if finite:
param_v=np.array(next(self.fstream).split(), dtype=float)
else:
param_v=[-state.unlimited, state.unlimited]
swap = next(self.fstream).strip() != '0'
center = center + z_axis*param_v[0]
h = param_v[1] - param_v[0]
result = SurfaceFactory.cylinder(r=r, center=center, xaxis=x_axis, axis=z_axis, h=h)
result.reparam(param_u, param_v)
if swap:
result.swap()
return result
def disc(self):
dim = int( self.read_next_non_whitespace().strip())
center = np.array(next(self.fstream).split(), dtype=float)
r = float( next(self.fstream).strip())
z_axis = np.array(next(self.fstream).split(), dtype=float)
x_axis = np.array(next(self.fstream).split(), dtype=float)
degen = next(self.fstream).strip() != '0'
angles = [float( next(self.fstream).strip()) for i in range(4)]
param_u = np.array(next(self.fstream).split(), dtype=float)
param_v = np.array(next(self.fstream).split(), dtype=float)
swap = next(self.fstream).strip() != '0'
if degen:
result = SurfaceFactory.disc(r=r, center=center, xaxis=x_axis, normal=z_axis, type='radial')
else:
if not(np.allclose(np.diff(angles), pi/2, atol=1e-10)):
raise RuntimeError('Unknown square parametrization of disc elementary surface')
result = SurfaceFactory.disc(r=r, center=center, xaxis=x_axis, normal=z_axis, type='square')
result.reparam(param_u, param_v)
if swap:
result.swap()
return result
def plane(self):
dim = int( self.read_next_non_whitespace().strip())
center = np.array(next(self.fstream).split(), dtype=float)
normal = np.array(next(self.fstream).split(), dtype=float)
x_axis = np.array(next(self.fstream).split(), dtype=float)
finite = next(self.fstream).strip() != '0'
if finite:
param_u= np.array(next(self.fstream).split(), dtype=float)
param_v= np.array(next(self.fstream).split(), dtype=float)
else:
param_u= [-state.unlimited, +state.unlimited]
param_v= [-state.unlimited, +state.unlimited]
swap = next(self.fstream).strip() != '0'
result = Surface() * [param_u[1]-param_u[0], param_v[1]-param_v[0]] + [param_u[0],param_v[0]]
result.rotate(rotate_local_x_axis(x_axis, normal))
result = flip_and_move_plane_geometry(result,center,normal)
result.reparam(param_u, param_v)
if(swap):
result.swap()
return result
def torus(self):
dim = int( self.read_next_non_whitespace().strip())
r2 = float( next(self.fstream).strip())
r1 = float( next(self.fstream).strip())
center = np.array(next(self.fstream).split(), dtype=float)
z_axis = np.array(next(self.fstream).split(), dtype=float)
x_axis = np.array(next(self.fstream).split(), dtype=float)
select_out= next(self.fstream).strip() != '0' # I have no idea what this does :(
param_u = np.array(next(self.fstream).split(), dtype=float)
param_v = np.array(next(self.fstream).split(), dtype=float)
swap = next(self.fstream).strip() != '0'
result = SurfaceFactory.torus(minor_r=r1, major_r=r2, center=center, normal=z_axis, xaxis=x_axis)
result.reparam(param_u, param_v)
if(swap):
result.swap()
return result
def sphere(self):
dim = int( self.read_next_non_whitespace().strip())
r = float( next(self.fstream).strip())
center = np.array(next(self.fstream).split(), dtype=float)
z_axis = np.array(next(self.fstream).split(), dtype=float)
x_axis = np.array(next(self.fstream).split(), dtype=float)
param_u = np.array(next(self.fstream).split(), dtype=float)
param_v = np.array(next(self.fstream).split(), dtype=float)
swap = next(self.fstream).strip() != '0'
result = SurfaceFactory.sphere(r=r, center=center, xaxis=x_axis, zaxis=z_axis).swap()
if swap:
result.swap()
result.reparam(param_u, param_v)
return result
def splines(self, pardim):
cls = G2.classes[pardim-1]
_, rational = self.read_next_non_whitespace().strip().split()
rational = bool(int(rational))
bases = [self.read_basis() for _ in range(pardim)]
ncps = 1
for b in bases:
ncps *= b.num_functions()
cps = [tuple(map(float, next(self.fstream).split()))
for _ in range(ncps)]
args = bases + [cps, rational]
return cls(*args)
def surface_of_linear_extrusion(self):
dim = int( self.read_next_non_whitespace().strip())
crv = self.splines(1)
normal = np.array(self.read_next_non_whitespace().split(), dtype=float)
finite = next(self.fstream).strip() != '0'
param_u = np.array(next(self.fstream).split(), dtype=float)
if finite:
param_v=np.array(next(self.fstream).split(), dtype=float)
else:
param_v=[-state.unlimited, +state.unlimited]
swap = next(self.fstream).strip() != '0'
result = SurfaceFactory.extrude(crv + normal*param_v[0], normal*(param_v[1]-param_v[0]))
result.reparam(param_u, param_v)
if swap:
result.swap()
return result
def bounded_surface(self):
objtype = int( next(self.fstream).strip() )
# create the underlying surface which all trimming curves are to be applied
if objtype in G2.g2_generators:
constructor = getattr(self, G2.g2_generators[objtype].__name__)
surface = constructor()
elif objtype == 200:
surface = self.splines(2)
else:
raise IOError('Unsopported trimmed surface or malformed input file')
# for all trimming loops
numb_loops = int( self.read_next_non_whitespace() )
all_loops = []
for i in range(numb_loops):
# for all cuve pieces of that loop
numb_crvs, space_epsilon = next(self.fstream).split()
state.parametric_absolute_tolerance = float(space_epsilon)
one_loop = []
for j in range(int(numb_crvs)):
# read a physical and parametric representation of the same curve
_, parameter_curve_type, space_curve_type = map(int, self.read_next_non_whitespace().split())
two_curves = []
for crv_type in [parameter_curve_type, space_curve_type]:
if crv_type in G2.g2_generators:
constructor = getattr(self, G2.g2_generators[crv_type].__name__)
crv = constructor()
elif crv_type == 100:
crv = self.splines(1)
else:
raise IOError('Unsopported trimming curve or malformed input file')
two_curves.append(crv)
# only keep the parametric version (re-generate physical one if we need it)
one_loop.append(two_curves[0])
self.trimming_curves.append(two_curves[1])
all_loops.append(one_loop)
return TrimmedSurface(surface.bases[0], surface.bases[1], surface.controlpoints, surface.rational, all_loops, raw=True)
g2_type = [100, 200, 700] # curve, surface, volume identifiers
classes = [Curve, Surface, Volume]
g2_generators = {120:line, 130:circle, 140:ellipse,
260:cylinder, 292:disc, 270:sphere, 290:torus, 250:plane,
210:bounded_surface, 261:surface_of_linear_extrusion} #, 280:cone
def __init__(self, filename):
if filename[-3:] != '.g2':
filename += '.g2'
self.filename = filename
self.trimming_curves = []
def __enter__(self):
return self
def write(self, obj):
if not hasattr(self, 'fstream'):
self.onlywrite = True
self.fstream = open(self.filename, 'w')
if not self.onlywrite:
raise IOError('Could not write to file %s' % (self.filename))
"""Write the object in GoTools format. """
if isinstance(obj[0], SplineObject): # input SplineModel or list
for o in obj:
self.write(o)
return
for i in range(obj.pardim):
if obj.periodic(i):
obj = obj.split(obj.start(i), i)
self.fstream.write('{} 1 0 0\n'.format(G2.g2_type[obj.pardim-1]))
self.fstream.write('{} {}\n'.format(obj.dimension, int(obj.rational)))
for b in obj.bases:
self.fstream.write('%i %i\n' % (len(b.knots) - b.order, b.order))
self.fstream.write(' '.join('%.16g' % k for k in b.knots))
self.fstream.write('\n')
for cp in obj:
self.fstream.write(' '.join('%.16g' % x for x in cp))
self.fstream.write('\n')
def read(self):
if not hasattr(self, 'fstream'):
self.onlywrite = False
self.fstream = open(self.filename, 'r')
if self.onlywrite:
raise IOError('Could not read from file %s' % (self.filename))
result = []
for line in self.fstream:
line = line.strip()
if not line:
continue
# read object type
objtype, major, minor, patch = map(int, line.split())
if (major, minor, patch) != (1, 0, 0):
raise IOError('Unknown G2 format')
# if obj type is in factory methods (cicle, torus etc), create it now
if objtype in G2.g2_generators:
constructor = getattr(self, G2.g2_generators[objtype].__name__)
result.append( constructor() )
continue
# for "normal" splines (Curves, Surfaces, Volumes) create it now
pardim = [i for i in range(len(G2.g2_type)) if G2.g2_type[i] == objtype]
if not pardim:
raise IOError('Unknown G2 object type {}'.format(objtype))
pardim = pardim[0] + 1
result.append(self.splines(pardim))
return result
def read_basis(self):
ncps, order = map(int, next(self.fstream).split())
kts = list(map(float, next(self.fstream).split()))
return BSplineBasis(order, kts, -1)
def __exit__(self, exc_type, exc_value, traceback):
if hasattr(self, 'fstream'):
self.fstream.close()