polaroid-pp

Schlieren and contour plot tool
git clone https://git.0xfab.ch/polaroid-pp.git
Log | Files | Refs | Submodules | README | LICENSE

commit 1b42cb706511384eaa06459804ddfc89d7d44f02
parent 9a34d7f8b39b6aeb5d1ffd26ca0aff6fba794bfd
Author: Fabian Wermelinger <fabianw@mavt.ethz.ch>
Date:   Fri, 29 Apr 2016 17:07:45 +0200

added various cartridges

Diffstat:
Aapps/polaroidCamera/BoundedLogNormalizerCartridge.h | 51+++++++++++++++++++++++++++++++++++++++++++++++++++
Aapps/polaroidCamera/BoundedNormalizerCartridge.h | 44++++++++++++++++++++++++++++++++++++++++++++
Aapps/polaroidCamera/BoundedTransmissionCartridge.h | 42++++++++++++++++++++++++++++++++++++++++++
Mapps/polaroidCamera/Cartridges.h | 5+++++
Aapps/polaroidCamera/LogNormalizerCartridge.h | 49+++++++++++++++++++++++++++++++++++++++++++++++++
Mapps/polaroidCamera/NormalizerCartridge.h | 41++++++++++++++++++++++++++++-------------
Mapps/polaroidCamera/SceneProcessor.cpp | 16+++++++++++++---
Aapps/polaroidCamera/SchlierenCartridge.cpp | 67+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Mapps/polaroidCamera/SchlierenCartridge.h | 55+++++++++++++++++++++++++++++++++++++------------------
Mapps/polaroidCamera/TransmissionCartridge.h | 2+-
Minclude/Cartridge.h | 2+-
Minclude/Types.h | 12++++++++----
Msrc/Polaroid.cpp | 6++++++
13 files changed, 352 insertions(+), 40 deletions(-)

