37#include <boost/format.hpp>
43#include "DGtal/base/Common.h"
44#include "DGtal/helpers/StdDefs.h"
45#include "DGtal/io/readers/GenericReader.h"
46#include "DGtal/io/writers/GenericWriter.h"
200using namespace DGtal;
202int main(
int argc,
char* argv[] )
210 string inpainting_mask;
224 string scv {
"0xff0000"};
229 stringstream ssDescr;
230 ssDescr <<
"Computes a piecewise smooth approximation of a grey-level or color image, by optimizing the Ambrosio-Tortorelli functional (with u a 2-form and v a 0-form).";
231 ssDescr <<
"Usage: " << argv[0] <<
" -i toto.pgm\n"
232 <<
"Computes the Ambrosio-Tortorelli reconstruction/segmentation of an input image."
233 <<
"It outputs 2 or 3 images (of basename given by option --output) giving the"
234 <<
" reconstructed image u, and other images superposing u and the discontinuities v."
238 <<
" | a.(u-g)^2 + v^2 |grad u|^2 + le.|grad v|^2 + (l/4e).(1-v)^2 "
242 <<
"Discretized as (u 2-form, v 0-form, A vertex-edge bdry, B edge-face bdy, M vertex-edge average)" << endl
243 <<
"E(u,v) = a(u-g)^t (u-g) + u^t B diag(M v)^2 B^t u + l e v^t A^t A v + l/(4e) (1-v)^t (1-v)" << endl
245 <<
"Example: ./at-u2-v0 -i ../Images/cerclesTriangle64b02.pgm -o tmp -a 0.05 -e 1 --lambda-1 0.1 --lambda-2 0.00001 -g"
248 app.description(ssDescr.str());
251 app.add_option(
"-i,--input,1", f1,
"the input image PPM filename." )
253 ->check(CLI::ExistingFile);
254 app.add_option(
"--inpainting-mask,-m", inpainting_mask,
"the input inpainting mask filename." );
255 app.add_option(
"--output,-o", f2,
"the output image basename.");
256 auto lambdaOpt = app.add_option(
"--lambda,-l",l,
"the parameter lambda.");
257 app.add_flag(
"--metric-average,-M",metric,
"use metric average to smooth L1-metric." );
259 app.add_option(
"--lambda-1,-1",l1,
"the initial parameter lambda (l1).");
260 app.add_option(
"--lambda-2,-2",l2,
"the final parameter lambda (l2)." );
261 app.add_option(
"--lambda-ratio,-q",lr,
"the division ratio for lambda from l1 to l2.");
262 app.add_option(
"--alpha,-a",a,
"the parameter alpha.");
263 auto epsOpt = app.add_option(
"--epsilon,-e",
"the initial and final parameter epsilon of AT functional at the same time.");
265 app.add_option(
"--epsilon-1",e1,
"the initial parameter epsilon.");
266 app.add_option(
"--epsilon-2",e2,
"the final parameter epsilon.");
267 app.add_option(
"--epsilon-r",er,
"sets the ratio between two consecutive epsilon values of AT functional.");
269 app.add_option(
"--nbiter,-n",nbiter,
"the maximum number of iterations." );
270 auto snrOpt = app.add_option(
"--image-snr", isnr,
"the input image without deterioration if you wish to compute the SNR.");
271 app.add_option(
"--pixel-size,-p", pix_sz,
"the pixel size for outputing images (useful when one wants to see the discontinuities v on top of u).");
272 app.add_option(
"--color-v,-c",scv,
"the color chosen for displaying the singularities v (e.g. red is 0xff0000)." );
273 app.add_option(
"--verbose,-v", verb,
"the verbose level (0: silent, 1: less silent, etc)." );
275 app.get_formatter()->column_width(40);
276 CLI11_PARSE(app, argc, argv);
281 Color color_v( (
unsigned int) std::stoul( scv,
nullptr, 16 ), 255 );
282 if ( lambdaOpt->count()) l1 = l2 = l;
283 if ( l2 > l1 ) l2 = l1;
284 if ( lr <= 1.0 ) lr = sqrt(2);
285 if ( epsOpt->count() > 0 ){
289 bool snr = snrOpt->count() > 0;
292 bool color_image = f1.size() > 4 && f1.compare( f1.size() - 4, 4,
".ppm" ) == 0;
293 bool grey_image = f1.size() > 4 && f1.compare( f1.size() - 4, 4,
".pgm" ) == 0;
294 if ( ! color_image && ! grey_image )
296 trace.error() <<
"Input image file must be either a PGM (grey-level) or a PPM (color) image with these extensions."
304 AT.setMetricAverage( metric );
307 typedef ImageContainerBySTLVector<Domain, Color> ColorImage;
308 typedef ImageContainerBySTLVector<Domain, unsigned char> GreyLevelImage;
312 trace.beginBlock(
"Reading PPM image");
313 ColorImage image = PPMReader<ColorImage>::importPPM( f1 );
315 trace.beginBlock(
"Building AT");
316 domain = image.domain();
317 K.init( domain.lowerBound(), domain.upperBound(),
true );
319 AT.addInput( image, [] ( Color c ) ->
double {
return ((
double) c.red()) / 255.0; } );
320 AT.addInput( image, [] ( Color c ) ->
double {
return ((
double) c.green()) / 255.0; } );
321 AT.addInput( image, [] ( Color c ) ->
double {
return ((
double) c.blue()) / 255.0; } );
324 else if ( grey_image )
326 trace.beginBlock(
"Reading PGM image");
327 GreyLevelImage image = PGMReader<GreyLevelImage>::importPGM( f1 );
329 trace.beginBlock(
"Building AT");
330 domain = image.domain();
331 K.init( domain.lowerBound(), domain.upperBound(),
true );
333 AT.addInput( image, [] (
unsigned char c ) {
return ((
double) c) / 255.0; } );
338 if ( snr && color_image )
340 trace.beginBlock(
"Reading ideal PPM image");
341 ColorImage image = PPMReader<ColorImage>::importPPM( isnr );
343 AT.addInput( image, [] ( Color c ) ->
double {
return ((
double) c.red()) / 255.0; }, true );
344 AT.addInput( image, [] ( Color c ) ->
double {
return ((
double) c.green()) / 255.0; }, true );
345 AT.addInput( image, [] ( Color c ) ->
double {
return ((
double) c.blue()) / 255.0; }, true );
347 else if ( snr && grey_image )
349 trace.beginBlock(
"Reading ideal PGM image");
350 GreyLevelImage image = PGMReader<GreyLevelImage>::importPGM( isnr );
352 AT.addInput( image, [] (
unsigned char c ) {
return ((
double) c) / 255.0; }, true );
357 Domain out_domain( pix_sz * domain.lowerBound(),
358 pix_sz * domain.upperBound() + Point::diagonal( pix_sz ) );
361 double g_snr = snr ? AT.computeSNR() : 0.0;
363 if ( inpainting_mask.size() > 0 )
365 string fm = inpainting_mask;
366 trace.beginBlock(
"Reading inpainting mask");
367 GreyLevelImage mask = GenericReader<GreyLevelImage>::import( fm );
369 Calculus::PrimalForm2 m( AT.calculus );
370 for ( Calculus::Index index = 0; index < m.myContainer.rows(); index++)
372 auto cell = m.getSCell( index );
373 double col = ((double) mask( K.sCoords( cell ) )) / 255.0;
374 m.myContainer( index ) = col > 0.0 ? 1.0 : 0.0;
377 AT.setUFromInputAndMask();
381 ossGM << boost::format(
"%s-g-mask.pgm") %f2;
382 GreyLevelImage image_mg( domain );
385 ( AT.calculus, mg, image_mg, 0.0, 1.0, 1 );
386 PGMWriter<GreyLevelImage>::exportPGM( ossGM.str(), image_mg );
388 else if ( color_image )
391 ossGM << boost::format(
"%s-g-mask.ppm") %f2;
392 ColorImage image_mg( domain );
397 ( AT.calculus, mg0, mg1, mg2, image_mg, 0.0, 1.0, 1 );
398 PPMWriter<ColorImage, functors::Identity >::exportPPM( ossGM.str(), image_mg );
404 trace.info() << AT << std::endl;
409 trace.info() <<
"************ lambda = " << l1 <<
" **************" << endl;
411 for ( eps = e1; eps >= e2; eps /= er )
413 trace.info() <<
" ======= epsilon = " << eps <<
" ========" << endl;
414 AT.setEpsilon( eps );
417 trace.progressBar( n, nbiter );
421 n_v = AT.computeVariation();
422 }
while ( ( n_v > 0.0001 ) && ( ++n < nbiter ) );
423 trace.progressBar( n, nbiter );
424 trace.info() <<
"[#### last variation = " << n_v <<
" " << endl;
428 if ( verb > 0 ) trace.beginBlock(
"Writing u[0] as PGM image");
429 ostringstream ossU, ossV, ossW;
430 ossU << boost::format(
"%s-a%.5f-l%.7f-u.pgm") % f2 % a % l1;
431 ossV << boost::format(
"%s-a%.5f-l%.7f-u-v.pgm") % f2 % a % l1;
432 ossW << boost::format(
"%s-a%.5f-l%.7f-u-v.ppm") % f2 % a % l1;
433 const Calculus::PrimalForm2 u = AT.getU( 0 );
434 const Calculus::PrimalForm1 v = AT.M01 * AT.getV();
436 GreyLevelImage image_u( domain );
438 ( AT.calculus, u, image_u, 0.0, 1.0, 1 );
439 PGMWriter<GreyLevelImage>::exportPGM( ossU.str(), image_u );
441 GreyLevelImage image_uv( out_domain );
443 ( AT.calculus, u, image_uv, 0.0, 1.0, pix_sz );
445 ( AT.calculus, v, image_uv, 0.0, 1.0, pix_sz );
446 PGMWriter<GreyLevelImage>::exportPGM( ossV.str(), image_uv );
448 ColorImage cimage( out_domain );
450 ( AT.calculus, u, u, u, cimage, 0.0, 1.0, pix_sz );
452 ( AT.calculus, v, cimage, color_v, 0.0, 1.0, pix_sz );
453 PPMWriter<ColorImage, functors::Identity >::exportPPM( ossW.str(), cimage );
454 if ( verb > 0 ) trace.endBlock();
456 else if ( color_image )
458 if ( verb > 0 ) trace.beginBlock(
"Writing u[0,1,2] as PGM image");
459 ostringstream ossU, ossV;
460 ossU << boost::format(
"%s-a%.5f-l%.7f-u.ppm") % f2 % a % l1;
461 ossV << boost::format(
"%s-a%.5f-l%.7f-u-v.ppm") % f2 % a % l1;
462 const Calculus::PrimalForm2 u0 = AT.getU( 0 );
463 const Calculus::PrimalForm2 u1 = AT.getU( 1 );
464 const Calculus::PrimalForm2 u2 = AT.getU( 2 );
465 const Calculus::PrimalForm1 v = AT.M01 * AT.getV();
467 ColorImage image_u( domain );
469 ( AT.calculus, u0, u1, u2, image_u, 0.0, 1.0, 1 );
470 PPMWriter<ColorImage, functors::Identity >::exportPPM( ossU.str(), image_u );
471 ColorImage image_uv( out_domain );
473 ( AT.calculus, u0, u1, u2, image_uv, 0.0, 1.0, pix_sz );
475 ( AT.calculus, v, image_uv, color_v, 0.0, 1.0, pix_sz );
476 PPMWriter<ColorImage, functors::Identity >::exportPPM( ossV.str(), image_uv );
477 if ( verb > 0 ) trace.endBlock();
482 double u_snr = AT.computeSNR();
483 trace.info() <<
"- SNR of u = " << u_snr <<
" SNR of g = " << g_snr << endl;
void primalForm1ToGreyLevelImage(const Calculus &calculus, const typename Calculus::PrimalForm1 &v, Image &image, double cut_low=0.0, double cut_up=1.0, int pixel_size=1)
void form2ToGreyLevelImage(const Calculus &calculus, const AnyForm2 &u, Image &image, double cut_low=0.0, double cut_up=1.0, int pixel_size=1)
void threeForms2ToRGBColorImage(const Calculus &calculus, const AnyForm2 &u0, const AnyForm2 &u1, const AnyForm2 &u2, Image &image, double cut_low=0.0, double cut_up=1.0, int pixel_size=1)
void primalForm1ToRGBColorImage(const Calculus &calculus, const typename Calculus::PrimalForm1 &v, Image &image, Color color, double cut_low=0.0, double cut_up=1.0, int pixel_size=1)
DGtal::LinearOperator< Calculus, dim, duality, dim, duality > diagonal(const DGtal::KForm< Calculus, dim, duality > &kform)
Aim: This class solves Ambrosio-Tortorelli functional in a plane for u a (vector of) 2-form(s) and v ...
DiscreteExteriorCalculus< 2, 2, LinearAlgebra > Calculus