SCM Repository
Annotation of /branches/pure-cfg/test/vr-curv-quad.diderot
Parent Directory
|
Revision Log
Revision 1025 - (view) (download)
1 : | glk | 866 | // vr-curv |
2 : | // | ||
3 : | glk | 911 | // Demonstration of curvature-based transfer functions in volume rendering |
4 : | // | ||
5 : | // process output with: | ||
6 : | // unu reshape -i mip.txt -s 4 480 345 | overrgb -i - -b 0.1 0.15 0.2 -g 1.2 -o - | unu quantize -b 8 -min 0 -max 255 -o vr-curv.png | ||
7 : | glk | 866 | |
8 : | glk | 1025 | //image(3)[] img = load("../data/quad-patches-pad.nrrd"); |
9 : | //field#2(3)[] F = img ⊛ bspln3; | ||
10 : | //image(2)[3] imgRGB = load("../data/hue.nhdr"); | ||
11 : | //field#0(2)[3] RGB = imgRGB ⊛ tent; | ||
12 : | glk | 911 | |
13 : | glk | 1025 | // HEY (bug) : trying to do load and convolve in one line, twice, |
14 : | // to produce effect of four lines above, but get: | ||
15 : | // /usr/bin/gcc -std=gnu99 -m64 -c -DNDEBUG -g -O2 -Wformat -Wreturn-type -Wuninitialized -D_THREAD_SAFE -I/Users/gk/diderot/diderot/pure-cfg/src/include -I/Users/gk/diderot/teem/arch/darwin.64/include vr-curv-quad.c | ||
16 : | // vr-curv-quad.c:30: error: conflicting types for ‘G__t’ | ||
17 : | // vr-curv-quad.c:29: error: previous declaration of ‘G__t’ was here | ||
18 : | // vr-curv-quad.c: In function ‘Diderot_InitGlobals’: | ||
19 : | field#2(3)[] F = load("../data/quad-patches-pad.nrrd") ⊛ bspln3; | ||
20 : | field#0(2)[3] RGB = load("../data/hue.nhdr") ⊛ tent; | ||
21 : | |||
22 : | glk | 913 | // set camera, image, and rendering parameters |
23 : | glk | 866 | vec3 camEye = [183.021, -411.857, 686.458]; |
24 : | vec3 camAt = [55.4037, 40.2407, 4.32732]; | ||
25 : | vec3 camUp = [0.0, 0.0, 1.0]; | ||
26 : | real camNear = -42.0; | ||
27 : | real camFar = 42.0; | ||
28 : | real camFOV = 7.15; | ||
29 : | int imgResU = 480; | ||
30 : | int imgResV = 345; | ||
31 : | glk | 911 | real rayStep = 0.2; |
32 : | real valOpacMin = -0.25; | ||
33 : | real valOpacMax = -0.10; | ||
34 : | glk | 913 | vec3 lightVspDir = [0.9, -2.5, -2.5]; |
35 : | vec3 lightRGB = [1.0, 1.0, 1.0]; | ||
36 : | real phongKa = 0.1; | ||
37 : | real phongKd = 0.7; | ||
38 : | real phongKs = 0.3; | ||
39 : | real phongSp = 45.0; | ||
40 : | glk | 866 | |
41 : | glk | 913 | // (boilerplate) computation of camera and light info |
42 : | glk | 866 | real camDist = |camAt - camEye|; |
43 : | real camVspNear = camNear + camDist; | ||
44 : | real camVspFar = camFar + camDist; | ||
45 : | vec3 camN = normalize(camAt - camEye); | ||
46 : | vec3 camU = normalize(camN × camUp); | ||
47 : | vec3 camV = camN × camU; | ||
48 : | real camVmax = tan(camFOV*3.1415926536/360.0)*camDist; | ||
49 : | real camUmax = camVmax*real(imgResU)/real(imgResV); | ||
50 : | glk | 911 | vec3 lightDir = normalize(lightVspDir[0]*camU + |
51 : | lightVspDir[1]*camV + | ||
52 : | lightVspDir[2]*camN); | ||
53 : | glk | 866 | |
54 : | glk | 911 | strand RayCast (int ui, int vi) { |
55 : | glk | 866 | real rayU = lerp(-camUmax, camUmax, -0.5, real(ui), real(imgResU)-0.5); |
56 : | real rayV = lerp(-camVmax, camVmax, -0.5, real(vi), real(imgResV)-0.5); | ||
57 : | vec3 rayVec = (camDist*camN + rayU*camU + rayV*camV)/camDist; | ||
58 : | |||
59 : | real rayN = camVspNear; | ||
60 : | glk | 911 | real rayTransp = 1.0; |
61 : | vec3 rayRGB = [0.0, 0.0, 0.0]; | ||
62 : | output vec4 outRGBA = [0.0, 0.0, 0.0, 0.0]; | ||
63 : | glk | 866 | |
64 : | glk | 911 | update { |
65 : | glk | 866 | vec3 pos = camEye + rayN*rayVec; |
66 : | if (inside (pos,F)) { | ||
67 : | glk | 911 | real val = F(pos); |
68 : | glk | 866 | if (val > valOpacMin) { // we have some opacity |
69 : | glk | 1025 | vec3 grad = -∇F(pos); |
70 : | vec3 norm = normalize(grad); | ||
71 : | glk | 911 | // begin curvature computation |
72 : | tensor[3,3] H = -∇(∇F)(pos); | ||
73 : | glk | 866 | tensor[3,3] P = identity[3] - norm⊗norm; |
74 : | tensor[3,3] G = (P•H•P)/|grad|; | ||
75 : | real disc = max(0.0, sqrt(2.0*|G|^2 - trace(G)^2)); | ||
76 : | glk | 911 | real k1 = (trace(G) + disc)/2.0; |
77 : | real k2 = (trace(G) - disc)/2.0; | ||
78 : | // finished curvature computation; begin finding sample RGBA | ||
79 : | k1 = max(-1.0, min(1.0, 4.2*k1)); | ||
80 : | k2 = max(-1.0, min(1.0, 4.2*k2)); | ||
81 : | vec3 matRGB = RGB([k1,k2]) if inside([k1,k2],RGB) | ||
82 : | else [1.0,0.0,1.0]; | ||
83 : | real alpha = min(1.0, lerp(0.0, 1.0, | ||
84 : | valOpacMin, val, valOpacMax)); | ||
85 : | real ld = max(0.0, norm • lightDir); | ||
86 : | real hd = max(0.0, norm • normalize(lightDir + | ||
87 : | normalize(camEye - pos))); | ||
88 : | // Phong shading | ||
89 : | vec3 pntRGB = (phongKa*matRGB | ||
90 : | + phongKd*ld*modulate(matRGB, lightRGB) | ||
91 : | + phongKs*hd^phongSp*lightRGB); | ||
92 : | // composite with existing ray color and transparency | ||
93 : | rayRGB = rayRGB + rayTransp*alpha*pntRGB; | ||
94 : | rayTransp = rayTransp*(1.0 - alpha); | ||
95 : | glk | 866 | } |
96 : | } | ||
97 : | glk | 911 | if (rayTransp < 0.01) { // early ray termination |
98 : | rayTransp = 0.0; | ||
99 : | outRGBA = [rayRGB[0], rayRGB[1], rayRGB[2], 1.0 - rayTransp]; | ||
100 : | glk | 866 | stabilize; |
101 : | } | ||
102 : | if (rayN > camVspFar) { | ||
103 : | glk | 911 | outRGBA = [rayRGB[0], rayRGB[1], rayRGB[2], 1.0 - rayTransp]; |
104 : | glk | 866 | stabilize; |
105 : | } | ||
106 : | rayN = rayN + rayStep; | ||
107 : | } | ||
108 : | |||
109 : | /* render: output maxval */ | ||
110 : | } | ||
111 : | |||
112 : | initially [ RayCast(ui, vi) | vi in 0..(imgResV-1), ui in 0..(imgResU-1) ]; |
root@smlnj-gforge.cs.uchicago.edu | ViewVC Help |
Powered by ViewVC 1.0.0 |