diff --git a/apps/polaroidCamera/BoundedLogNormalizerCartridge.h b/apps/polaroidCamera/BoundedLogNormalizerCartridge.h @@ -0,0 +1,51 @@ +// File : BoundedLogNormalizerCartridge.h +// Date : Fri 29 Apr 2016 09:39:39 AM CEST +// Author : Fabian Wermelinger +// Description: Bounded Log Data Normalizer Cartridge +// Copyright 2016 ETH Zurich. All Rights Reserved. +#ifndef BOUNDEDLOGNORMALIZERCARTRIDGE_H_QF9H4ZPF +#define BOUNDEDLOGNORMALIZERCARTRIDGE_H_QF9H4ZPF + +#include <cassert> +#include <algorithm> +#include <cmath> +#include "Cartridge.h" + +class BoundedLogNormalizerCartridge : public Cartridge +{ +public: + BoundedLogNormalizerCartridge(ArgumentParser& parser) : Cartridge(parser) {} + + virtual void capture(PhotoPaper& photo, Slice& data) + { + const Real upper = m_parser("-upper_bound").asDouble(10.0); + const Real lower = m_parser("-lower_bound").asDouble(1.0); + + photo.resize(data.width(), data.height()); + + // set description + string desc("2D_Bounded_Log_Normalized"); + photo.set_description(desc.c_str()); + + // compute min/max for shader + const Real dataMinInv = 1.0/lower; + const Real fac = 1.0/log(upper*dataMinInv); + assert(!isnan(dataMinInv)); + assert(!isnan(fac)); + + // pixel shader + for (int h=0; h < data.height(); ++h) + for (int w=0; w < data.width(); ++w) + { + const Real bound = std::max(static_cast<Real>(0.0), std::min(static_cast<Real>(1.0), log(data(w,h)*dataMinInv)*fac)); + assert(!isnan(bound)); + photo.set_pixel(w, h, bound); + } + + photo.write(); + } + + virtual void compute(Slice& data) { } +}; + +#endif /* BOUNDEDLOGNORMALIZERCARTRIDGE_H_QF9H4ZPF */ diff --git a/apps/polaroidCamera/BoundedNormalizerCartridge.h b/apps/polaroidCamera/BoundedNormalizerCartridge.h @@ -0,0 +1,44 @@ +// File : BoundedNormalizerCartridge.h +// Date : Fri 29 Apr 2016 09:33:09 AM CEST +// Author : Fabian Wermelinger +// Description: Bounded Normalizer Cartridge +// Copyright 2016 ETH Zurich. All Rights Reserved. +#ifndef BOUNDEDNORMALIZERCARTRIDGE_H_DYS5RKGA +#define BOUNDEDNORMALIZERCARTRIDGE_H_DYS5RKGA + +#include <algorithm> +#include "Cartridge.h" + +class BoundedNormalizerCartridge : public Cartridge +{ +public: + BoundedNormalizerCartridge(ArgumentParser& parser) : Cartridge(parser) {} + + virtual void capture(PhotoPaper& photo, Slice& data) + { + const Real upper = m_parser("-upper_bound").asDouble(1.0); + const Real lower = m_parser("-lower_bound").asDouble(0.0); + + photo.resize(data.width(), data.height()); + + // set description + string desc("2D_Bounded_Normalized"); + photo.set_description(desc.c_str()); + + const Real data_normInv = 1.0 / (upper - lower); + + // pixel shader + for (int h=0; h < data.height(); ++h) + for (int w=0; w < data.width(); ++w) + { + const Real bound = std::max(static_cast<Real>(0.0),std::min(static_cast<Real>(1.0), (data(w,h)-lower)*data_normInv)); + photo.set_pixel(w, h, bound); + } + + photo.write(); + } + + virtual void compute(Slice& data) { } +}; + +#endif /* BOUNDEDNORMALIZERCARTRIDGE_H_DYS5RKGA */ diff --git a/apps/polaroidCamera/BoundedTransmissionCartridge.h b/apps/polaroidCamera/BoundedTransmissionCartridge.h @@ -0,0 +1,42 @@ +// File : BoundedTransmissionCartridge.h +// Date : Fri 29 Apr 2016 09:22:10 AM CEST +// Author : Fabian Wermelinger +// Description: Data BoundedTransmission Cartridge +// Copyright 2016 ETH Zurich. All Rights Reserved. +#ifndef BOUNDEDTRANSMISSIONCARTRIDGE_H_YDXLKPJU +#define BOUNDEDTRANSMISSIONCARTRIDGE_H_YDXLKPJU + +#include <algorithm> +#include "Cartridge.h" + +class BoundedTransmissionCartridge : public Cartridge +{ +public: + BoundedTransmissionCartridge(ArgumentParser& parser) : Cartridge(parser) {} + + virtual void capture(PhotoPaper& photo, Slice& data) + { + const Real upper = m_parser("-upper_bound").asDouble(1.0); + const Real lower = m_parser("-lower_bound").asDouble(0.0); + + photo.resize(data.width(), data.height()); + + // set description + string desc("2D_Bounded_Transmission"); + photo.set_description(desc.c_str()); + + // pixel shader + for (int h=0; h < data.height(); ++h) + for (int w=0; w < data.width(); ++w) + { + const Real bound = std::max(lower,std::min(upper,data(w,h))); + photo.set_pixel(w, h, bound); + } + + photo.write(); + } + + virtual void compute(Slice& data) { } +}; + +#endif /* BOUNDEDTRANSMISSIONCARTRIDGE_H_YDXLKPJU */ diff --git a/apps/polaroidCamera/Cartridges.h b/apps/polaroidCamera/Cartridges.h @@ -8,5 +8,10 @@ #include "TransmissionCartridge.h" #include "NormalizerCartridge.h" +#include "LogNormalizerCartridge.h" +#include "BoundedTransmissionCartridge.h" +#include "BoundedNormalizerCartridge.h" +#include "BoundedLogNormalizerCartridge.h" +#include "SchlierenCartridge.h" #endif /* CARTRIDGES_H_LTWKNANR */ diff --git a/apps/polaroidCamera/LogNormalizerCartridge.h b/apps/polaroidCamera/LogNormalizerCartridge.h @@ -0,0 +1,49 @@ +// File : LogNormalizerCartridge.h +// Date : Fri 29 Apr 2016 09:05:30 AM CEST +// Author : Fabian Wermelinger +// Description: Log Normalizer Cartridge +// Copyright 2016 ETH Zurich. All Rights Reserved. +#ifndef LOGNORMALIZERCARTRIDGE_H_GEKZAQWI +#define LOGNORMALIZERCARTRIDGE_H_GEKZAQWI + +#include <cassert> +#include "NormalizerCartridge.h" + +class LogNormalizerCartridge : public NormalizerCartridge +{ +public: + LogNormalizerCartridge(ArgumentParser& parser) : NormalizerCartridge(parser) {} + + virtual void capture(PhotoPaper& photo, Slice& data) + { + photo.resize(data.width(), data.height()); + + // set description + string desc("2D_Log_Normalized"); + photo.set_description(desc.c_str()); + + // compute min/max for shader + if (!m_bComputed) + { + m_dataMin = data.min(); + m_dataMax = data.max(); + } + const Real dataMinInv = 1.0/m_dataMin; + const Real fac = 1.0/log(m_dataMax*dataMinInv); + assert(!isnan(dataMinInv)); + assert(!isnan(fac)); + + // pixel shader + for (int h=0; h < data.height(); ++h) + for (int w=0; w < data.width(); ++w) + { + const Real logData = log(data(w,h)*dataMinInv); + assert(!isnan(logData)); + photo.set_pixel(w, h, fac*logData); + } + + photo.write(); + } +}; + +#endif /* LOGNORMALIZERCARTRIDGE_H_GEKZAQWI */ diff --git a/apps/polaroidCamera/NormalizerCartridge.h b/apps/polaroidCamera/NormalizerCartridge.h @@ -6,41 +6,56 @@ #ifndef NORMALIZERCARTRIDGE_H_ONAUDSVX #define NORMALIZERCARTRIDGE_H_ONAUDSVX +#include <algorithm> +#include <cmath> #include "Cartridge.h" class NormalizerCartridge : public Cartridge { +protected: + Real m_dataMin; + Real m_dataMax; + public: - NormalizerCartridge(ArgumentParser& parser) : Cartridge(parser) {} + NormalizerCartridge(ArgumentParser& parser) : Cartridge(parser), m_dataMin(HUGE_VAL), m_dataMax(-HUGE_VAL) {} virtual void capture(PhotoPaper& photo, Slice& data) { photo.resize(data.width(), data.height()); // set description - string desc("2D Normalized Slice Data"); + string desc("2D_Normalized"); photo.set_description(desc.c_str()); // compute min/max for shader - Real dataMin = data(0,0); - Real dataMax = data(0,0); - for (int h=0; h < data.height(); ++h) - for (int w=0; w < data.width(); ++w) - { - dataMin = min(dataMin, data(w,h)); - dataMax = max(dataMax, data(w,h)); - } - const Real data_normInv = 1.0 / (dataMax - dataMin); + if (!m_bComputed) + { + m_dataMin = data.min(); + m_dataMax = data.max(); + } + const Real data_normInv = 1.0 / (m_dataMax - m_dataMin); // pixel shader for (int h=0; h < data.height(); ++h) for (int w=0; w < data.width(); ++w) - photo.set_pixel(w, h, (data(w,h) - dataMin) * data_normInv); + photo.set_pixel(w, h, (data(w,h) - m_dataMin) * data_normInv); photo.write(); } - virtual void compute(Slice& data) { } + virtual void compute(Slice& data) + { + m_dataMin = std::min(m_dataMin, data.min()); + m_dataMax = std::max(m_dataMax, data.max()); + m_bComputed = true; + } + + virtual void reset() + { + m_bComputed = false; + m_dataMin = HUGE_VAL; + m_dataMax = -HUGE_VAL; + } }; #endif /* NORMALIZERCARTRIDGE_H_ONAUDSVX */ diff --git a/apps/polaroidCamera/SceneProcessor.cpp b/apps/polaroidCamera/SceneProcessor.cpp @@ -18,6 +18,16 @@ void SceneProcessor::_prepare_cam() m_cartridge = new TransmissionCartridge(m_parser); else if (cart == "normalizer") m_cartridge = new NormalizerCartridge(m_parser); + else if (cart == "log_normalizer") + m_cartridge = new LogNormalizerCartridge(m_parser); + else if (cart == "bounded_transmission") + m_cartridge = new BoundedTransmissionCartridge(m_parser); + else if (cart == "bounded_normalizer") + m_cartridge = new BoundedNormalizerCartridge(m_parser); + else if (cart == "bounded_log_normalizer") + m_cartridge = new BoundedLogNormalizerCartridge(m_parser); + else if (cart == "schlieren") + m_cartridge = new SchlierenCartridge(m_parser); else { if (m_isroot) @@ -78,7 +88,7 @@ void SceneProcessor::process1212(const vector<char*>& scenes) mycam.capture(*m_photo); if (m_isroot) - printf("[Progress %3.1f %%]\n", static_cast<double>(i)/scenes.size()*100.0); + printf("[Progress %3.1f %%]\n", static_cast<double>(i+1)/scenes.size()*100.0); } } @@ -95,7 +105,7 @@ void SceneProcessor::process1122(const vector<char*>& scenes) cam.insertCartridge(m_cartridge); cam.computeScene(); if (m_isroot) - printf("[Scene Progress %3.1f %%]\n", static_cast<double>(i)/scenes.size()*100.0); + printf("[Scene Progress %3.1f %%]\n", static_cast<double>(i+1)/scenes.size()*100.0); } // 2.)apply pixel shader and write photo @@ -111,6 +121,6 @@ void SceneProcessor::process1122(const vector<char*>& scenes) cam.capture(*m_photo); if (m_isroot) - printf("[Shader Progress %3.1f %%]\n", static_cast<double>(i)/scenes.size()*100.0); + printf("[Shader Progress %3.1f %%]\n", static_cast<double>(i+1)/scenes.size()*100.0); } } diff --git a/apps/polaroidCamera/SchlierenCartridge.cpp b/apps/polaroidCamera/SchlierenCartridge.cpp @@ -0,0 +1,67 @@ +// File : SchlierenCartridge.cpp +// Date : Fri 29 Apr 2016 11:30:07 AM CEST +// Author : Fabian Wermelinger +// Description: Schlieren shader implementation +// Copyright 2016 ETH Zurich. All Rights Reserved. +#include "SchlierenCartridge.h" + +void SchlierenCartridge::_compute(Slice& data) +{ + m_gradX.resize(data.width(), data.height()); + m_gradY.resize(data.width(), data.height()); + _gradX(data); + _gradY(data); + for (int h=0; h < data.height(); ++h) + for (int w=0; w < data.width(); ++w) + { + const Real IphiI2 = m_gradX(w,h)*m_gradX(w,h) + m_gradY(w,h)*m_gradY(w,h); + data(w,h) = sqrt(IphiI2); + } +} + +// 4th-order Finite Differences +void SchlierenCartridge::_gradX(const Slice& data) +{ + // this assumes + // 1.) uniform grid-spacing + // 2.) maximum extend = 1 + const Real fac = data.get_max3Ddim()/12.0; + + // left + for (int h=0; h < data.height(); ++h) + for (int w=0; w < 2; ++w) + m_gradX(w,h) = fac*(-25.0*data(w,h) + 48.0*data(w+1,h) - 36.0*data(w+2,h) + 16.0*data(w+3,h) - 3.0*data(w+4,h)); + + // right + for (int h=0; h < data.height(); ++h) + for (int w=data.width()-2; w < data.width(); ++w) + m_gradX(w,h) = fac*(-25.0*data(w,h) + 48.0*data(w-1,h) - 36.0*data(w-2,h) + 16.0*data(w-3,h) - 3.0*data(w-4,h)); + + // interior + for (int h=0; h < data.height(); ++h) + for (int w=2; w < data.width()-2; ++w) + m_gradX(w,h) = fac*(-data(w+2,h) + 8.0*data(w+1,h) - 8.0*data(w-1,h) + data(w-2,h)); +} + +void SchlierenCartridge::_gradY(const Slice& data) +{ + // this assumes + // 1.) uniform grid-spacing + // 2.) maximum extend = 1 + const Real fac = data.get_max3Ddim()/12.0; + + // left + for (int h=0; h < 2; ++h) + for (int w=0; w < data.width(); ++w) + m_gradY(w,h) = fac*(-25.0*data(w,h) + 48.0*data(w,h+1) - 36.0*data(w,h+2) + 16.0*data(w,h+3) - 3.0*data(w,h+4)); + + // right + for (int h=data.height()-2; h < data.height(); ++h) + for (int w=0; w < data.width(); ++w) + m_gradY(w,h) = fac*(-25.0*data(w,h) + 48.0*data(w,h-1) - 36.0*data(w,h-2) + 16.0*data(w,h-3) - 3.0*data(w,h-4)); + + // interior + for (int h=2; h < data.height()-2; ++h) + for (int w=0; w < data.width(); ++w) + m_gradY(w,h) = fac*(-data(w,h+2) + 8.0*data(w,h+1) - 8.0*data(w,h-1) + data(w,h-2)); +} diff --git a/apps/polaroidCamera/SchlierenCartridge.h b/apps/polaroidCamera/SchlierenCartridge.h @@ -6,39 +6,58 @@ #ifndef SCHLIERENCARTRIDGE_H_AUHX7ILM #define SCHLIERENCARTRIDGE_H_AUHX7ILM -#include "Cartridge.h" +#include "NormalizerCartridge.h" -class SchlierenCartridge : public Cartridge +class SchlierenCartridge : public NormalizerCartridge { + Slice m_gradX; + Slice m_gradY; + + // helper + void _compute(Slice& data); + void _gradX(const Slice& data); + void _gradY(const Slice& data); + public: - SchlierenCartridge(ArgumentParser& parser) : Cartridge(parser) {} + SchlierenCartridge(ArgumentParser& parser) : NormalizerCartridge(parser) {} - virtual void capture(PhotoPaper& photo, Slice& data) const + virtual void capture(PhotoPaper& photo, Slice& data) { + const Real ka = m_parser("-ka").asDouble(1.0); + const Real k0 = m_parser("-k0").asDouble(0.0); + const Real k1 = m_parser("-k1").asDouble(1.0); + photo.resize(data.width(), data.height()); // set description - string desc("2D Schlieren Data"); + string desc("2D_Schlieren"); photo.set_description(desc.c_str()); - // compute min/max - Real data_min = data(0,0); - Real data_max = data(0,0); + // compute schlieren + if (!m_bComputed) + { + _compute(data); + m_dataMin = data.min(); + m_dataMax = data.max(); + } + const Real dataMaxInv = 1.0/m_dataMax; + + // apply shader + const Real fac = -ka/(k1 - k0); for (int h=0; h < data.height(); ++h) for (int w=0; w < data.width(); ++w) - { - data_min = min(data_min, data(w,h)); - data_max = max(data_max, data(w,h)); - } - const Real data_normInv = 1.0 / (data_max - data_min); - - // compute pixel - for (int h=0; h < data.height(); ++h) - for (int w=0; w < data.width(); ++w) - photo.set_pixel(w, h, (data(w,h) - data_min) * data_normInv); + photo.set_pixel(w, h, std::exp(fac*(data(w,h)*dataMaxInv - k0))); photo.write(); } + + virtual void compute(Slice& data) + { + _compute(data); + m_dataMin = std::min(m_dataMin, data.min()); + m_dataMax = std::max(m_dataMax, data.max()); + m_bComputed = true; + } }; #endif /* SCHLIERENCARTRIDGE_H_AUHX7ILM */ diff --git a/apps/polaroidCamera/TransmissionCartridge.h b/apps/polaroidCamera/TransmissionCartridge.h @@ -18,7 +18,7 @@ public: photo.resize(data.width(), data.height()); // set description - string desc("2D Slice Data Propagation"); + string desc("2D_Transmission"); photo.set_description(desc.c_str()); // put pixel diff --git a/include/Cartridge.h b/include/Cartridge.h @@ -22,7 +22,7 @@ public: virtual void capture(PhotoPaper& photo, Slice& data) = 0; virtual void compute(Slice& data) = 0; - inline void reset() { m_bComputed = false; } + virtual void reset() { m_bComputed = false; } }; #endif /* CARTRIDGE_H_IRUVKOCT */ diff --git a/include/Types.h b/include/Types.h @@ -16,13 +16,15 @@ class Slice2D int m_width, m_height, m_N; Real* m_data; + int m_max3Ddim; + public: Slice2D() : m_width(0), m_height(0), m_N(0), m_data(nullptr) {} Slice2D(const int width, const int height) : m_width(width), m_height(height), m_N((width+_EX-_SX)*(height+_EY-_SY)) { m_data = new Real[m_N]; - for (int i = 0; i < m_N; ++i) - m_data[i] = 0.0; + // for (int i = 0; i < m_N; ++i) + // m_data[i] = 0.0; } virtual ~Slice2D() { delete [] m_data; } @@ -49,8 +51,8 @@ public: m_N = (width+_EX-_SX)*(height+_EY-_SY); if (m_data) delete [] m_data; m_data = new Real[m_N]; - for (int i = 0; i < m_N; ++i) - m_data[i] = 0.0; + // for (int i = 0; i < m_N; ++i) + // m_data[i] = 0.0; } inline int width() const { return m_width; } @@ -69,6 +71,8 @@ public: f = std::max(f, m_data[i]); return f; } + inline void set_max3Ddim(const int d) { m_max3Ddim = d; } + inline int get_max3Ddim() const { return m_max3Ddim; } }; typedef Slice2D<0,0,0,0> Slice; diff --git a/src/Polaroid.cpp b/src/Polaroid.cpp @@ -4,6 +4,7 @@ // Description: Polaroid Implementation // Copyright 2016 ETH Zurich. All Rights Reserved. #include <cmath> +#include <algorithm> #include <string> #include <sstream> #include <cstdlib> @@ -46,6 +47,9 @@ void Polaroid::load_hdf5(const char* filename, const double fraction, const int dataspace_id = H5Screate_simple(rank, dims, NULL); int status = H5Dread(dataset_id, HDF_PRECISION, dataspace_id, file_dataspace_id, H5P_DEFAULT, data); + const int Nmax = std::max(maxDim[0], std::max(maxDim[1], maxDim[2])); + m_data.set_max3Ddim(Nmax); + /* release stuff */ delete [] dims; status = H5Dclose(dataset_id); @@ -102,6 +106,8 @@ void Polaroid::load_wavelet(const char* filename, const double fraction, const i const int NBY = myreader.yblocks(); const int NBZ = myreader.zblocks(); const int maxDim[3] = {NBX*_BLOCKSIZE_, NBY*_BLOCKSIZE_, NBZ*_BLOCKSIZE_}; + const int Nmax = std::max(maxDim[0], std::max(maxDim[1], maxDim[2])); + m_data.set_max3Ddim(Nmax); Real blockdata[_BLOCKSIZE_][_BLOCKSIZE_][_BLOCKSIZE_];