fill.cpp 21 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704
  1. /****************************************************************************
  2. *
  3. * MODULE: r.terraflow
  4. *
  5. * COPYRIGHT (C) 2007 Laura Toma
  6. *
  7. * This program is free software; you can redistribute it and/or modify
  8. * it under the terms of the GNU General Public License as published by
  9. * the Free Software Foundation; either version 2 of the License, or
  10. * (at your option) any later version.
  11. *
  12. * This program is distributed in the hope that it will be useful,
  13. * but WITHOUT ANY WARRANTY; without even the implied warranty of
  14. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  15. * GNU General Public License for more details.
  16. *
  17. *****************************************************************************/
  18. #include <stdlib.h>
  19. #include <ctype.h>
  20. #include <time.h>
  21. #include <string>
  22. #include "fill.h"
  23. #include "common.h"
  24. #include "water.h"
  25. #include "sortutils.h"
  26. #include "streamutils.h"
  27. #include "filldepr.h"
  28. #include "grid.h"
  29. /* globals in common.H
  30. extern statsRecorder *stats; stats file
  31. extern userOptions *opt; command-line options
  32. extern struct Cell_head *region; header of the region
  33. extern dimension_type nrows, ncols;
  34. */
  35. #define FILL_DEBUG if(0)
  36. #define FILL_SAVEALL if(0)
  37. /* defined in this module */
  38. void recordWatersheds(AMI_STREAM<labelElevType> *labeledWater);
  39. long assignDirections(AMI_STREAM<plateauStats> *statstr,
  40. AMI_STREAM<plateauType> *platstr,
  41. AMI_STREAM<waterType> *waterstr);
  42. void assignFinalDirections(AMI_STREAM<plateauStats> *statstr,
  43. AMI_STREAM<plateauType> *platstr,
  44. AMI_STREAM<waterType> *waterstr);
  45. template<class T1, class T2, class T3, class T4, class FUN>
  46. void mergeStreamGridGrid(AMI_STREAM<T1> *grid1,
  47. AMI_STREAM<T2> *grid2,
  48. dimension_type rows, dimension_type cols,
  49. AMI_STREAM<T3> *str,
  50. FUN fo,
  51. AMI_STREAM<T4> *outStream);
  52. void merge2waterBase(AMI_STREAM<waterType> *unsortedWaterStr,
  53. AMI_STREAM<direction_type> *dirStr,
  54. AMI_STREAM<elevation_type> *elStr,
  55. AMI_STREAM<waterWindowBaseType> *merge);
  56. AMI_STREAM<waterGridType>*
  57. merge2waterGrid(AMI_STREAM<waterType> *unsortedWaterStr,
  58. AMI_STREAM<direction_type> *dirStr,
  59. AMI_STREAM<elevation_type> *elStr);
  60. AMI_STREAM<boundaryType> *
  61. findBoundariesMain(AMI_STREAM<labelElevType> *labeledWater);
  62. /* ********************************************************************** */
  63. /* some helper classes */
  64. /* ********************************************************************** */
  65. class printElevation {
  66. public:
  67. char *operator()(const elevation_type &p) {
  68. static char buf[20];
  69. sprintf(buf, "%.1f", (float)p);
  70. return buf;
  71. }
  72. };
  73. class printDirection {
  74. public:
  75. char *operator()(const direction_type &p) {
  76. static char buf[20];
  77. sprintf(buf, "%3d", p);
  78. return buf;
  79. }
  80. char *operator()(const waterWindowBaseType &p) {
  81. static char buf[20];
  82. sprintf(buf, "%3d", p.dir);
  83. return buf;
  84. }
  85. #if(0)
  86. char *operator()(const waterWindowBaseType &p) {
  87. static char buf[3];
  88. buf[0] = directionSymbol(p.dir);
  89. buf[1] = '\0';
  90. return buf;
  91. }
  92. #endif
  93. };
  94. class printLabel {
  95. public:
  96. char *operator()(const labelElevType&p) {
  97. static char buf[8];
  98. sprintf(buf, CCLABEL_FMT, p.getLabel());
  99. return buf;
  100. }
  101. char *operator()(const waterGridType &p) {
  102. static char buf[8];
  103. sprintf(buf, CCLABEL_FMT, p.getLabel());
  104. return buf;
  105. }
  106. char *operator()(const waterType &p) {
  107. static char buf[8];
  108. sprintf(buf, CCLABEL_FMT, p.getLabel());
  109. return buf;
  110. }
  111. };
  112. class printDepth {
  113. public:
  114. char *operator()(const waterGridType &p) {
  115. static char buf[3];
  116. sprintf(buf, "%1u", p.depth);
  117. return buf;
  118. }
  119. };
  120. char *
  121. verbosedir(std::string s) {
  122. static char buf[BUFSIZ];
  123. sprintf(buf, "dump/%s", s.c_str());
  124. return buf;
  125. }
  126. /* ---------------------------------------------------------------------- */
  127. /* fill the terrain (if necessary) and compute flow direction stream;
  128. elstr must exist and contain grid data before call, filledstr and
  129. dirstr are created; elstr is deleted and replaced with the
  130. classified elstr, which has boundary nodata distinguished from
  131. inner nodata */
  132. AMI_STREAM<waterWindowBaseType>*
  133. computeFlowDirections(AMI_STREAM<elevation_type>*& elstr,
  134. AMI_STREAM<elevation_type>*& filledstr,
  135. AMI_STREAM<direction_type>*& dirstr,
  136. AMI_STREAM<labelElevType> *& labeledWater) {
  137. Rtimer rt, rtTotal;
  138. AMI_STREAM<elevation_type> *elstr_reclass=NULL;
  139. AMI_STREAM<ElevationWindow > *winstr=NULL;
  140. AMI_STREAM<plateauStats> *statstr=NULL;
  141. AMI_STREAM<plateauType> *platstr=NULL;
  142. AMI_STREAM<waterType> *waterstr=NULL;
  143. AMI_STREAM<waterGridType> *mergedWaterStr=NULL;
  144. AMI_STREAM<boundaryType> *boundaryStr=NULL;
  145. AMI_STREAM<waterWindowType> *waterWindows=NULL;
  146. rt_start(rtTotal);
  147. assert(elstr && filledstr == NULL && dirstr == NULL && labeledWater == NULL);
  148. if (stats) {
  149. stats->comment("------------------------------");
  150. stats->comment("COMPUTING FLOW DIRECTIONS");
  151. /* adjust nodata -- boundary nodata distinguished from inner
  152. nodata */
  153. stats->comment("classifying nodata (inner & boundary)");
  154. }
  155. elstr_reclass = classifyNodata(elstr);
  156. delete elstr;
  157. elstr = elstr_reclass;
  158. /* ---------------------------------------------------------------------- */
  159. /* find the plateaus. */
  160. /* ---------------------------------------------------------------------- */
  161. if (stats) {
  162. stats->comment("----------", opt->verbose);
  163. stats->comment("assigning preliminary directions");
  164. }
  165. rt_start(rt);
  166. dirstr = new AMI_STREAM<direction_type>;
  167. winstr = new AMI_STREAM<ElevationWindow>();
  168. statstr = new AMI_STREAM<plateauStats>;
  169. platstr = findPlateaus(elstr, nrows, ncols, nodataType::ELEVATION_NODATA,
  170. winstr, dirstr, statstr);
  171. delete winstr; /* not used; not made */
  172. rt_stop(rt);
  173. if (stats) {
  174. stats->recordTime("findingPlateaus", rt);
  175. stats->recordLength("plateaus", platstr);
  176. stats->recordLength("plateau stats", statstr);
  177. }
  178. FILL_SAVEALL {
  179. /* printStream(*stats, statstr); */
  180. FILL_DEBUG cout << "sort plateauStr (by ij): ";
  181. AMI_STREAM<plateauType> *tmp = sort(platstr, ijCmpPlateauType());
  182. printStream2Grid(tmp, nrows, ncols,
  183. verbosedir("label1.asc"), plateauType::printLabel);
  184. delete tmp;
  185. }
  186. /* ---------------------------------------------------------------------- */
  187. /* assign labels and directions & BFS ordering. depressions have
  188. labels, but no direction information.
  189. */
  190. /* ---------------------------------------------------------------------- */
  191. rt_start(rt);
  192. waterstr = new AMI_STREAM<waterType>();
  193. assignDirections(statstr, platstr, waterstr);
  194. delete platstr;
  195. delete statstr;
  196. rt_stop(rt);
  197. if (stats) {
  198. stats->recordTime("assigning directions", rt);
  199. *stats << "maxWatershedCount=" << labelFactory::getLabelCount() << endl;
  200. }
  201. rt_start(rt);
  202. mergedWaterStr = merge2waterGrid(waterstr, dirstr, elstr);
  203. delete dirstr;
  204. delete waterstr;
  205. rt_stop(rt);
  206. if (stats) {
  207. stats->recordTime("merging", rt);
  208. stats->recordLength("mergedWaterStr", mergedWaterStr);
  209. }
  210. /* ---------------------------------------------------------------------- */
  211. /* watershed analysis */
  212. /* IN: mergedWaterStr, labelFactory::... */
  213. /* ---------------------------------------------------------------------- */
  214. if (stats) {
  215. stats->comment("--------------", opt->verbose);
  216. stats->comment("generating watersheds and watershed graph");
  217. }
  218. rt_start(rt);
  219. waterWindows = new AMI_STREAM<waterWindowType>();
  220. createWaterWindows(mergedWaterStr, nrows, ncols, waterWindows);
  221. delete mergedWaterStr;
  222. rt_stop(rt);
  223. if (stats) {
  224. stats->recordTime("creating windows", rt);
  225. stats->recordLength("waterWindows", waterWindows);
  226. }
  227. /* ---------------------------------------------------------------------- */
  228. rt_start(rt);
  229. labeledWater = new AMI_STREAM<labelElevType>();
  230. boundaryStr = new AMI_STREAM<boundaryType>();
  231. generateWatersheds(&waterWindows, nrows, ncols, labeledWater, boundaryStr);
  232. /* do we need to make boundaries here?? */
  233. delete waterWindows;
  234. /* cout << "bogus boundary length = " << boundaryStr->stream_len() << endl;*/
  235. assert(boundaryStr->stream_len() == 0);
  236. delete boundaryStr;
  237. assert(labeledWater->stream_len() == nrows * ncols);
  238. rt_stop(rt);
  239. if (stats)
  240. stats->recordTime("generating watersheds", rt);
  241. /* ---------------------------------------------------------------------- */
  242. /* find boundaries */
  243. /* ---------------------------------------------------------------------- */
  244. FILL_DEBUG cerr << "sort labeledWater (by ij):";
  245. sort(&labeledWater, ijCmpLabelElevType());
  246. #ifdef SAVE_ASCII
  247. cerr << "saving WATERSHED grid as watershed_grid\n";
  248. printStream2Grid(labeledWater, nrows, ncols,
  249. "watershed.asc", labelElevType::printLabel);
  250. #endif
  251. boundaryStr = findBoundariesMain(labeledWater);
  252. /* ---------------------------------------------------------------------- */
  253. /* filling */
  254. /* valid streams are: boundaryStr, labeledWater */
  255. /* ---------------------------------------------------------------------- */
  256. rt_start(rt);
  257. elevation_type *raise;
  258. /*find the raise elevations */
  259. FILL_DEBUG cerr << "sort boundaryStr (by elev): ";
  260. sort(&boundaryStr, elevCmpBoundaryType());
  261. raise = fill_depression(boundaryStr, labelFactory::getLabelCount());
  262. delete boundaryStr;
  263. rt_stop(rt);
  264. if (stats)
  265. stats->recordTime("filling depressions", rt);
  266. /*fill the terrain*/
  267. rt_start(rt);
  268. filledstr = new AMI_STREAM<elevation_type>();
  269. commit_fill(labeledWater, raise, labelFactory::getLabelCount(), filledstr);
  270. assert(filledstr->stream_len() == nrows * ncols);
  271. delete [] raise;
  272. rt_stop(rt);
  273. if (stats) {
  274. stats->recordTime("updating filled grid", rt);
  275. stats->recordLength("filledstr", filledstr);
  276. /* ---------------------------------------------------------------------- */
  277. /* find plateaus again and reassign directions */
  278. /* ---------------------------------------------------------------------- */
  279. stats->comment("------------------------------");
  280. stats->comment("REASSIGNING DIRECTIONS");
  281. }
  282. rt_start(rt);
  283. winstr = NULL;
  284. dirstr = new AMI_STREAM<direction_type>();
  285. statstr = new AMI_STREAM<plateauStats>();
  286. platstr = findPlateaus(filledstr, nrows, ncols, nodataType::ELEVATION_NODATA,
  287. winstr, dirstr, statstr);
  288. rt_stop(rt);
  289. if (stats) {
  290. stats->recordTime("findingPlateaus2", rt);
  291. stats->recordLength("final plateaus", platstr);
  292. stats->recordLength("final plateau stats", statstr);
  293. }
  294. FILL_SAVEALL {
  295. FILL_DEBUG cout << "sort plateauStr (by ij): ";
  296. AMI_STREAM<plateauType> *tmp = sort(platstr, ijCmpPlateauType());
  297. printStream2Grid(tmp, nrows, ncols,
  298. verbosedir("plateaus.asc"), plateauType::printLabel);
  299. delete tmp;
  300. }
  301. /* assign final directions */
  302. rt_start(rt);
  303. waterstr = new AMI_STREAM<waterType>();
  304. long dc = assignDirections(statstr, platstr, waterstr);
  305. if(dc && stats) {
  306. *stats << "WARNING: " << dc << " depressions (islands) detected\n";
  307. }
  308. delete platstr;
  309. delete statstr;
  310. rt_stop(rt);
  311. if (stats)
  312. stats->recordTime("final directions", rt);
  313. /* merge */
  314. rt_start(rt);
  315. AMI_STREAM<waterWindowBaseType> *flowStream;
  316. /*STREAM_TO_OPTION(flowStream, "flowStream");*/
  317. char path[BUFSIZ];
  318. char* base_dir = getenv(STREAM_TMPDIR);
  319. assert(base_dir);
  320. sprintf(path, "%s/flowStream", base_dir);
  321. flowStream = new AMI_STREAM<waterWindowBaseType>(path);
  322. /*flowStream->persist(PERSIST_PERSISTENT); */
  323. if (stats)
  324. stats->comment("creating flowStream: ", flowStream->sprint());
  325. merge2waterBase(waterstr, dirstr, filledstr, flowStream);
  326. delete waterstr;
  327. rt_stop(rt);
  328. if (stats)
  329. stats->recordTime("merge water,dir,elev to flow", rt);
  330. rt_stop(rtTotal);
  331. #ifdef SAVE_ASCII
  332. /*write grids as ascii file */
  333. printGridStream(filledstr, nrows, ncols,
  334. "filled_elev.asc", printElevation());
  335. printGridStream(flowStream, nrows, ncols,
  336. "direction.asc", printDirection());
  337. #endif
  338. if (stats) {
  339. stats->recordTime("Total compute flow direction running time", rtTotal);
  340. stats->comment("compute flow directions done.");
  341. }
  342. return flowStream;
  343. }
  344. /* ---------------------------------------------------------------------- */
  345. void
  346. recordWatersheds(AMI_STREAM<labelElevType> *labeledWater) {
  347. AMI_err ae;
  348. long labelCount = 0;
  349. AMI_STREAM<labelElevType> *tmp;
  350. if (stats)
  351. *stats << "--- watershed stats" << endl;
  352. FILL_DEBUG cout << "sort labeledWater (by wat label): ";
  353. tmp = sort(labeledWater, labelCmpLabelElevType());
  354. labelElevType *p;
  355. cclabel_type prev(LABEL_UNDEF);
  356. tmp->seek(0);
  357. while((ae = tmp->read_item(&p)) == AMI_ERROR_NO_ERROR) {
  358. if(p->getLabel() != prev) {
  359. labelCount++;
  360. prev = p->getLabel();
  361. }
  362. }
  363. assert(ae == AMI_ERROR_END_OF_STREAM);
  364. if (stats) {
  365. *stats << "watershed count = " << labelCount << endl;
  366. *stats << "---" << endl;
  367. stats->flush();
  368. }
  369. delete tmp;
  370. }
  371. /* ********************************************************************** */
  372. /* assign directions to plateaus that have sinks;
  373. * reassign labels to depressions (don't drain out).
  374. * all plateaus are written out to the waterstr. */
  375. long
  376. assignDirections(AMI_STREAM<plateauStats> *statstr,
  377. AMI_STREAM<plateauType> *platstr,
  378. AMI_STREAM<waterType> *waterstr) {
  379. size_t fmem;
  380. AMI_err ae;
  381. plateauStats *ps;
  382. if (stats) {
  383. stats->comment("----------", opt->verbose);
  384. stats->comment("assigning directions on plateaus");
  385. }
  386. labelFactory::reset(); /* we are relabeling now */
  387. statstr->seek(0);
  388. platstr->seek(0);
  389. fmem = getAvailableMemory();
  390. long depressionCount=0;
  391. long spillCount=0;
  392. while((ae = statstr->read_item(&ps)) == AMI_ERROR_NO_ERROR) {
  393. if(ps->size*sizeof(gridElement) > fmem) {
  394. cerr << "WARNING: grid larger than memory (ignored)" << endl;
  395. }
  396. assert(ps->label != LABEL_NODATA);
  397. if(ps->hasSpill) {
  398. spillCount++;
  399. grid *platGrid = new grid(ps->iMin, ps->jMin, ps->iMax, ps->jMax,
  400. ps->size, ps->label);
  401. platGrid->load(*platstr);
  402. platGrid->assignDirections(opt->d8 ? 1 : 0);
  403. platGrid->save(*waterstr); /* this doesn't save labels */
  404. delete platGrid;
  405. } else {
  406. /* depression - just give contiguous labels only */
  407. depressionCount++;
  408. cclabel_type label = labelFactory::getNewLabel();
  409. for(int i=0; i<ps->size; i++) {
  410. plateauType *pt;
  411. platstr->read_item(&pt);
  412. pt->cclabel = label; /* assign new label */
  413. waterType wt(*pt); /* write it out */
  414. ae = waterstr->write_item(wt);
  415. assert(ae == AMI_ERROR_NO_ERROR);
  416. }
  417. }
  418. }
  419. if (stats) {
  420. *stats << "depression count = " << depressionCount << endl;
  421. *stats << "spill count = " << spillCount << endl;
  422. }
  423. return depressionCount;
  424. }
  425. /* ********************************************************************** */
  426. /* assign directions to plateaus that have sinks;
  427. * check that there are no depressions.
  428. */
  429. void
  430. assignFinalDirections(AMI_STREAM<plateauStats> *statstr,
  431. AMI_STREAM<plateauType> *platstr,
  432. AMI_STREAM<waterType> *waterstr) {
  433. AMI_err ae;
  434. plateauStats *ps;
  435. if (stats)
  436. stats->comment("assigning final directions");
  437. statstr->seek(0);
  438. platstr->seek(0);
  439. while((ae = statstr->read_item(&ps)) == AMI_ERROR_NO_ERROR) {
  440. if(ps->hasSpill) {
  441. grid *platGrid = new grid(ps->iMin, ps->jMin, ps->iMax, ps->jMax,
  442. ps->size, ps->label);
  443. platGrid->load(*platstr);
  444. platGrid->assignDirections(opt->d8 ? 1 : 0);
  445. platGrid->save(*waterstr); /* this doesn't save labels */
  446. delete platGrid;
  447. } else {
  448. /* could be legitimate */
  449. cerr << "WARNING: depression detected: " << *ps << endl;
  450. continue;
  451. }
  452. }
  453. };
  454. /* ********************************************************************** */
  455. class directionElevationMerger {
  456. public:
  457. waterGridType operator()(elevation_type el, direction_type dir,
  458. waterType p) {
  459. /* check that no (boundary) nodata values got in here */
  460. assert(el != nodataType::ELEVATION_BOUNDARY);
  461. assert(!is_nodata(el)); /* p should be a valid grid cell */
  462. return waterGridType(el, p.getDirection(), p.getLabel(), p.getDepth());
  463. }
  464. waterGridType operator()(elevation_type el, direction_type dir) {
  465. waterGridType wg(el, dir);
  466. if(el == nodataType::ELEVATION_BOUNDARY) { /* hack XXX (approved) */
  467. wg.setLabel(LABEL_BOUNDARY);
  468. }
  469. /* nodata => boundary or undef */
  470. assert(!is_nodata(el) ||
  471. (wg.getLabel()==LABEL_BOUNDARY || wg.getLabel()==LABEL_UNDEF));
  472. return wg;
  473. }
  474. };
  475. /* ---------------------------------------------------------------------- */
  476. AMI_STREAM<waterGridType> *
  477. merge2waterGrid(AMI_STREAM<waterType> *unsortedWaterStr,
  478. AMI_STREAM<direction_type> *dirStr,
  479. AMI_STREAM<elevation_type> *elStr) {
  480. AMI_STREAM<waterType> *waterStr;
  481. FILL_DEBUG cout << "sort waterStr (by ij): ";
  482. waterStr = sort(unsortedWaterStr, ijCmpWaterType());
  483. FILL_SAVEALL printStream2Grid(waterStr, nrows, ncols,
  484. verbosedir("platlabels.asc"), printLabel());
  485. AMI_STREAM<waterGridType> *mergedWaterStr = new AMI_STREAM<waterGridType>();
  486. mergeStreamGridGrid(elStr, dirStr,
  487. nrows, ncols,
  488. waterStr,
  489. directionElevationMerger(),
  490. mergedWaterStr);
  491. delete waterStr;
  492. FILL_SAVEALL printGridStream(mergedWaterStr, nrows, ncols,
  493. verbosedir("mergedlabels.asc"), printLabel());
  494. assert(mergedWaterStr->stream_len());
  495. return mergedWaterStr;
  496. }
  497. /* ---------------------------------------------------------------------- */
  498. void
  499. merge2waterBase(AMI_STREAM<waterType> *unsortedWaterStr,
  500. AMI_STREAM<direction_type> *dirStr,
  501. AMI_STREAM<elevation_type> *elStr,
  502. AMI_STREAM<waterWindowBaseType> *merge) {
  503. AMI_STREAM<waterType> *waterStr;
  504. FILL_DEBUG cout << "sort waterStr (by ij): ";
  505. waterStr = sort(unsortedWaterStr, ijCmpWaterType());
  506. mergeStreamGridGrid(elStr, dirStr,
  507. nrows, ncols,
  508. waterStr,
  509. directionElevationMerger(),
  510. merge);
  511. delete waterStr;
  512. }
  513. /* ---------------------------------------------------------------------- */
  514. /*
  515. * merge 2 grids and a stream together to form a new grid.
  516. * str should be sorted in ij order
  517. */
  518. template<class T1, class T2, class T3, class T4, class FUN>
  519. void
  520. mergeStreamGridGrid(AMI_STREAM<T1> *grid1,
  521. AMI_STREAM<T2> *grid2,
  522. dimension_type rows, dimension_type cols,
  523. AMI_STREAM<T3> *str,
  524. FUN fo,
  525. AMI_STREAM<T4> *outStream) {
  526. T1 *t1p;
  527. T2 *t2p;
  528. T3 *t3p;
  529. AMI_err aeUpd, ae;
  530. grid1->seek(0);
  531. grid2->seek(0);
  532. str->seek(0);
  533. aeUpd = str->read_item(&t3p);
  534. assert(aeUpd == AMI_ERROR_NO_ERROR || aeUpd == AMI_ERROR_END_OF_STREAM);
  535. for(dimension_type row = 0; row < rows; row++) {
  536. for(dimension_type col = 0; col < cols; col++) {
  537. ae = grid1->read_item(&t1p);
  538. assert(ae == AMI_ERROR_NO_ERROR);
  539. ae = grid2->read_item(&t2p);
  540. assert(ae == AMI_ERROR_NO_ERROR);
  541. T4 t4;
  542. if(aeUpd == AMI_ERROR_NO_ERROR && t3p->i == row && t3p->j == col) {
  543. /* cell present in stream */
  544. t4 = fo(*t1p, *t2p, *t3p);
  545. aeUpd = str->read_item(&t3p);
  546. assert(aeUpd == AMI_ERROR_NO_ERROR ||
  547. aeUpd == AMI_ERROR_END_OF_STREAM);
  548. } else {
  549. /* not in stream */
  550. t4 = fo(*t1p, *t2p);
  551. }
  552. ae = outStream->write_item(t4);
  553. assert(ae == AMI_ERROR_NO_ERROR);
  554. }
  555. /*assert(outStream->stream_len() == (row+1) * cols); */
  556. }
  557. assert(outStream->stream_len() == rows * cols);
  558. return;
  559. }
  560. /* ---------------------------------------------------------------------- */
  561. /* make boundaryStr from labeledWater */
  562. AMI_STREAM<boundaryType> *
  563. findBoundariesMain(AMI_STREAM<labelElevType> *labeledWater) {
  564. AMI_STREAM<boundaryType> *boundaryStr;
  565. Rtimer rt;
  566. rt_start(rt);
  567. boundaryStr = new AMI_STREAM<boundaryType>();
  568. FILL_SAVEALL printGridStream(labeledWater, nrows, ncols,
  569. verbosedir("labels.asc"), printLabel());
  570. findBoundaries(labeledWater, nrows, ncols, boundaryStr);
  571. if (stats)
  572. stats->recordLength("all boundaries", boundaryStr);
  573. FILL_SAVEALL {
  574. FILL_DEBUG cout << "sort boundaryStr (by ij): ";
  575. sort(&boundaryStr, ijCmpBoundaryType());
  576. removeDuplicatesEx(&boundaryStr, ijCmpBoundaryType());
  577. printStream2Grid(boundaryStr, nrows, ncols,
  578. verbosedir("boundary.asc"), boundaryType::print);
  579. }
  580. FILL_DEBUG cout << "sort boundaryStr (by wat label): ";
  581. sort(&boundaryStr, waterCmpBoundaryType());
  582. removeDuplicatesEx(&boundaryStr, boundaryCmpBoundaryType());
  583. rt_stop(rt);
  584. if (stats) {
  585. stats->recordTime("generating boundaries", rt);
  586. stats->recordLength("boundary stream", boundaryStr);
  587. }
  588. /*if(GETOPT("veryfillverbose")) printStream(cout, boundaryStr);
  589. SAVE_ON_OPTION(boundaryStr, "saveBoundaryStream");
  590. */
  591. return boundaryStr;
  592. }