Changeset 1616
- Timestamp:
- 02/11/10 11:40:12 (6 months ago)
- Location:
- trunk/yt
- Files:
-
- 5 modified
-
commands.py (modified) (3 diffs)
-
extensions/volume_rendering/software_sampler.py (modified) (2 diffs)
-
lagos/HDF5LightReader.c (modified) (1 diff)
-
_amr_utils/FixedInterpolator.c (modified) (1 diff)
-
_amr_utils/VolumeIntegrator.pyx (modified) (12 diffs)
Legend:
- Unmodified
- Added
- Removed
-
trunk/yt/commands.py
r1539 r1616 28 28 from yt.recipes import _fix_pf 29 29 import yt.cmdln as cmdln 30 import optparse, os, os.path, math, sys 30 import optparse, os, os.path, math, sys, time, subprocess 31 31 32 32 _common_options = dict( … … 206 206 return arg_iterate 207 207 208 def _get_vcs_type(path): 209 if os.path.exists(os.path.join(path, ".hg")): 210 return "hg" 211 if os.path.exists(os.path.join(path, ".svn")): 212 return "svn" 213 return None 214 215 def _get_svn_version(path): 216 p = subprocess.Popen(["svn", "info", path], stdout = subprocess.PIPE, 217 stderr = subprocess.STDOUT) 218 stdout, stderr = p.communicate() 219 return stdout 220 221 def _update_svn(path): 222 f = open(os.path.join(path, "yt_updater.log"), "a") 223 f.write("\n\nUPDATE PROCESS: %s\n\n" % (time.asctime())) 224 p = subprocess.Popen(["svn", "up", path], stdout = subprocess.PIPE, 225 stderr = subprocess.STDOUT) 226 stdout, stderr = p.communicate() 227 f.write(stdout) 228 f.write("\n\n") 229 if p.returncode: 230 print "BROKEN: See %s" % (os.path.join(path, "yt_updater.log")) 231 sys.exit(1) 232 f.write("Rebuilding modules\n\n") 233 p = subprocess.Popen([sys.executable, "setup.py", "build_ext", "-i"], cwd=path, 234 stdout = subprocess.PIPE, stderr = subprocess.STDOUT) 235 stdout, stderr = p.communicate() 236 f.write(stdout) 237 f.write("\n\n") 238 if p.returncode: 239 print "BROKEN: See %s" % (os.path.join(path, "yt_updater.log")) 240 sys.exit(1) 241 f.write("Successful!\n") 242 243 def _update_hg(path): 244 from mercurial import hg, ui, commands 245 f = open(os.path.join(path, "yt_updater.log"), "a") 246 u = ui.ui() 247 u.pushbuffer() 248 config_fn = os.path.join(path, ".hg", "hgrc") 249 print "Reading configuration from ", config_fn 250 u.readconfig(config_fn) 251 repo = hg.repository(u, path) 252 commands.pull(u, repo) 253 f.write(u.popbuffer()) 254 f.write("\n\n") 255 u.pushbuffer() 256 commands.identify(u, repo) 257 if "+" in u.popbuffer(): 258 print "Can't rebuild modules by myself." 259 print "You will have to do this yourself. Here's a sample commands:" 260 print 261 print " $ cd %s" % (path) 262 print " $ hg up" 263 print " $ %s setup.py develop" % (sys.executable) 264 sys.exit(1) 265 f.write("Rebuilding modules\n\n") 266 p = subprocess.Popen([sys.executable, "setup.py", "build_ext", "-i"], cwd=path, 267 stdout = subprocess.PIPE, stderr = subprocess.STDOUT) 268 stdout, stderr = p.communicate() 269 f.write(stdout) 270 f.write("\n\n") 271 if p.returncode: 272 print "BROKEN: See %s" % (os.path.join(path, "yt_updater.log")) 273 sys.exit(1) 274 f.write("Successful!\n") 275 276 def _get_hg_version(path): 277 from mercurial import hg, ui, commands 278 u = ui.ui() 279 u.pushbuffer() 280 repo = hg.repository(u, path) 281 commands.identify(u, repo) 282 return u.popbuffer() 283 284 _vcs_identifier = dict(svn = _get_svn_version, 285 hg = _get_hg_version) 286 _vcs_updater = dict(svn = _update_svn, 287 hg = _update_hg) 288 208 289 class YTCommands(cmdln.Cmdln): 209 290 name="yt" … … 220 301 """ 221 302 self.cmdloop() 303 304 @cmdln.option("-u", "--update-source", action="store_true", 305 default = False, 306 help="Update the yt installation, if able") 307 def do_instinfo(self, subcmd, opts): 308 """ 309 ${cmd_name}: Get some information about the yt installation 310 311 ${cmd_usage} 312 ${cmd_option_list} 313 """ 314 import pkg_resources 315 yt_provider = pkg_resources.get_provider("yt") 316 path = os.path.dirname(yt_provider.module_path) 317 print 318 print "yt module located at:" 319 print " %s" % (path) 320 if "site-packages" not in path: 321 vc_type = _get_vcs_type(path) 322 vstring = _vcs_identifier[vc_type](path) 323 print 324 print "The current version of the code is:" 325 print 326 print "---" 327 print vstring 328 print "---" 329 print 330 print "This installation CAN be automatically updated." 331 if opts.update_source: 332 _vcs_updater[vc_type](path) 333 print "Updated successfully." 222 334 223 335 @add_cmd_options(['outputfn','bn','thresh','dm_only','skip']) -
trunk/yt/extensions/volume_rendering/software_sampler.py
r1607 r1616 75 75 hit = 0 76 76 tnow = time.time() 77 every = na.ceil(len(partitioned_grids) / 100.0)78 77 79 78 vp = VectorPlane(vectors, norm_vec, back_center, … … 90 89 tfp.ns = nsamples 91 90 92 pbar = get_pbar("Ray casting ", len(partitioned_grids)) 91 92 total_cells = sum(na.prod(g.my_data.shape) for g in partitioned_grids) 93 pbar = get_pbar("Ray casting ", total_cells) 94 total_cells = 0 93 95 for i,g in enumerate(partitioned_grids[ind]): 94 if (i % every) == 0: 95 pbar.update(i) 96 pbar.update(total_cells) 96 97 pos = g.cast_plane(tfp, vp) 98 total_cells += na.prod(g.my_data.shape) 97 99 pbar.finish() 98 100 -
trunk/yt/lagos/HDF5LightReader.c
r1602 r1616 34 34 35 35 #include "numpy/ndarrayobject.h" 36 37 #ifndef npy_float128 38 #define npy_float128 npy_longdouble 39 #endif 36 40 37 41 -
trunk/yt/_amr_utils/FixedInterpolator.c
r1579 r1616 98 98 normval += grad[i]*grad[i]; 99 99 } 100 normval = sqrt(normval); 101 for (i = 0; i < 3; i++) grad[i] /= -normval; 102 //fprintf(stderr, "Normval: %0.3lf %0.3lf %0.3lf %0.3lf\n", 103 // normval, grad[0], grad[1], grad[2]); 100 if (normval != 0.0){ 101 normval = sqrt(normval); 102 for (i = 0; i < 3; i++) grad[i] /= -normval; 103 //fprintf(stderr, "Normval: %0.3lf %0.3lf %0.3lf %0.3lf\n", 104 // normval, grad[0], grad[1], grad[2]); 105 }else{ 106 grad[0]=grad[1]=grad[2]=0.0; 107 } 104 108 } -
trunk/yt/_amr_utils/VolumeIntegrator.pyx
r1607 r1616 46 46 47 47 cdef inline int iclip(int i, int a, int b): 48 return imin(imax(i, a), b) 48 if i < a: return a 49 if i > b: return b 50 return i 49 51 50 52 cdef inline np.float64_t fclip(np.float64_t f, … … 77 79 cdef int nbins 78 80 cdef public int ns 79 cdef np.float64_t dbin 81 cdef np.float64_t dbin, idbin 80 82 cdef np.float64_t light_color[3] 81 83 cdef np.float64_t light_dir[3] … … 97 99 self.nbins = tf_obj.nbins 98 100 self.dbin = (self.x_bounds[1] - self.x_bounds[0])/self.nbins 101 self.idbin = 1.0/self.dbin 99 102 self.light_color[0] = tf_obj.light_color[0] 100 103 self.light_color[1] = tf_obj.light_color[1] … … 109 112 self.use_light = tf_obj.use_light 110 113 114 @cython.boundscheck(False) 115 @cython.wraparound(False) 116 cdef void interpolate(self, np.float64_t dv, np.float64_t *trgba): 117 cdef int bin_id, channel 118 cdef np.float64_t bv, dy, dd, tf 119 bin_id = <int> ((dv - self.x_bounds[0]) * self.idbin) 120 # Recall that linear interpolation is y0 + (x-x0) * dx/dy 121 dd = dv-(self.x_bounds[0] + bin_id * self.dbin) # x - x0 122 for channel in range(4): 123 bv = self.vs[channel][bin_id] # This is x0 124 dy = self.vs[channel][bin_id+1]-bv # dy 125 # This is our final value for transfer function on the entering face 126 trgba[channel] = bv+dd*dy*self.idbin 127 128 @cython.boundscheck(False) 129 @cython.wraparound(False) 111 130 cdef void eval_transfer(self, np.float64_t dt, np.float64_t dv, 112 131 np.float64_t *rgba, np.float64_t *grad): 113 132 cdef int i 114 cdef int bin_id 115 cdef np.float64_t tf, trgba[4], bv, dx, dy, dd, ta, dot_prod 116 dx = self.dbin 117 133 cdef np.float64_t ta, tf, trgba[4], dot_prod 134 self.interpolate(dv, trgba) 118 135 # get source alpha first 119 136 # First locate our points 120 bin_id = iclip(<int> floor((dv - self.x_bounds[0]) / dx),121 0, self.nbins-2)122 # Recall that linear interpolation is y0 + (x-x0) * dx/dy123 bv = self.vs[3][bin_id] # This is x0124 dy = self.vs[3][bin_id+1]-bv # dy125 dd = dv-(self.x_bounds[0] + bin_id * dx) # x - x0126 # This is our final value for transfer function on the entering face127 tf = bv+dd*(dy/dx)128 ta = tf # Store the source alpha129 137 dot_prod = 0.0 130 for i in range(3): 131 dot_prod += self.light_dir[i] * grad[i] 132 #print dot_prod, grad[0], grad[1], grad[2] 133 dot_prod = fmax(0.0, dot_prod) 134 for i in range(3): 135 # Recall that linear interpolation is y0 + (x-x0) * dx/dy 136 bv = self.vs[i][bin_id] # This is x0 137 dy = self.vs[i][bin_id+1]-bv # dy 138 dd = dv-(self.x_bounds[0] + bin_id * dx) # x - x0 139 # This is our final value for transfer function on the entering face 140 tf = bv+dd*(dy/dx) + dot_prod * self.light_color[i] 141 # alpha blending 142 rgba[i] += (1. - rgba[3])*ta*tf*dt 143 #update alpha 144 rgba[3] += (1. - rgba[3])*ta*dt 145 # We should really do some alpha blending. 146 # Front to back blending is defined as: 147 # dst.rgb = dst.rgb + (1 - dst.a) * src.a * src.rgb 148 # dst.a = dst.a + (1 - dst.a) * src.a 138 if self.use_light: 139 for i in range(3): 140 dot_prod += self.light_dir[i] * grad[i] 141 dot_prod = fmax(0.0, dot_prod) 142 for i in range(3): 143 trgba[i] += dot_prod*self.light_color[i] 144 # alpha blending 145 ta = (1.0 - rgba[3])*dt*trgba[3] 146 for i in range(4): 147 rgba[i] += ta*trgba[i] 149 148 150 149 cdef class VectorPlane: … … 182 181 self.pdy = (self.bounds[3] - self.bounds[2])/self.nv 183 182 183 @cython.boundscheck(False) 184 @cython.wraparound(False) 184 185 cdef void get_start_stop(self, np.float64_t *ex, int *rv): 185 186 # Extrema need to be re-centered 186 187 cdef np.float64_t cx, cy 188 cdef int i 187 189 cx = cy = 0.0 188 190 for i in range(3): … … 216 218 cdef np.float64_t right_edge[3] 217 219 cdef np.float64_t dds[3] 220 cdef np.float64_t idds[3] 218 221 cdef public np.float64_t min_dds 219 222 cdef int dims[3] … … 235 238 self.dims[i] = dims[i] 236 239 self.dds[i] = (self.right_edge[i] - self.left_edge[i])/dims[i] 240 self.idds[i] = 1.0/self.dds[i] 237 241 self.my_data = data 238 242 self.data = <np.float64_t*> data.data … … 294 298 cdef np.float64_t intersect_t = 1.0 295 299 cdef np.float64_t intersect[3], tmax[3], tdelta[3] 296 cdef np.float64_t enter_t, dist, alpha, dt 300 cdef np.float64_t enter_t, dist, alpha, dt, exit_t 297 301 cdef np.float64_t tr, tl, temp_x, temp_y, dv 298 302 for i in range(3): … … 328 332 cur_ind[i] = <int> floor((intersect[i] + 329 333 step[i]*1e-8*self.dds[i] - 330 self.left_edge[i]) /self.dds[i])334 self.left_edge[i])*self.idds[i]) 331 335 tmax[i] = (((cur_ind[i]+step[i])*self.dds[i])+ 332 336 self.left_edge[i]-v_pos[i])/v_dir[i] 337 # This deals with the asymmetry in having our indices refer to the 338 # left edge of a cell, but the right edge of the brick being one 339 # extra zone out. 333 340 if cur_ind[i] == self.dims[i] and step[i] < 0: 334 341 cur_ind[i] = self.dims[i] - 1 … … 354 361 if tmax[0] < tmax[1]: 355 362 if tmax[0] < tmax[2]: 356 self.sample_values(v_pos, v_dir, enter_t, tmax[0], cur_ind, 363 exit_t = fmin(tmax[0], 1.0) 364 self.sample_values(v_pos, v_dir, enter_t, exit_t, cur_ind, 357 365 rgba, tf) 358 366 cur_ind[0] += step[0] 359 dt = fmin(tmax[0], 1.0) - enter_t360 367 enter_t = tmax[0] 361 368 tmax[0] += tdelta[0] 362 369 else: 363 self.sample_values(v_pos, v_dir, enter_t, tmax[2], cur_ind, 370 exit_t = fmin(tmax[2], 1.0) 371 self.sample_values(v_pos, v_dir, enter_t, exit_t, cur_ind, 364 372 rgba, tf) 365 373 cur_ind[2] += step[2] 366 dt = fmin(tmax[2], 1.0) - enter_t367 374 enter_t = tmax[2] 368 375 tmax[2] += tdelta[2] 369 376 else: 370 377 if tmax[1] < tmax[2]: 371 self.sample_values(v_pos, v_dir, enter_t, tmax[1], cur_ind, 378 exit_t = fmin(tmax[1], 1.0) 379 self.sample_values(v_pos, v_dir, enter_t, exit_t, cur_ind, 372 380 rgba, tf) 373 381 cur_ind[1] += step[1] 374 dt = fmin(tmax[1], 1.0) - enter_t375 382 enter_t = tmax[1] 376 383 tmax[1] += tdelta[1] 377 384 else: 378 self.sample_values(v_pos, v_dir, enter_t, tmax[2], cur_ind, 385 exit_t = fmin(tmax[2], 1.0) 386 self.sample_values(v_pos, v_dir, enter_t, exit_t, cur_ind, 379 387 rgba, tf) 380 388 cur_ind[2] += step[2] 381 dt = fmin(tmax[2], 1.0) - enter_t382 389 enter_t = tmax[2] 383 390 tmax[2] += tdelta[2] … … 385 392 return hit 386 393 394 @cython.boundscheck(False) 395 @cython.wraparound(False) 387 396 cdef void sample_values(self, 388 397 np.float64_t v_pos[3], … … 394 403 TransferFunctionProxy tf): 395 404 cdef np.float64_t cp[3], dp[3], temp, dt, t, dv 396 cdef np.float64_t grad[3] 405 cdef np.float64_t grad[3], ds[3] 406 grad[0] = grad[1] = grad[2] = 0.0 397 407 cdef int dti, i 398 dt = (exit_t - enter_t) / (tf.ns-1) # five samples, so divide by four 399 for dti in range(tf.ns - 1): 400 t = enter_t + dt * dti 408 dt = (exit_t - enter_t) / (tf.ns) # 4 samples should be dt=0.25 409 for i in range(3): 410 temp = ci[i] * self.dds[i] + self.left_edge[i] 411 dp[i] = (enter_t + 0.5 * dt) * v_dir[i] + v_pos[i] - temp 412 dp[i] *= self.idds[i] 413 ds[i] = v_dir[i] * self.idds[i] * dt 414 for dti in range(tf.ns): 401 415 for i in range(3): 402 cp[i] = v_pos[i] + t * v_dir[i] 403 dp[i] = fclip(fmod(cp[i], self.dds[i])/self.dds[i], 0, 1.0) 416 dp[i] += ds[i] 404 417 dv = trilinear_interpolate(self.dims, ci, dp, self.data) 418 if not ((dv > tf.x_bounds[0]) and (dv < tf.x_bounds[1])): 419 continue 405 420 if tf.use_light == 1: 406 421 eval_gradient(self.dims, ci, dp, self.data, grad)
