SCM Repository
[diderot] / branches / pure-cfg / test / vr-lite-cam.diderot |
View of /branches/pure-cfg/test/vr-lite-cam.diderot
Parent Directory
|
Revision Log
Revision 1302 -
(download)
(annotate)
Fri Jun 10 00:00:05 2011 UTC (9 years, 10 months ago) by jhr
File size: 4515 byte(s)
Fri Jun 10 00:00:05 2011 UTC (9 years, 10 months ago) by jhr
File size: 4515 byte(s)
Minor tweaks
// vr-lite-cam // // simple volume renderer, now with a single directional light, // Phong shading, and the new camera. // // The main thing missing from this is a more general transfer function // // process output with: // 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-lite-cam.png //string dataFile = "../data/txs-pad3.nrrd"; //vec3 camEye = [25.0, 15.0, 10.0]; //vec3 camAt = [3.0, 3.0, 3.0]; //vec3 camUp = [0.0, 0.0, 1.0]; //real camNear = -5.0; //real camFar = 5.0; //real camFOV = 16.0; //int imgResU = 480; //int imgResV = 345; //real rayStep = 0.1; //real valOpacMin = 0.15; // highest value with opacity 0.0 //real valOpacMax = 0.20; // lowest value with opacity 1.0 string dataFile = "../data/vfrhand-nohip.nhdr"; vec3 camEye = [127.331, -1322.05, 272.53]; vec3 camAt = [63.0, 82.6536, 98.0]; vec3 camUp = [0.9987, 0.0459166, -0.0221267]; real camNear = -78.0; real camFar = 78.0; real camFOV = 5.0; int imgResU = 480; int imgResV = 345; real rayStep = 0.1; real valOpacMin = 400.0; // 400.0 for skin, 1150.0 for bone real valOpacMax = 700.0; // 700.0 for skin, 1450.0 for bone real camDist = |camAt - camEye|; real camVspNear = camDist + camNear; real camVspFar = camDist + camFar; vec3 camN = normalize(camAt - camEye); vec3 camU = normalize(camN × camUp); vec3 camV = camN × camU; real camVmax = tan(camFOV*π/360.0)*camDist; real camUmax = camVmax*real(imgResU)/real(imgResV); vec3 lightVspDir = [0.9, -1.0, -2.5]; vec3 lightDir = normalize(lightVspDir[0]*camU + lightVspDir[1]*camV + lightVspDir[2]*camN); //vec3 lightRGB = [1.0, 1.0, 1.0]; real phongKa = 0.05; real phongKd = 0.65; real phongKs = 0.45; real phongSp = 50.0; //field#4(3)[] F = bspln5 ⊛ load(dataFile); field#2(3)[] F = bspln3 ⊛ load(dataFile); //field#1(3)[] F = ctmr ⊛ load(dataFile); //field#1(3)[] F = c1tent ⊛ load(dataFile); strand RayCast (int ui, int vi) { real rayU = lerp(-camUmax, camUmax, -0.5, real(ui), real(imgResU)-0.5); real rayV = lerp(-camVmax, camVmax, -0.5, real(vi), real(imgResV)-0.5); vec3 rayVec = (camDist*camN + rayU*camU + rayV*camV)/camDist; vec3 toEye = normalize(-rayVec); real rayN = camVspNear; // ########## BEGIN per-ray initialization real rayTransp = 1.0; // vec3 rayRGB = [0.0, 0.0, 0.0]; real rayGrey = 0.0; output vec4 outRGBA = [0.0, 0.0, 0.0, 0.0]; // ########## END per-ray initialization update { vec3 rayPos = camEye + rayN*rayVec; if (inside (rayPos,F)) { // ########## BEGIN per-sample code real val = F(rayPos); if (val > valOpacMin) { // we have some opacity vec3 grad = ∇F(rayPos); // (here as easy target for optimization) vec3 norm = normalize(-∇F(rayPos)); real alpha = min(1.0, lerp(0.0, 1.0, valOpacMin, val, valOpacMax)); real ld = max(0.0, norm • lightDir); real hd = max(0.0, norm • normalize(lightDir + toEye)); // contrived assignment of RGB from XYZ, only sensible // for txs-pad3.nrrd, not for vfrhand-nohip.nhdr!! // vec3 matRGB = [lerp(0.2, 1.0, 1.0, rayPos[0], 8.0), // lerp(0.2, 1.0, 1.0, rayPos[1], 8.0), // lerp(0.2, 1.0, 1.0, rayPos[2], 8.0)]; // vec3 pntRGB = (phongKa*matRGB // + phongKd*ld*modulate(matRGB, lightRGB) // + phongKs*hd^phongSp*lightRGB); // rayRGB += rayTransp*alpha*pntRGB; real mat = lerp(0.2, 1.0, 1.0, rayPos[0], 8.0); real pnt = phongKa * mat + phongKd * ld * mat + phongKs * hd^phongSp; rayGrey += rayTransp*alpha*pnt; rayTransp = rayTransp*(1.0 - alpha); } // ########## END per-sample code } if (rayTransp < 0.01) { // early ray termination rayTransp = 0.0; // outRGBA = [rayRGB[0], rayRGB[1], rayRGB[2], 1.0-rayTransp]; /* FIXME */ outRGBA = [rayGrey, rayGrey, rayGrey, 1.0-rayTransp]; /* FIXME */ stabilize; } if (rayN > camVspFar) { // outRGBA = [rayRGB[0], rayRGB[1], rayRGB[2], 1.0-rayTransp]; /* FIXME */ outRGBA = [rayGrey, rayGrey, rayGrey, 1.0-rayTransp]; /* FIXME */ stabilize; } rayN += rayStep; } stabilize { // outRGBA = [rayRGB[0], rayRGB[1], rayRGB[2], 1.0-rayTransp]; } } 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 |