-
Notifications
You must be signed in to change notification settings - Fork 161
/
gaussian_curvature.rs
executable file
·589 lines (518 loc) · 26.1 KB
/
gaussian_curvature.rs
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
/*
This tool is part of the WhiteboxTools geospatial analysis library.
Authors: Dr. John Lindsay
Created: 12/01/2022
Last Modified: 15/01/2022
License: MIT
*/
use whitebox_raster::*;
use crate::tools::*;
use num_cpus;
use std::env;
use std::f64;
use std::io::{Error, ErrorKind};
use std::path;
use std::sync::mpsc;
use std::sync::Arc;
use std::thread;
use whitebox_common::utils::{
get_formatted_elapsed_time,
haversine_distance,
vincenty_distance
};
/// This tool calculates the Gaussian curvature from a digital elevation model (DEM). Gaussian curvature
/// is the product of maximal and minimal curvatures, and retains values in each point of the topographic
/// surface after its bending without breaking, stretching, and compressing (Florinsky, 2017). Gaussian
/// curvature is measured in units of m<sup>-2</sup>.
///
/// ![](../../doc_img/GaussianCurvature.png)
///
/// The user must specify the name of the input DEM (`--dem`) and the output raster (`--output`). The
/// The Z conversion factor (`--zfactor`) is only important when the vertical and horizontal units are not the
/// same in the DEM. When this is the case, the algorithm will multiply each elevation in the DEM by the
/// Z Conversion Factor. Curvature values are often very small and as such the user may opt to log-transform
/// the output raster (`--log`). Transforming the values applies the equation by Shary et al. (2002):
///
/// *Θ*' = sign(*Θ*) ln(1 + 10<sup>*n*</sup>|*Θ*|)
///
/// where *Θ* is the parameter value and *n* is dependent on the grid cell size.
///
/// For DEMs in projected coordinate systems, the tool uses the 3rd-order bivariate
/// Taylor polynomial method described by Florinsky (2016). Based on a polynomial fit
/// of the elevations within the 5x5 neighbourhood surrounding each cell, this method is considered more
/// robust against outlier elevations (noise) than other methods. For DEMs in geographic coordinate systems
/// (i.e. angular units), the tool uses the 3x3 polynomial fitting method for equal angle grids also
/// described by Florinsky (2016).
///
/// # References
/// Florinsky, I. (2016). Digital terrain analysis in soil science and geology. Academic Press.
///
/// Florinsky, I. V. (2017). An illustrated introduction to general geomorphometry. Progress in Physical
/// Geography, 41(6), 723-752.
///
/// Shary P. A., Sharaya L. S. and Mitusov A. V. (2002) Fundamental quantitative methods of land surface analysis.
/// Geoderma 107: 1–32.
///
/// `TangentialCurvature`, `ProfileCurvature`, `PlanCurvature`, `MeanCurvature`, `MinimalCurvature`, `MaximalCurvature`
pub struct GaussianCurvature {
name: String,
description: String,
toolbox: String,
parameters: Vec<ToolParameter>,
example_usage: String,
}
impl GaussianCurvature {
pub fn new() -> GaussianCurvature {
// public constructor
let name = "GaussianCurvature".to_string();
let toolbox = "Geomorphometric Analysis".to_string();
let description = "Calculates a mean curvature raster from an input DEM.".to_string();
let mut parameters = vec![];
parameters.push(ToolParameter {
name: "Input DEM File".to_owned(),
flags: vec!["-i".to_owned(), "--dem".to_owned()],
description: "Input raster DEM file.".to_owned(),
parameter_type: ParameterType::ExistingFile(ParameterFileType::Raster),
default_value: None,
optional: false,
});
parameters.push(ToolParameter {
name: "Output File".to_owned(),
flags: vec!["-o".to_owned(), "--output".to_owned()],
description: "Output raster file.".to_owned(),
parameter_type: ParameterType::NewFile(ParameterFileType::Raster),
default_value: None,
optional: false,
});
parameters.push(ToolParameter {
name: "Log-transform the output?".to_owned(),
flags: vec!["--log".to_owned()],
description: "Display output values using a log-scale."
.to_owned(),
parameter_type: ParameterType::Boolean,
default_value: Some("false".to_string()),
optional: true,
});
parameters.push(ToolParameter {
name: "Z Conversion Factor".to_owned(),
flags: vec!["--zfactor".to_owned()],
description:
"Optional multiplier for when the vertical and horizontal units are not the same."
.to_owned(),
parameter_type: ParameterType::Float,
default_value: None,
optional: true,
});
let sep: String = path::MAIN_SEPARATOR.to_string();
let e = format!("{}", env::current_exe().unwrap().display());
let mut parent = env::current_exe().unwrap();
parent.pop();
let p = format!("{}", parent.display());
let mut short_exe = e
.replace(&p, "")
.replace(".exe", "")
.replace(".", "")
.replace(&sep, "");
if e.contains(".exe") {
short_exe += ".exe";
}
let usage = format!(
">>.*{} -r={} -v --wd=\"*path*to*data*\" --dem=DEM.tif -o=output.tif",
short_exe, name
)
.replace("*", &sep);
GaussianCurvature {
name: name,
description: description,
toolbox: toolbox,
parameters: parameters,
example_usage: usage,
}
}
}
impl WhiteboxTool for GaussianCurvature {
fn get_source_file(&self) -> String {
String::from(file!())
}
fn get_tool_name(&self) -> String {
self.name.clone()
}
fn get_tool_description(&self) -> String {
self.description.clone()
}
fn get_tool_parameters(&self) -> String {
let mut s = String::from("{\"parameters\": [");
for i in 0..self.parameters.len() {
if i < self.parameters.len() - 1 {
s.push_str(&(self.parameters[i].to_string()));
s.push_str(",");
} else {
s.push_str(&(self.parameters[i].to_string()));
}
}
s.push_str("]}");
s
}
fn get_example_usage(&self) -> String {
self.example_usage.clone()
}
fn get_toolbox(&self) -> String {
self.toolbox.clone()
}
fn run<'a>(
&self,
args: Vec<String>,
working_directory: &'a str,
verbose: bool,
) -> Result<(), Error> {
let tool_name = self.get_tool_name();
let sep: String = path::MAIN_SEPARATOR.to_string();
// Read in the environment variables and get the necessary values
let configs = whitebox_common::configs::get_configs()?;
let max_procs = configs.max_procs;
// read the arguments
let mut input_file: String = String::new();
let mut output_file: String = String::new();
let mut log_transform = false;
let mut z_factor = 1f64;
if args.len() <= 1 {
return Err(Error::new(
ErrorKind::InvalidInput,
"Tool run with too few parameters.",
));
}
for i in 0..args.len() {
let mut arg = args[i].replace("\"", "");
arg = arg.replace("\'", "");
let cmd = arg.split("="); // in case an equals sign was used
let vec = cmd.collect::<Vec<&str>>();
let mut keyval = false;
if vec.len() > 1 {
keyval = true;
}
let flag_val = vec[0].to_lowercase().replace("--", "-");
if flag_val == "-i" || flag_val == "-input" || flag_val == "-dem" {
input_file = if keyval {
vec[1].to_string()
} else {
args[i + 1].to_string()
};
} else if flag_val == "-o" || flag_val == "-output" {
output_file = if keyval {
vec[1].to_string()
} else {
args[i + 1].to_string()
};
} else if flag_val == "-zfactor" {
z_factor = if keyval {
vec[1]
.to_string()
.parse::<f64>()
.expect(&format!("Error parsing {}", flag_val))
} else {
args[i + 1]
.to_string()
.parse::<f64>()
.expect(&format!("Error parsing {}", flag_val))
};
} else if flag_val == "-log" {
if vec.len() == 1 || !vec[1].to_string().to_lowercase().contains("false") {
log_transform = true;
}
}
}
if verbose {
let welcome_len = format!("* Welcome to {} *", tool_name).len().max(28);
// 28 = length of the 'Powered by' by statement.
println!("{}", "*".repeat(welcome_len));
println!("* Welcome to {} {}*", tool_name, " ".repeat(welcome_len - 15 - tool_name.len()));
println!("* Powered by WhiteboxTools {}*", " ".repeat(welcome_len - 28));
println!("* www.whiteboxgeo.com {}*", " ".repeat(welcome_len - 23));
println!("{}", "*".repeat(welcome_len));
}
let mut progress: usize;
let mut old_progress: usize = 1;
let start = Instant::now();
if !input_file.contains(&sep) && !input_file.contains("/") {
input_file = format!("{}{}", working_directory, input_file);
}
if !output_file.contains(&sep) && !output_file.contains("/") {
output_file = format!("{}{}", working_directory, output_file);
}
// Read in the input raster
let input = Arc::new(Raster::new(&input_file, "r")?);
let rows = input.configs.rows as isize;
let columns = input.configs.columns as isize;
let nodata = input.configs.nodata;
let resx = input.configs.resolution_x;
let resy = input.configs.resolution_y;
let res = (resx + resy) / 2.;
let mut num_procs = num_cpus::get() as isize;
if max_procs > 0 && max_procs < num_procs {
num_procs = max_procs;
}
let (tx, rx) = mpsc::channel();
if !input.is_in_geographic_coordinates() {
// Based on Florinsky (2016) pg. 246
let log_multiplier = match res {
x if x >= 0. && x < 1. => { 10f64.powi(2) },
x if x >= 1. && x < 10. => { 10f64.powi(3) },
x if x >= 10. && x < 100. => { 10f64.powi(4) },
x if x >= 100. && x < 1000. => { 10f64.powi(5) },
x if x >= 1000. && x < 5000. => { 10f64.powi(6) },
x if x >= 5000. && x < 10000. => { 10f64.powi(7) },
x if x >= 10000. && x < 75000. => { 10f64.powi(8) },
_ => { 10f64.powi(9) },
};
for tid in 0..num_procs {
let input = input.clone();
let tx = tx.clone();
thread::spawn(move || {
let mut z12: f64;
let mut p: f64;
let mut q: f64;
let mut r: f64;
let mut s: f64;
let mut t: f64;
let mut gaussian_curv: f64;
let offsets = [
[-2, -2], [-1, -2], [0, -2], [1, -2], [2, -2],
[-2, -1], [-1, -1], [0, -1], [1, -1], [2, -1],
[-2, 0], [-1, 0], [0, 0], [1, 0], [2, 0],
[-2, 1], [-1, 1], [0, 1], [1, 1], [2, 1],
[-2, 2], [-1, 2], [0, 2], [1, 2], [2, 2]
];
let mut z = [0f64; 25];
for row in (0..rows).filter(|r| r % num_procs == tid) {
let mut data = vec![nodata; columns as usize];
for col in 0..columns {
z12 = input.get_value(row, col);
if z12 != nodata {
for n in 0..25 {
z[n] = input.get_value(row + offsets[n][1], col + offsets[n][0]);
if z[n] != nodata {
z[n] *= z_factor;
} else {
z[n] = z12 * z_factor;
}
}
/*
The following equations have been taken from Florinsky (2016) Principles and Methods
of Digital Terrain Modelling, Chapter 4, pg. 117. Note that I believe Florinsky reversed
the equations for q and p.
*/
r = 1f64 / (35f64 * res * res) * (2. * (z[0] + z[4] + z[5] + z[9] + z[10] + z[14] + z[15] + z[19] + z[20] + z[24])
- 2. * (z[2] + z[7] + z[12] + z[17] + z[22]) - z[1] - z[3] - z[6] - z[8]
- z[11] - z[13] - z[16] - z[18] - z[21] - z[23]);
t = 1f64 / (35f64 * res * res) * (2. * (z[0] + z[1] + z[2] + z[3] + z[4] + z[20] + z[21] + z[22] + z[23] + z[24])
- 2. * (z[10] + z[11] + z[12] + z[13] + z[14]) - z[5] - z[6] - z[7] - z[8]
- z[9] - z[15] - z[16] - z[17] - z[18] - z[19]);
s = 1. / (100. * res * res) * (z[8] + z[16] - z[6] - z[18] + 4. * (z[4] + z[20] - z[0] - z[24])
+ 2. * (z[3] + z[9] + z[15] + z[21] - z[1] - z[5] - z[19] - z[23]));
p = 1. / (420. * res) * (44. * (z[3] + z[23] - z[1] - z[21]) + 31. * (z[0] + z[20] - z[4] - z[24]
+ 2. * (z[8] + z[18] - z[6] - z[16])) + 17. * (z[14] - z[10] + 4. * (z[13] - z[11]))
+ 5. * (z[9] + z[19] - z[5] - z[15]));
q = 1. / (420. * res) * (44. * (z[5] + z[9] - z[15] - z[19]) + 31. * (z[20] + z[24] - z[0] - z[4]
+ 2. * (z[6] + z[8] - z[16] - z[18])) + 17. * (z[2] - z[22] + 4. * (z[7] - z[17]))
+ 5. * (z[1] + z[3] - z[21] - z[23]));
/*
The following equation has been taken from Florinsky (2016) Principles and Methods
of Digital Terrain Modelling, Chapter 2, pg. 18.
*/
gaussian_curv = (r * t - s * s) / (1. + p * p + q * q).powi(2);
if log_transform {
// Based on Florinsky (2016) pg. 244 eq. 8.1
gaussian_curv = gaussian_curv.signum() * (1. + log_multiplier * gaussian_curv.abs()).ln();
}
data[col as usize] = gaussian_curv;
}
}
tx.send((row, data)).unwrap();
}
});
}
} else { // geographic coordinates
let phi1 = input.get_y_from_row(0);
let lambda1 = input.get_x_from_column(0);
let phi2 = phi1;
let lambda2 = input.get_x_from_column(-1);
let linear_res = vincenty_distance((phi1, lambda1), (phi2, lambda2));
let lr2 = haversine_distance((phi1, lambda1), (phi2, lambda2));
let diff = 100. * (linear_res - lr2).abs() / linear_res;
let use_haversine = diff < 0.5; // if the difference is less than 0.5%, use the faster haversine method to calculate distances.
// Based on Florinsky (2016) pg. 246
let log_multiplier = match linear_res {
x if x >= 0. && x < 1. => { 10f64.powi(2) },
x if x >= 1. && x < 10. => { 10f64.powi(3) },
x if x >= 10. && x < 100. => { 10f64.powi(4) },
x if x >= 100. && x < 1000. => { 10f64.powi(5) },
x if x >= 1000. && x < 5000. => { 10f64.powi(6) },
x if x >= 5000. && x < 10000. => { 10f64.powi(7) },
x if x >= 10000. && x < 75000. => { 10f64.powi(8) },
_ => { 10f64.powi(9) },
};
for tid in 0..num_procs {
let input = input.clone();
let tx = tx.clone();
thread::spawn(move || {
let mut z4: f64;
let mut p: f64;
let mut q: f64;
let mut r: f64;
let mut s: f64;
let mut t: f64;
let mut a: f64;
let mut b: f64;
let mut c: f64;
let mut d: f64;
let mut e: f64;
let mut phi1: f64;
let mut lambda1: f64;
let mut phi2: f64;
let mut lambda2: f64;
let mut gaussian_curv: f64;
let offsets = [
[-1, -1], [0, -1], [1, -1],
[-1, 0], [0, 0], [1, 0],
[-1, 1], [0, 1], [1, 1]
];
let mut z = [0f64; 25];
for row in (0..rows).filter(|r| r % num_procs == tid) {
let mut data = vec![nodata; columns as usize];
for col in 0..columns {
z4 = input.get_value(row, col);
if z4 != nodata {
for n in 0..9 {
z[n] = input.get_value(row + offsets[n][1], col + offsets[n][0]);
if z[n] != nodata {
z[n] *= z_factor;
} else {
z[n] = z4 * z_factor;
}
}
// Calculate a, b, c, d, and e.
phi1 = input.get_y_from_row(row);
lambda1 = input.get_x_from_column(col);
phi2 = phi1;
lambda2 = input.get_x_from_column(col-1);
b = if use_haversine {
haversine_distance((phi1, lambda1), (phi2, lambda2))
} else {
vincenty_distance((phi1, lambda1), (phi2, lambda2))
};
phi2 = input.get_y_from_row(row+1);
lambda2 = lambda1;
d = if use_haversine {
haversine_distance((phi1, lambda1), (phi2, lambda2))
} else {
vincenty_distance((phi1, lambda1), (phi2, lambda2))
};
phi2 = input.get_y_from_row(row-1);
lambda2 = lambda1;
e = if use_haversine {
haversine_distance((phi1, lambda1), (phi2, lambda2))
} else {
vincenty_distance((phi1, lambda1), (phi2, lambda2))
};
phi1 = input.get_y_from_row(row+1);
lambda1 = input.get_x_from_column(col);
phi2 = phi1;
lambda2 = input.get_x_from_column(col-1);
a = if use_haversine {
haversine_distance((phi1, lambda1), (phi2, lambda2))
} else {
vincenty_distance((phi1, lambda1), (phi2, lambda2))
};
phi1 = input.get_y_from_row(row-1);
lambda1 = input.get_x_from_column(col);
phi2 = phi1;
lambda2 = input.get_x_from_column(col-1);
c = if use_haversine {
haversine_distance((phi1, lambda1), (phi2, lambda2))
} else {
vincenty_distance((phi1, lambda1), (phi2, lambda2))
};
/*
The following equations have been taken from Florinsky (2016) Principles and Methods
of Digital Terrain Modelling, Chapter 4, pg. 117.
*/
r = (c * c * (z[0] + z[2] - 2. * z[1]) + b * b * (z[3] + z[5] - 2. * z[4]) + a * a * (z[6] + z[8] - 2. * z[7]))
/ (a.powi(4) + b.powi(4) + c.powi(4));
t = 2. / (3. * d * e * (d + e) * (a.powi(4) + b.powi(4) + c.powi(4)))
* ((d * (a.powi(4) + b.powi(4) + b * b * c * c) - c * c * e * (a * a - b * b)) * (z[0] + z[2])
- (d * (a.powi(4) + c.powi(4) + b * b * c * c) + e * (a.powi(4) + c.powi(4) + a * a * b * b)) * (z[3] + z[5])
+ (e * (b.powi(4) + c.powi(4) + a * a * b * b) + a * a * d * (b * b - c * c)) * (z[6] + z[8])
+ d * (b.powi(4) * (z[1] - 3. * z[4]) + c.powi(4) * (3. * z[1] - z[4]) + (a.powi(4) - 2. * b * b * c * c) * (z[1] - z[4]))
+ e * (a.powi(4) * (3. * z[7] - z[4]) + b.powi(4) * (z[7] - 3. * z[4]) + (c.powi(4) - 2. * a * a * b * b) * (z[7] - z[4]))
- 2. * (a * a * d * (b * b - c * c) * z[7] - c * c * e * (a * a - b * b) * z[1]));
s = (c * (a * a * (d + e) + b * b * e) * (z[2] - z[0]) - b * (a * a * d - c * c * e) * (z[3] - z[5]) + a * (c * c * (d + e) + b * b * d) * (z[6] - z[8]))
/ (2. * (a * a * c * c * (d + e).powi(2) + b * b * (a * a * d * d + c * c * e * e)));
p = (a * a * c * d * (d + e) * (z[2] - z[0]) + b * (a * a * d * d + c * c * e * e) * (z[5] - z[3]) + a * c * c * e * (d + e) * (z[8] - z[6]))
/ (2. * (a * a * c * c * (d + e).powi(2) + b * b * (a * a * d * d + c * c * e * e)));
q = 1. / (3. * d * e * (d + e) * (a.powi(4) + b.powi(4) + c.powi(4)))
* ((d * d * (a.powi(4) + b.powi(4) + b * b * c * c) + c * c * e * e * (a * a - b * b)) * (z[0] + z[2])
- (d * d * (a.powi(4) + c.powi(4) + b * b * c * c) - e * e * (a.powi(4) + c.powi(4) + a * a * b * b)) * (z[3] + z[5])
- (e * e * (b.powi(4) + c.powi(4) + a * a * b * b) - a * a * d * d * (b * b - c * c)) * (z[6] + z[8])
+ d * d * (b.powi(4) * (z[1] - 3. * z[4]) + c.powi(4) * (3. * z[1] - z[4]) + (a.powi(4) - 2. * b * b * c * c) * (z[1] - z[4]))
+ e * e * (a.powi(4) * (z[4] - 3. * z[7]) + b.powi(4) * (3. * z[4] - z[7]) + (c.powi(4) - 2. * a * a * b * b) * (z[4] - z[7]))
- 2. * (a * a * d * d * (b * b - c * c) * z[7] + c * c * e * e * (a * a - b * b) * z[1]));
/*
The following equation has been taken from Florinsky (2016) Principles and Methods
of Digital Terrain Modelling, Chapter 2, pg. 18.
*/
gaussian_curv = (r * t - s * s) / (1. + p * p + q * q).powi(2);
if log_transform {
// Based on Florinsky (2016) pg. 244 eq. 8.1
gaussian_curv = gaussian_curv.signum() * (1. + log_multiplier * gaussian_curv.abs()).ln();
}
data[col as usize] = gaussian_curv;
}
}
tx.send((row, data)).unwrap();
}
});
}
}
let mut output = Raster::initialize_using_file(&output_file, &input);
output.configs.data_type = DataType::F32;
for row in 0..rows {
let (r, data) = rx.recv().expect("Error receiving data from thread.");
output.set_row_data(r, data);
if verbose {
progress = (100.0_f64 * row as f64 / (rows - 1) as f64) as usize;
if progress != old_progress {
println!("Progress: {}%", progress);
old_progress = progress;
}
}
}
//////////////////////
// Output the image //
//////////////////////
if verbose {
println!("Saving data...")
};
let elapsed_time = get_formatted_elapsed_time(start);
output.add_metadata_entry(format!(
"Created by whitebox_tools\' {} tool",
tool_name
));
output.add_metadata_entry(format!("Input file: {}", input_file));
output.add_metadata_entry(format!("Elapsed Time (excluding I/O): {}", elapsed_time));
let _ = match output.write() {
Ok(_) => {
if verbose {
println!("Output file written")
}
}
Err(e) => return Err(e),
};
if verbose {
println!(
"{}",
&format!("Elapsed Time (excluding I/O): {}", elapsed_time)
);
}
Ok(())
}
}