CoastalME (Coastal Modelling Environment)
Simulates the long-term behaviour of complex coastlines
Loading...
Searching...
No Matches
gis_utils.cpp
Go to the documentation of this file.
1
22
23/* ===============================================================================================================================
24 This file is part of CoastalME, the Coastal Modelling Environment. CoastalME is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.
25
26 This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
27
28 You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
29===============================================================================================================================*/
30#include <assert.h>
31
32#include <cstdio>
33
34#include <vector>
35using std::vector;
36
37#include <iostream>
38using std::cerr;
39using std::endl;
40using std::ios;
41
42#include <cmath>
43using std::atan2;
44
45#include <cfloat>
46#include <climits>
47#include <cstdint>
48
49#include <cstring>
50using std::strstr;
51
52#include <gdal.h>
53#include <gdal_priv.h>
54#include <cpl_string.h>
55
56#include "cme.h"
57#include "coast.h"
58#include "raster_grid.h"
59#include "2d_point.h"
60#include "2di_point.h"
61
62//===============================================================================================================================
64//===============================================================================================================================
65double CSimulation::dGridCentroidXToExtCRSX(int const nGridX) const
66{
67 // TODO 064
68 return (m_dGeoTransform[0] + (nGridX * m_dGeoTransform[1]) + (m_dGeoTransform[1] / 2));
69}
70
71//===============================================================================================================================
73//===============================================================================================================================
74double CSimulation::dGridCentroidYToExtCRSY(int const nGridY) const
75{
76 // TODO 064
77 return (m_dGeoTransform[3] + (nGridY * m_dGeoTransform[5]) + (m_dGeoTransform[5] / 2));
78}
79
80//===============================================================================================================================
82//===============================================================================================================================
84{
85 // TODO 064
86 double const dX = m_dGeoTransform[0] + (pPtiIn->nGetX() * m_dGeoTransform[1]) + (m_dGeoTransform[1] / 2);
87 double const dY = m_dGeoTransform[3] + (pPtiIn->nGetY() * m_dGeoTransform[5]) + (m_dGeoTransform[5] / 2);
88
89 return CGeom2DPoint(dX, dY);
90}
91
92//===============================================================================================================================
94//===============================================================================================================================
95double CSimulation::dGridXToExtCRSX(double const dGridX) const
96{
97 // TODO 064 Xgeo = GT(0) + Xpixel*GT(1) + Yline*GT(2)
98 return m_dGeoTransform[0] + (dGridX * m_dGeoTransform[1]) - 1;
99}
100
101//===============================================================================================================================
103//===============================================================================================================================
104double CSimulation::dGridYToExtCRSY(double const dGridY) const
105{
106 // TODO 064 Ygeo = GT(3) + Xpixel*GT(4) + Yline*GT(5)
107 return m_dGeoTransform[3] + (dGridY * m_dGeoTransform[5]) - 1;
108}
109
110//===============================================================================================================================
112//===============================================================================================================================
113double CSimulation::dExtCRSXToGridX(double const dExtCRSX) const
114{
115 // TODO 064
116 return ((dExtCRSX - m_dGeoTransform[0]) / m_dGeoTransform[1]) - 1;
117}
118
119//===============================================================================================================================
121//===============================================================================================================================
122double CSimulation::dExtCRSYToGridY(double const dExtCRSY) const
123{
124 // TODO 064
125 return ((dExtCRSY - m_dGeoTransform[3]) / m_dGeoTransform[5]) - 1;
126}
127
128//===============================================================================================================================
130//===============================================================================================================================
132{
133 // TODO 064
134 int const nX = nRound(((pPtIn->dGetX() - m_dGeoTransform[0]) / m_dGeoTransform[1]) - 1);
135 int const nY = nRound(((pPtIn->dGetY() - m_dGeoTransform[3]) / m_dGeoTransform[5]) - 1);
136
137 return CGeom2DIPoint(nX, nY);
138}
139
140//===============================================================================================================================
142//===============================================================================================================================
144{
145 double const dXDist = Pt1->dGetX() - Pt2->dGetX();
146 double const dYDist = Pt1->dGetY() - Pt2->dGetY();
147
148 return hypot(dXDist, dYDist);
149}
150
151//===============================================================================================================================
153//===============================================================================================================================
155{
156 double const dXDist = Pti1->nGetX() - Pti2->nGetX();
157 double const dYDist = Pti1->nGetY() - Pti2->nGetY();
158
159 return hypot(dXDist, dYDist);
160}
161
162//===============================================================================================================================
164//===============================================================================================================================
165double CSimulation::dGetDistanceBetween(double const dX1, double const dY1, double const dX2, double const dY2)
166{
167 double const dXDist = dX1 - dX2;
168 double const dYDist = dY1 - dY2;
169
170 return hypot(dXDist, dYDist);
171}
172
173//===============================================================================================================================
175//===============================================================================================================================
176double CSimulation::dTriangleAreax2(CGeom2DPoint const* pPtA, CGeom2DPoint const* pPtB, CGeom2DPoint const* pPtC)
177{
178 return (pPtB->dGetX() - pPtA->dGetX()) * (pPtC->dGetY() - pPtA->dGetY()) - (pPtB->dGetY() - pPtA->dGetY()) * (pPtC->dGetX() - pPtA->dGetX());
179}
180
181//===============================================================================================================================
183//===============================================================================================================================
184bool CSimulation::bIsWithinValidGrid(int const nX, int const nY) const
185{
186 if ((nX < 0) || (nX >= m_nXGridSize))
187 return false;
188
189 if ((nY < 0) || (nY >= m_nYGridSize))
190 return false;
191
192 if (m_pRasterGrid->m_Cell[nX][nY].bBasementElevIsMissingValue())
193 return false;
194
195 return true;
196}
197
198//===============================================================================================================================
200//===============================================================================================================================
202{
203 int const nX = Pti->nGetX();
204 int const nY = Pti->nGetY();
205
206 return this->bIsWithinValidGrid(nX, nY);
207}
208
209//===============================================================================================================================
211//===============================================================================================================================
212void CSimulation::KeepWithinValidGrid(int& nX, int& nY) const
213{
214 nX = tMax(nX, 0);
215 nX = tMin(nX, m_nXGridSize - 1);
216
217 nY = tMax(nY, 0);
218 nY = tMin(nY, m_nYGridSize - 1);
219}
220
221//===============================================================================================================================
223//===============================================================================================================================
225{
226 KeepWithinValidGrid(Pti0->nGetX(), Pti0->nGetY(), *Pti1->pnGetX(), *Pti1->pnGetY());
227}
228
229//===============================================================================================================================
232//===============================================================================================================================
233void CSimulation::KeepWithinValidGrid(int nX0, int nY0, int& nX1, int& nY1) const
234{
235 // Safety check: make sure that the first point is within the valid grid
236 if (nX0 >= m_nXGridSize)
237 nX0 = m_nXGridSize - 1;
238
239 else if (nX0 < 0)
240 nX0 = 0;
241
242 if (nY0 >= m_nYGridSize)
243 nY0 = m_nYGridSize - 1;
244
245 else if (nY0 < 0)
246 nY0 = 0;
247
248 // OK let's go
249 int const nDiffX = nX0 - nX1;
250 int const nDiffY = nY0 - nY1;
251
252 if (nDiffX == 0)
253 {
254 // The two points have the same x coordinates, so we just need to constrain the y co-ord
255 if (nY1 < nY0)
256 {
257 nY1 = -1;
258
259 do
260 {
261 nY1++;
262 } while (m_pRasterGrid->m_Cell[nX1][nY1].bBasementElevIsMissingValue());
263
264 return;
265 }
266 else
267 {
268 nY1 = m_nYGridSize;
269
270 do
271 {
272 nY1--;
273 } while (m_pRasterGrid->m_Cell[nX1][nY1].bBasementElevIsMissingValue());
274
275 return;
276 }
277 }
278 else if (nDiffY == 0)
279 {
280 // The two points have the same y coordinates, so we just need to constrain the x co-ord
281 if (nX1 < nX0)
282 {
283 nX1 = -1;
284
285 do
286 {
287 nX1++;
288 } while (m_pRasterGrid->m_Cell[nX1][nY1].bBasementElevIsMissingValue());
289
290 return;
291 }
292 else
293 {
294 nX1 = m_nXGridSize;
295
296 do
297 {
298 nX1--;
299 } while (m_pRasterGrid->m_Cell[nX1][nY1].bBasementElevIsMissingValue());
300
301 return;
302 }
303 }
304 else
305 {
306 // The two points have different x coordinates and different y coordinates, so we have to work harder. First find which of the coordinates is the greatest distance outside the grid, and constrain that co-ord for efficiency (since this will reduce the number of times round the loop). Note that both may be inside the grid, if the incorrect co-ord is in the invalid margin, in which case arbitrarily contrain the x co-ord
307 int nXDistanceOutside = 0;
308 int nYDistanceOutside = 0;
309
310 if (nX1 < 0)
311 nXDistanceOutside = -nX1;
312 else if (nX1 >= m_nXGridSize)
313 nXDistanceOutside = nX1 - m_nXGridSize + 1;
314
315 if (nY1 < 0)
316 nYDistanceOutside = -nY1;
317 else if (nY1 >= m_nYGridSize)
318 nXDistanceOutside = nY1 - m_nYGridSize + 1;
319
320 if (nXDistanceOutside >= nYDistanceOutside)
321 {
322 // Constrain the x co-ord
323 if (nX1 < nX0)
324 {
325 // The incorrect x co-ord is less than the correct x co-ord: constrain it and find the y co-ord
326 nX1 = -1;
327
328 do
329 {
330 nX1++;
331
332 nY1 = nY0 + nRound(((nX1 - nX0) * nDiffY) / static_cast<double>(nDiffX));
333 } while ((nY1 < 0) || (nY1 >= m_nYGridSize) || (m_pRasterGrid->m_Cell[nX1][nY1].bBasementElevIsMissingValue()));
334
335 return;
336 }
337 else
338 {
339 // The incorrect x co-ord is greater than the correct x-co-ord: constrain it and find the y co-ord
340 nX1 = m_nXGridSize;
341
342 do
343 {
344 nX1--;
345
346 nY1 = nY0 + nRound(((nX1 - nX0) * nDiffY) / static_cast<double>(nDiffX));
347 } while ((nY1 < 0) || (nY1 >= m_nYGridSize) || (m_pRasterGrid->m_Cell[nX1][nY1].bBasementElevIsMissingValue()));
348
349 return;
350 }
351 }
352 else
353 {
354 // Constrain the y co-ord
355 if (nY1 < nY0)
356 {
357 // The incorrect y co-ord is less than the correct y-co-ord: constrain it and find the x co-ord
358 nY1 = -1;
359
360 do
361 {
362 nY1++;
363
364 nX1 = nX0 + nRound(((nY1 - nY0) * nDiffX) / static_cast<double>(nDiffY));
365 } while ((nX1 < 0) || (nX1 >= m_nXGridSize) || (m_pRasterGrid->m_Cell[nX1][nY1].bBasementElevIsMissingValue()));
366
367 return;
368 }
369 else
370 {
371 // The incorrect y co-ord is greater than the correct y co-ord: constrain it and find the x co-ord
372 nY1 = m_nYGridSize;
373
374 do
375 {
376 nY1--;
377
378 nX1 = nX0 +
379 nRound(((nY1 - nY0) * nDiffX) / static_cast<double>(nDiffY));
380 } while ((nX1 < 0) || (nX1 >= m_nXGridSize) || (m_pRasterGrid->m_Cell[nX1][nY1].bBasementElevIsMissingValue()));
381
382 return;
383 }
384 }
385 }
386}
387
388//===============================================================================================================================
390//===============================================================================================================================
391double CSimulation::dKeepWithin360(double const dAngle)
392{
393 double dNewAngle = dAngle;
394
395 // Sort out -ve angles
396 while (dNewAngle < 0)
397 dNewAngle += 360;
398
399 // Sort out angles > 360
400 while (dNewAngle > 360)
401 dNewAngle -= 360;
402
403 return dNewAngle;
404}
405
406//===============================================================================================================================
408//===============================================================================================================================
410{
411 double const dPt1X = pPt1->dGetX();
412 double const dPt1Y = pPt1->dGetY();
413 double const dPt2X = pPt2->dGetX();
414 double const dPt2Y = pPt2->dGetY();
415 double const dPtAvgX = (dPt1X + dPt2X) / 2;
416 double const dPtAvgY = (dPt1Y + dPt2Y) / 2;
417
418 return CGeom2DPoint(dPtAvgX, dPtAvgY);
419}
420
421// //===============================================================================================================================
422// //! Returns an integer point (grid CRS) which is the approximate average of (i.e. is midway between) two other grid CRS integer points
423// //===============================================================================================================================
424// CGeom2DIPoint CSimulation::PtiAverage(CGeom2DIPoint const* pPti1, CGeom2DIPoint const* pPti2)
425// {
426// int const nPti1X = pPti1->nGetX();
427// int const nPti1Y = pPti1->nGetY();
428// int const nPti2X = pPti2->nGetX();
429// int const nPti2Y = pPti2->nGetY();
430// int const nPtiAvgX = (nPti1X + nPti2X) / 2;
431// int const nPtiAvgY = (nPti1Y + nPti2Y) / 2;
432//
433// return CGeom2DIPoint(nPtiAvgX, nPtiAvgY);
434// }
435
436//===============================================================================================================================
438//===============================================================================================================================
439CGeom2DIPoint CSimulation::PtiWeightedAverage(CGeom2DIPoint const* pPti1, CGeom2DIPoint const* pPti2, double const dWeight)
440{
441 int const nPti1X = pPti1->nGetX();
442 int const nPti1Y = pPti1->nGetY();
443 int const nPti2X = pPti2->nGetX();
444 int const nPti2Y = pPti2->nGetY();
445 double const dOtherWeight = 1.0 - dWeight;
446
447 int const nPtiWeightAvgX = nRound((dWeight * nPti2X) + (dOtherWeight * nPti1X));
448 int const nPtiWeightAvgY = nRound((dWeight * nPti2Y) + (dOtherWeight * nPti1Y));
449
450 return CGeom2DIPoint(nPtiWeightAvgX, nPtiWeightAvgY);
451}
452
453//===============================================================================================================================
455//===============================================================================================================================
456CGeom2DPoint CSimulation::PtAverage(vector<CGeom2DPoint>* pVIn)
457{
458 int const nSize = static_cast<int>(pVIn->size());
459
460 if (nSize == 0)
462
463 double dAvgX = 0;
464 double dAvgY = 0;
465
466 for (int n = 0; n < nSize; n++)
467 {
468 dAvgX += pVIn->at(n).dGetX();
469 dAvgY += pVIn->at(n).dGetY();
470 }
471
472 dAvgX /= nSize;
473 dAvgY /= nSize;
474
475 return CGeom2DPoint(dAvgX, dAvgY);
476}
477
478// //===============================================================================================================================
479// //! Returns a point (grid CRS) which is the average of a vector of grid CRS points
480// //===============================================================================================================================
481// CGeom2DIPoint CSimulation::PtiAverage(vector<CGeom2DIPoint>* pVIn)
482// {
483// int nSize = static_cast<int>(pVIn->size());
484// if (nSize == 0)
485// return CGeom2DIPoint(INT_NODATA, INT_NODATA);
486//
487// double dAvgX = 0;
488// double dAvgY = 0;
489//
490// for (int n = 0; n < nSize; n++)
491// {
492// dAvgX += pVIn->at(n).nGetX();
493// dAvgY += pVIn->at(n).nGetY();
494// }
495//
496// dAvgX /= nSize;
497// dAvgY /= nSize;
498//
499// return CGeom2DIPoint(nRound(dAvgX), nRound(dAvgY));
500// }
501
502//===============================================================================================================================
504//===============================================================================================================================
506{
507 CGeom2DIPoint PtiCentroid(0, 0);
508 int const nSize = static_cast<int>(pVIn->size());
509 int nX0 = 0; // Current vertex X
510 int nY0 = 0; // Current vertex Y
511 int nX1 = 0; // Next vertex X
512 int nY1 = 0; // Next vertex Y
513
514 double dA = 0; // Partial signed area
515 double dSignedArea = 0.0;
516
517 // For all vertices except last
518 for (int i = 0; i < nSize - 1; ++i)
519 {
520 nX0 = pVIn->at(i).nGetX();
521 nY0 = pVIn->at(i).nGetY();
522 nX1 = pVIn->at(i + 1).nGetX();
523 nY1 = pVIn->at(i + 1).nGetY();
524
525 dA = (nX0 * nY1) - (nX1 * nY0);
526 dSignedArea += dA;
527 PtiCentroid.AddXAddY((nX0 + nX1) * dA, (nY0 + nY1) * dA);
528 }
529
530 // Do last vertex separately to avoid performing an expensive modulus operation in each iteration
531 nX0 = pVIn->at(nSize - 1).nGetX();
532 nY0 = pVIn->at(nSize - 1).nGetY();
533 nX1 = pVIn->at(0).nGetX();
534 nY1 = pVIn->at(0).nGetY();
535
536 dA = (nX0 * nY1) - (nX1 * nY0);
537 dSignedArea += dA;
538 PtiCentroid.AddXAddY((nX0 + nX1) * dA, (nY0 + nY1) * dA);
539
540 dSignedArea *= 0.5;
541 PtiCentroid.DivXDivY(6.0 * dSignedArea, 6.0 * dSignedArea);
542
543 return PtiCentroid;
544}
545
546//===============================================================================================================================
547// Returns a vector which is perpendicular to an existing vector
548//===============================================================================================================================
549// vector<CGeom2DPoint> CSimulation::VGetPerpendicular(CGeom2DPoint const*
550// PtStart, CGeom2DPoint const* PtNext, double const dDesiredLength, int const
551// nHandedness)
552// {
553// // Returns a two-point vector which passes through PtStart with a scaled
554// length
555// double dXLen = PtNext->dGetX() - PtStart->dGetX();
556// double dYLen = PtNext->dGetY() - PtStart->dGetY();
557//
558// double dLength = hypot(dXLen, dYLen);
559// double dScaleFactor = dDesiredLength / dLength;
560//
561// // The difference vector is (dXLen, dYLen), so the perpendicular
562// difference vector is (-dYLen, dXLen) or (dYLen, -dXLen)
563// CGeom2DPoint EndPt;
564// if (nHandedness == RIGHT_HANDED)
565// {
566// EndPt.SetX(PtStart->dGetX() + (dScaleFactor * dYLen));
567// EndPt.SetY(PtStart->dGetY() - (dScaleFactor * dXLen));
568// }
569// else
570// {
571// EndPt.SetX(PtStart->dGetX() - (dScaleFactor * dYLen));
572// EndPt.SetY(PtStart->dGetY() + (dScaleFactor * dXLen));
573// }
574//
575// vector<CGeom2DPoint> VNew;
576// VNew.push_back(*PtStart);
577// VNew.push_back(EndPt);
578// return VNew;
579// }
580
581// //===============================================================================================================================
582// //! Returns a CGeom2DPoint which is the 'other' point of a two-point vector passing through PtStart, and which is perpendicular to the two-point vector from PtStart to PtNext
583// //===============================================================================================================================
584// CGeom2DPoint CSimulation::PtGetPerpendicular(CGeom2DPoint const* PtStart, CGeom2DPoint const* PtNext, double const dDesiredLength, int const nHandedness)
585// {
586// double const dXLen = PtNext->dGetX() - PtStart->dGetX();
587// double const dYLen = PtNext->dGetY() - PtStart->dGetY();
588// double dLength;
589//
590// if (bFPIsEqual(dXLen, 0.0, TOLERANCE))
591// dLength = dYLen;
592// else if (bFPIsEqual(dYLen, 0.0, TOLERANCE))
593// dLength = dXLen;
594// else
595// dLength = hypot(dXLen, dYLen);
596//
597// double const dScaleFactor = dDesiredLength / dLength;
598//
599// // The difference vector is (dXLen, dYLen), so the perpendicular difference vector is (-dYLen, dXLen) or (dYLen, -dXLen)
600// CGeom2DPoint EndPt;
601//
602// if (nHandedness == RIGHT_HANDED)
603// {
604// EndPt.SetX(PtStart->dGetX() + (dScaleFactor * dYLen));
605// EndPt.SetY(PtStart->dGetY() - (dScaleFactor * dXLen));
606// }
607// else
608// {
609// EndPt.SetX(PtStart->dGetX() - (dScaleFactor * dYLen));
610// EndPt.SetY(PtStart->dGetY() + (dScaleFactor * dXLen));
611// }
612//
613// return EndPt;
614// }
615
616//===============================================================================================================================
618//===============================================================================================================================
619CGeom2DIPoint CSimulation::PtiGetPerpendicular(CGeom2DIPoint const* PtiStart, CGeom2DIPoint const* PtiNext, double const dDesiredLength, int const nHandedness)
620{
621 double const dXLen = PtiNext->nGetX() - PtiStart->nGetX();
622 double const dYLen = PtiNext->nGetY() - PtiStart->nGetY();
623 double dLength;
624
625 if (bFPIsEqual(dXLen, 0.0, TOLERANCE))
626 dLength = dYLen;
627
628 else if (bFPIsEqual(dYLen, 0.0, TOLERANCE))
629 dLength = dXLen;
630
631 else
632 dLength = hypot(dXLen, dYLen);
633
634 double const dScaleFactor = dDesiredLength / dLength;
635
636 // The difference vector is (dXLen, dYLen), so the perpendicular difference vector is (-dYLen, dXLen) or (dYLen, -dXLen)
637 CGeom2DIPoint EndPti;
638
639 if (nHandedness == RIGHT_HANDED)
640 {
641 EndPti.SetX(PtiStart->nGetX() + nRound(dScaleFactor * dYLen));
642 EndPti.SetY(PtiStart->nGetY() - nRound(dScaleFactor * dXLen));
643 }
644
645 else
646 {
647 EndPti.SetX(PtiStart->nGetX() - nRound(dScaleFactor * dYLen));
648 EndPti.SetY(PtiStart->nGetY() + nRound(dScaleFactor * dXLen));
649 }
650
651 return EndPti;
652}
653
654//===============================================================================================================================
656//===============================================================================================================================
657CGeom2DIPoint CSimulation::PtiGetPerpendicular(int const nStartX, int const nStartY, int const nNextX, int const nNextY, double const dDesiredLength, int const nHandedness)
658{
659 double const dXLen = nNextX - nStartX;
660 double const dYLen = nNextY - nStartY;
661 double dLength;
662
663 if (bFPIsEqual(dXLen, 0.0, TOLERANCE))
664 dLength = dYLen;
665
666 else if (bFPIsEqual(dYLen, 0.0, TOLERANCE))
667 dLength = dXLen;
668
669 else
670 dLength = hypot(dXLen, dYLen);
671
672 double const dScaleFactor = dDesiredLength / dLength;
673
674 // The difference vector is (dXLen, dYLen), so the perpendicular difference vector is (-dYLen, dXLen) or (dYLen, -dXLen)
675 CGeom2DIPoint EndPti;
676
677 if (nHandedness == RIGHT_HANDED)
678 {
679 EndPti.SetX(nStartX + nRound(dScaleFactor * dYLen));
680 EndPti.SetY(nStartY - nRound(dScaleFactor * dXLen));
681 }
682
683 else
684 {
685 EndPti.SetX(nStartX - nRound(dScaleFactor * dYLen));
686 EndPti.SetY(nStartY + nRound(dScaleFactor * dXLen));
687 }
688
689 return EndPti;
690}
691
692//===============================================================================================================================
694//===============================================================================================================================
695double CSimulation::dAngleSubtended(CGeom2DIPoint const* pPtiA, CGeom2DIPoint const* pPtiB, CGeom2DIPoint const* pPtiC)
696{
697 double const dXDistBtoA = pPtiB->nGetX() - pPtiA->nGetX();
698 double const dYDistBtoA = pPtiB->nGetY() - pPtiA->nGetY();
699 double const dXDistCtoA = pPtiC->nGetX() - pPtiA->nGetX();
700 double const dYDistCtoA = pPtiC->nGetY() - pPtiA->nGetY();
701 double const dDotProduct = dXDistBtoA * dXDistCtoA + dYDistBtoA * dYDistCtoA;
702 double const dPseudoCrossProduct = dXDistBtoA * dYDistCtoA - dYDistBtoA * dXDistCtoA;
703 double const dAngle = atan2(dPseudoCrossProduct, dDotProduct);
704
705 return dAngle;
706}
707
708//===============================================================================================================================
710//===============================================================================================================================
712{
713 // Register all available GDAL raster and vector drivers (GDAL 2)
714 GDALAllRegister();
715
716 // If the user hasn't specified a GIS output format, assume that we will use the same GIS format as the input basement DEM
717 if (m_strRasterGISOutFormat.empty())
719
720 // Load the raster GDAL driver
721 GDALDriver *pDriver = GetGDALDriverManager()->GetDriverByName(m_strRasterGISOutFormat.c_str());
722
723 if (NULL == pDriver)
724 {
725 // Can't load raster GDAL driver. Incorrectly specified?
726 cerr << ERR << "Unknown raster GIS output format '" << m_strRasterGISOutFormat << "'." << endl;
727 return false;
728 }
729
730 // Get the metadata for this raster driver
731 char **papszMetadata = pDriver->GetMetadata();
732
733 // for (int i = 0; papszMetadata[i] != NULL; i++)
734 // cout << papszMetadata[i] << endl;
735 // cout << endl;
736
737 // Need to test if this is a raster driver
738 if (!CSLFetchBoolean(papszMetadata, GDAL_DCAP_RASTER, false))
739 {
740 // This is not a raster driver
741 cerr << ERR << "GDAL driver '" << m_strRasterGISOutFormat << "' is not a raster driver. Choose another format." << endl;
742 return false;
743 }
744
745 // This driver is OK, so store its longname and the default file extension
746 string strTmp = CSLFetchNameValue(papszMetadata, "DMD_LONGNAME");
748 strTmp = CSLFetchNameValue(papszMetadata, "DMD_EXTENSIONS"); // Note DMD_EXTENSION (no S, is a single value) appears not to be implemented for newer drivers
749 strTmp = strTrim(&strTmp);
750
751 // We have a space-separated list of one or more file extensions: use the first extension in the list
752 long unsigned int const nPos = strTmp.find(SPACE);
753
754 if (nPos == string::npos)
755 {
756 // No space i.e. just one extension
758 }
759
760 else
761 {
762 // There's a space, so we must have more than one extension
763 m_strGDALRasterOutputDriverExtension = strTmp.substr(0, nPos);
764 }
765
766 // Set up any defaults for raster files that are created using this driver
768
769 // Now do various tests of the driver's capabilities
770 if (!CSLFetchBoolean(papszMetadata, GDAL_DCAP_CREATE, false))
771 {
772 // This raster driver does not support the Create() method, does it support CreateCopy()?
773 if (!CSLFetchBoolean(papszMetadata, GDAL_DCAP_CREATECOPY, false))
774 {
775 cerr << ERR << "Cannot write using raster GDAL driver '" << m_strRasterGISOutFormat << " since neither Create() or CreateCopy() are supported'. Choose another GDAL raster format." << endl;
776 return false;
777 }
778
779 // Can't use Create() but can use CreateCopy()
780 m_bGDALCanCreate = false;
781 }
782
783 // Next, test to see what data types the driver can write and from this, work out the largest int and float we can write
784 if (strstr(CSLFetchNameValue(papszMetadata, "DMD_CREATIONDATATYPES"), "Float"))
785 {
787 m_GDALWriteFloatDataType = GDT_Float32;
788 }
789
790 if (strstr(CSLFetchNameValue(papszMetadata, "DMD_CREATIONDATATYPES"), "UInt32"))
791 {
793
794 m_GDALWriteIntDataType = GDT_UInt32;
795 m_lGDALMaxCanWrite = UINT32_MAX;
797
799 m_GDALWriteFloatDataType = GDT_UInt32;
800
801 return true;
802 }
803
804 if (strstr(CSLFetchNameValue(papszMetadata, "DMD_CREATIONDATATYPES"), "Int32"))
805 {
807
808 m_GDALWriteIntDataType = GDT_Int32;
809 m_lGDALMaxCanWrite = INT32_MAX;
810 m_lGDALMinCanWrite = INT32_MIN;
811
813 m_GDALWriteFloatDataType = GDT_Int32;
814
815 return true;
816 }
817
818 if (strstr(CSLFetchNameValue(papszMetadata, "DMD_CREATIONDATATYPES"), "UInt16"))
819 {
820 m_bGDALCanWriteInt32 = false;
821
822 m_GDALWriteIntDataType = GDT_UInt16;
823 m_lGDALMaxCanWrite = UINT16_MAX;
825
827 m_GDALWriteFloatDataType = GDT_UInt16;
828
829 return true;
830 }
831
832 if (strstr(CSLFetchNameValue(papszMetadata, "DMD_CREATIONDATATYPES"), "Int16"))
833 {
834 m_bGDALCanWriteInt32 = false;
835
836 m_GDALWriteIntDataType = GDT_Int16;
837 m_lGDALMaxCanWrite = INT16_MAX;
838 m_lGDALMinCanWrite = INT16_MIN;
839
841 m_GDALWriteFloatDataType = GDT_Int16;
842
843 return true;
844 }
845
846 if (strstr(CSLFetchNameValue(papszMetadata, "DMD_CREATIONDATATYPES"), "Byte"))
847 {
848 m_bGDALCanWriteInt32 = false;
849
850 m_GDALWriteIntDataType = GDT_Byte;
851 m_lGDALMaxCanWrite = UINT8_MAX;
853
855 m_GDALWriteFloatDataType = GDT_Byte;
856
857 return true;
858 }
859
860 // This driver does not even support byte output
861 cerr << ERR << "Cannot write using raster GDAL driver '" << m_strRasterGISOutFormat << ", not even byte output is supported'. Choose another GIS raster format." << endl;
862 return false;
863}
864
865//===============================================================================================================================
867//===============================================================================================================================
869{
870 // Load the vector GDAL driver (this assumes that GDALAllRegister() has already been called)
871 GDALDriver *pDriver = GetGDALDriverManager()->GetDriverByName(m_strVectorGISOutFormat.c_str());
872
873 if (NULL == pDriver)
874 {
875 // Can't load vector GDAL driver. Incorrectly specified?
876 cerr << ERR << "Unknown vector GIS output format '" << m_strVectorGISOutFormat << "'." << endl;
877 return false;
878 }
879
880 // Get the metadata for this vector driver
881 char **papszMetadata = pDriver->GetMetadata();
882
883 // For GDAL2, need to test if this is a vector driver
884 if (!CSLFetchBoolean(papszMetadata, GDAL_DCAP_VECTOR, false))
885 {
886 // This is not a vector driver
887 cerr << ERR << "GDAL driver '" << m_strVectorGISOutFormat << "' is not a vector driver. Choose another format." << endl;
888 return false;
889 }
890
891 if (!CSLFetchBoolean(papszMetadata, GDAL_DCAP_CREATE, false))
892 {
893 // Driver does not support create() method
894 cerr << ERR << "Cannot write vector GIS files using GDAL driver '" << m_strRasterGISOutFormat << "'. Choose another format." << endl;
895 return false;
896 }
897
898 // Driver is OK, now set some options for individual drivers
899 if (m_strVectorGISOutFormat == "ESRI Shapefile")
900 {
901 // Set this, so that just a single dataset-with-one-layer shapefile is created, rather than a directory (see http://www.gdal.org/ogr/drv_shapefile.html)
903 }
904
905 else if (m_strVectorGISOutFormat == "geojson")
906 {
908 }
909
910 else if (m_strVectorGISOutFormat == "gpkg")
911 {
913 }
914
915 // TODO 033 Others
916
917 return true;
918}
919
920//===============================================================================================================================
922//===============================================================================================================================
924{
925 // Increment file number
926 m_nGISSave++;
927
928 // Set for next save
929 if (m_bSaveRegular)
930 {
932 }
933 else
934 {
935 if (m_nThisSave < m_nUSave - 1)
936 {
937 // Still have user-defined save times remaining
938 m_nThisSave++;
939 }
940 else
941 {
942 // Finished user-defined times, switch to regular interval using last value as interval
943 double dLastInterval;
944
945 if (m_nUSave > 1)
946 dLastInterval = m_dUSaveTime[m_nUSave - 1] - m_dUSaveTime[m_nUSave - 2];
947 else
948 dLastInterval = m_dUSaveTime[m_nUSave - 1];
949
950 m_dRegularSaveTime = m_dSimElapsed + dLastInterval;
951 m_dRegularSaveInterval = dLastInterval;
952 m_bSaveRegular = true;
953 }
954 }
955
958 return false;
959
962 return false;
963
964 if (m_bTalusSave)
966 return false;
967
970 return false;
971
974 return false;
975
976 if (m_bCliffToeSave)
978 return false;
979
980 if (m_bSeaDepthSave)
982 return false;
983
986 return false;
987
990 return false;
991
992 // Don't write platform erosion files if there is no platform erosion
994 {
997 return false;
998
1001 return false;
1002
1005 return false;
1006
1009 return false;
1010
1013 return false;
1014
1017 return false;
1018
1021 return false;
1022
1025 return false;
1026
1029 return false;
1030
1033 return false;
1034
1037 return false;
1038
1041 return false;
1042 }
1043
1044 if (m_bLandformSave)
1046 return false;
1047
1050 return false;
1051
1054 return false;
1055
1058 return false;
1059
1062 return false;
1063
1064 // Don't write suspended sediment files if there is no fine sediment
1066 {
1067 if (m_bSuspSedSave)
1069 return false;
1070
1073 return false;
1074 }
1075
1078 return false;
1079
1080 for (int nLayer = 0; nLayer < m_nLayers; nLayer++)
1081 {
1083 {
1085 return false;
1086 }
1087
1089 {
1091 return false;
1092 }
1093
1095 {
1097 return false;
1098 }
1099
1101 {
1103 return false;
1104 }
1105
1107 {
1109 return false;
1110 }
1111
1113 {
1115 return false;
1116 }
1117 }
1118
1119 if (m_bSliceSave)
1120 {
1121 for (int i = 0; i < static_cast<int>(m_VdSliceElev.size()); i++)
1122 {
1124 return false;
1125 }
1126 }
1127
1129 {
1131 return false;
1132 }
1133
1135 {
1137 return false;
1138 }
1139
1141 {
1143 return false;
1144 }
1145
1146 // Don't write cliff collapse files if we aren't considering cliff collapse
1148 {
1150 {
1152 {
1154 return false;
1155 }
1156
1158 {
1160 return false;
1161 }
1162
1164 {
1166 return false;
1167 }
1168
1170 {
1172 return false;
1173 }
1174
1175#ifdef _DEBUG
1177 {
1178 if (! bWriteRasterGISFile(RASTER_PLOT_CLIFF_COLLAPSE_TIMESTEP, &RASTER_PLOT_CLIFF_COLLAPSE_TIMESTEP_TITLE))
1179 return false;
1180 }
1181#endif
1182 }
1183
1185 {
1187 {
1189 return false;
1190 }
1191
1193 {
1195 return false;
1196 }
1197
1199 {
1201 return false;
1202 }
1203 }
1204
1206 {
1208 {
1210 return false;
1211 }
1212
1214 {
1216 return false;
1217 }
1218 }
1219
1221 {
1223 {
1225 return false;
1226 }
1227
1229 {
1231 return false;
1232 }
1233 }
1234 }
1235
1237 {
1239 return false;
1240 }
1241
1242 if (m_bSeaMaskSave)
1243 {
1245 return false;
1246 }
1247
1248 if (m_bBeachMaskSave)
1249 {
1251 return false;
1252 }
1253
1255 {
1257 return false;
1258 }
1259
1261 {
1263 return false;
1264 }
1265
1267 {
1269 return false;
1270
1272 return false;
1273 }
1274
1276 {
1278 return false;
1279 }
1280
1282 {
1284 return false;
1285 }
1286
1288 {
1290 return false;
1291 }
1292
1294 {
1296 return false;
1297 }
1298
1300 {
1302 return false;
1303 }
1304
1305 // if (m_bSetupSurgeFloodMaskSave)
1306 // {
1307 // if (! bWriteRasterGISFile(RASTER_PLOT_SETUP_SURGE_FLOOD_MASK, &RASTER_PLOT_SETUP_SURGE_FLOOD_MASK_TITLE))
1308 // return false;
1309 // }
1310
1311 // if (m_bSetupSurgeRunupFloodMaskSave)
1312 // {
1313 // if (! bWriteRasterGISFile(RASTER_PLOT_SETUP_SURGE_RUNUP_FLOOD_MASK, &RASTER_PLOT_SETUP_SURGE_RUNUP_FLOOD_MASK_TITLE))
1314 // return false;
1315 // }
1316 //
1317 return true;
1318}
1319
1320//===============================================================================================================================
1322//===============================================================================================================================
1324{
1325 // Always written
1326 if (m_bCoastSave)
1327 {
1329 return false;
1330
1332 return false;
1333
1335 return false;
1336 }
1337
1338 if (m_bCliffEdgeSave)
1339 {
1341 return false;
1342 }
1343
1344 if (m_bNormalsSave)
1345 {
1347 return false;
1348 }
1349
1351 {
1353 return false;
1354 }
1355
1357 {
1359 return false;
1360 }
1361
1363 {
1365 return false;
1366 }
1367
1369 {
1371 return false;
1372 }
1373
1375 {
1377 return false;
1378 }
1379
1381 {
1383 return false;
1384 }
1385
1387 {
1389 return false;
1390 }
1391
1393 {
1395 return false;
1396 }
1397
1399 {
1401 return false;
1402 }
1403
1405 {
1407 return false;
1408 }
1409
1411 {
1413 return false;
1414 }
1415
1417 {
1419 return false;
1420 }
1421
1423 {
1425 return false;
1426 }
1427
1429 {
1431 return false;
1432 }
1433
1434 // if (m_bWaveSetupSave)
1435 // {
1436 // if (! bWriteVectorGISFile(VECTOR_PLOT_WAVE_SETUP, &VECTOR_PLOT_WAVE_SETUP_TITLE))
1437 // return false;
1438 // }
1439 //
1440 // if (m_bStormSurgeSave)
1441 // {
1442 // if (! bWriteVectorGISFile(VECTOR_PLOT_STORM_SURGE, &VECTOR_PLOT_STORM_SURGE_TITLE))
1443 // return false;
1444 // }
1445 //
1446 // if (m_bRunUpSave)
1447 // {
1448 // if (! bWriteVectorGISFile(VECTOR_PLOT_RUN_UP, &VECTOR_PLOT_RUN_UP_TITLE))
1449 // return false;
1450 // }
1451 //
1452 // if (m_bRiverineFlooding && m_bVectorWaveFloodLineSave)
1453 // {
1454 // if (! bWriteVectorGISFile(VECTOR_PLOT_FLOOD_LINE, &VECTOR_PLOT_FLOOD_SWL_SETUP_LINE_TITLE))
1455 // return false;
1456 //
1457 // // if (! bWriteVectorGISFile(VECTOR_PLOT_FLOOD_SWL_SETUP_SURGE_LINE,
1458 // // &VECTOR_PLOT_FLOOD_SWL_SETUP_SURGE_LINE_TITLE)) return false;
1459 //
1460 // // if (! bWriteVectorGISFile(VECTOR_PLOT_FLOOD_SWL_SETUP_SURGE_RUNUP_LINE,
1461 // // &VECTOR_PLOT_FLOOD_SWL_SETUP_SURGE_RUNUP_LINE_TITLE)) return false;
1462 // }
1463
1464 return true;
1465}
1466
1467//===============================================================================================================================
1469//===============================================================================================================================
1470void CSimulation::GetRasterOutputMinMax(int const nDataItem, double& dMin, double& dMax, int const nLayer, double const dElev)
1471{
1472 // If this is a binary mask layer, we already know the max and min values
1474 (nDataItem == RASTER_PLOT_INUNDATION_MASK) ||
1475 (nDataItem == RASTER_PLOT_BEACH_MASK) ||
1476 (nDataItem == RASTER_PLOT_COAST) ||
1477 (nDataItem == RASTER_PLOT_NORMAL_PROFILE) ||
1478 (nDataItem == RASTER_PLOT_ACTIVE_ZONE) ||
1480 (nDataItem == RASTER_PLOT_SETUP_SURGE_FLOOD_MASK) ||
1482 (nDataItem == RASTER_PLOT_WAVE_FLOOD_LINE))
1483 {
1484 dMin = 0;
1485 dMax = 1;
1486
1487 return;
1488 }
1489
1490 // Not a binary mask layer, so we must find the max and min values
1491 dMin = DBL_MAX;
1492 dMax = DBL_MIN;
1493
1494 double dTmp = 0;
1495
1496 for (int nY = 0; nY < m_nYGridSize; nY++)
1497 {
1498 for (int nX = 0; nX < m_nXGridSize; nX++)
1499 {
1500 switch (nDataItem)
1501 {
1502 case (RASTER_PLOT_SLICE):
1503 dTmp = m_pRasterGrid->m_Cell[nX][nY].nGetLayerAtElev(dElev);
1504 break;
1505
1506 case (RASTER_PLOT_LANDFORM):
1507 dTmp = m_pRasterGrid->m_Cell[nX][nY].pGetLandform()->nGetLFCategory();
1508 break;
1509
1511 dTmp = INT_NODATA;
1512
1513 if (bIsInterventionCell(nX, nY))
1514 dTmp = m_pRasterGrid->m_Cell[nX][nY].pGetLandform()->nGetLFCategory();
1515
1516 break;
1517
1519 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetInterventionHeight();
1520 break;
1521
1522 case (RASTER_PLOT_POLYGON):
1523 dTmp = m_pRasterGrid->m_Cell[nX][nY].nGetPolygonID();
1524 break;
1525
1527 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetBasementElev();
1528 break;
1529
1531 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetAllSedTopElevIncTalus();
1532 break;
1533
1534 case (RASTER_PLOT_TALUS):
1535 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTalusDepth();
1536 break;
1537
1539 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTopElevIncSea();
1540 break;
1541
1543 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetConsSedSlope();
1544 break;
1545
1546 case (RASTER_PLOT_SEA_DEPTH):
1547 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetSeaDepth();
1548 break;
1549
1551 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotSeaDepth() / static_cast<double>(m_ulIter);
1552 break;
1553
1555 if (! m_pRasterGrid->m_Cell[nX][nY].bIsInContiguousSea())
1556 dTmp = m_dMissingValue;
1557 else
1558 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetWaveHeight();
1559
1560 break;
1561
1563 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotWaveHeight() / static_cast<double>(m_ulIter);
1564 break;
1565
1567 if (! m_pRasterGrid->m_Cell[nX][nY].bIsInContiguousSea())
1568 dTmp = m_dMissingValue;
1569
1570 else
1571 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetWaveAngle();
1572
1573 break;
1574
1576 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotWaveAngle() / static_cast<double>(m_ulIter);
1577 break;
1578
1580 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetBeachProtectionFactor();
1581
1582 if (! bFPIsEqual(dTmp, DBL_NODATA, TOLERANCE))
1583 dTmp = 1 - dTmp; // Output the inverse, seems more intuitive
1584
1585 break;
1586
1588 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetPotentialPlatformErosion();
1589 break;
1590
1592 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetActualPlatformErosion();
1593 break;
1594
1596 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotPotentialPlatformErosion();
1597 break;
1598
1600 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotActualPlatformErosion();
1601 break;
1602
1604 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetPotentialBeachErosion();
1605 break;
1606
1608 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetActualBeachErosion();
1609 break;
1610
1612 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotPotentialBeachErosion();
1613 break;
1614
1616 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotActualBeachErosion();
1617 break;
1618
1620 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetBeachDeposition();
1621 break;
1622
1624 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotBeachDeposition();
1625 break;
1626
1628 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetSuspendedSediment();
1629 break;
1630
1632 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotSuspendedSediment() /
1633 static_cast<double>(m_ulIter);
1634 break;
1635
1637 dTmp = m_pRasterGrid->m_Cell[nX][nY].pGetLayerAboveBasement(nLayer)->pGetUnconsolidatedSediment()->dGetFineDepth();
1638 break;
1639
1641 dTmp = m_pRasterGrid->m_Cell[nX][nY].pGetLayerAboveBasement(nLayer)->pGetUnconsolidatedSediment()->dGetSandDepth();
1642 break;
1643
1645 dTmp = m_pRasterGrid->m_Cell[nX][nY].pGetLayerAboveBasement(nLayer)->pGetUnconsolidatedSediment()->dGetCoarseDepth();
1646 break;
1647
1649 dTmp = m_pRasterGrid->m_Cell[nX][nY].pGetLayerAboveBasement(nLayer)->pGetConsolidatedSediment()->dGetFineDepth();
1650 break;
1651
1653 dTmp = m_pRasterGrid->m_Cell[nX][nY].pGetLayerAboveBasement(nLayer)->pGetConsolidatedSediment()->dGetSandDepth();
1654 break;
1655
1657 dTmp = m_pRasterGrid->m_Cell[nX][nY].pGetLayerAboveBasement(nLayer)->pGetConsolidatedSediment()->dGetCoarseDepth();
1658 break;
1659
1661 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetThisIterCliffCollapseErosionFine();
1662 break;
1663
1665 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetThisIterCliffCollapseErosionSand();
1666 break;
1667
1669 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetThisIterCliffCollapseErosionCoarse();
1670 break;
1671
1673 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotCliffCollapseFine();
1674 break;
1675
1677 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotCliffCollapseSand();
1678 break;
1679
1681 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotCliffCollapseCoarse();
1682 break;
1683
1685 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetThisIterCliffCollapseSandTalusDeposition();
1686 break;
1687
1689 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetThisIterCliffCollapseCoarseTalusDeposition();
1690 break;
1691
1693 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotSandTalusDeposition();
1694 break;
1695
1697 dTmp = m_pRasterGrid->m_Cell[nX][nY].dGetTotCoarseTalusDeposition();
1698 break;
1699
1701 dTmp = m_pRasterGrid->m_Cell[nX][nY].nGetShadowZoneNumber();
1702 break;
1703
1705 dTmp = m_pRasterGrid->m_Cell[nX][nY].nGetDownDriftZoneNumber();
1706 break;
1707
1709 dTmp = m_pRasterGrid->m_Cell[nX][nY].nGetDownDriftZoneNumber();
1710 break;
1711
1713 dTmp = m_pRasterGrid->m_Cell[nX][nY].nGetDownDriftZoneNumber();
1714 break;
1715
1717 int const nPoly = m_pRasterGrid->m_Cell[nX][nY].nGetPolygonID();
1718 int const nPolyCoast = m_pRasterGrid->m_Cell[nX][nY].nGetPolygonCoastID();
1719
1720 if (nPoly == INT_NODATA)
1721 dTmp = m_dMissingValue;
1722 else
1723 dTmp = m_VCoast[nPolyCoast].pGetPolygon(nPoly)->dGetBeachDepositionAndSuspensionAllUncons();
1724
1725 break;
1726 }
1727
1728 if (! bFPIsEqual(dTmp, DBL_NODATA, TOLERANCE))
1729 {
1730 if (dTmp > dMax)
1731 dMax = dTmp;
1732
1733 if (dTmp < dMin)
1734 dMin = dTmp;
1735 }
1736 }
1737 }
1738}
1739
1740//===============================================================================================================================
1742//===============================================================================================================================
1744{
1745 string const strDriver = strToLower(&m_strRasterGISOutFormat);
1746 string const strComment = "Created by " + PROGRAM_NAME + " for " + PLATFORM + " " + strGetBuild() + " running on " + strGetComputerName();
1747
1748 // TODO 034 Do these for all commonly-used file types
1749 if (strDriver == "aaigrid")
1750 {
1751 }
1752 else if (strDriver == "bmp")
1753 {
1754 }
1755 else if (strDriver == "gtiff")
1756 {
1757 if (m_bWorldFile)
1758 m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "TFW", "YES");
1759
1760 m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "NUM_THREADS", "ALL_CPUS");
1761 m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "COMPRESS", "LZW");
1762 }
1763 else if (strDriver == "hfa")
1764 {
1765 m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "NBITS", "4");
1766 }
1767 else if (strDriver == "jpeg")
1768 {
1769 m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "COMMENT", strComment.c_str());
1770 m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "QUALITY", "95");
1771 }
1772 else if (strDriver == "png")
1773 {
1774 if (m_bWorldFile)
1775 m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "WORLDFILE", "YES");
1776
1777 // m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "TITLE", "This is the title"); m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "DESCRIPTION", "This is a description");
1778 // m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "COPYRIGHT", "This is some copyright statement");
1779 m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "COMMENT", strComment.c_str());
1780 m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "NBITS", "4");
1781 }
1782 else if (strDriver == "rst")
1783 {
1784 m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "COMMENT", strComment.c_str());
1785 }
1786 else if (strDriver == "geojson")
1787 {
1788 // m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "COMMENT", strComment.c_str());
1789 }
1790 else if (strDriver == "gpkg")
1791 {
1792 // TODO 065 Does GDAL support overwriting raster gpkg files yet?
1793 m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "OVERWRITE", "YES");
1794 m_papszGDALRasterOptions = CSLSetNameValue(m_papszGDALRasterOptions, "USE_TILE_EXTENT", "YES");
1795 }
1796 else if (strDriver == "netcdf")
1797 {
1798 }
1799}
1800
1801//===============================================================================================================================
1803//===============================================================================================================================
1804int CSimulation::nGetOppositeDirection(int const nDirection)
1805{
1806 switch (nDirection)
1807 {
1808 case NORTH:
1809 return SOUTH;
1810
1811 case NORTH_EAST:
1812 return SOUTH - WEST;
1813
1814 case EAST:
1815 return WEST;
1816
1817 case SOUTH_EAST:
1818 return NORTH - WEST;
1819
1820 case SOUTH:
1821 return NORTH;
1822
1823 case SOUTH_WEST:
1824 return NORTH - EAST;
1825
1826 case WEST:
1827 return EAST;
1828
1829 case NORTH_WEST:
1830 return SOUTH - EAST;
1831 }
1832
1833 // Should never get here
1834 return NO_DIRECTION;
1835}
1836
1837// //===============================================================================================================================
1838// //! Given two integer points, calculates the slope and intercept of the line passing through the points
1839// //===============================================================================================================================
1840// void CSimulation::GetSlopeAndInterceptFromPoints(CGeom2DIPoint const* pPti1, CGeom2DIPoint const* pPti2, double& dSlope, double& dIntercept)
1841// {
1842// int nX1 = pPti1->nGetX();
1843// int nY1 = pPti1->nGetY();
1844// int nX2 = pPti2->nGetX();
1845// int nY2 = pPti2->nGetY();
1846//
1847// double dXDiff = nX1 - nX2;
1848// double dYDiff = nY1 - nY2;
1849//
1850// if (bFPIsEqual(dXDiff, 0.0, TOLERANCE))
1851// dSlope = 0;
1852// else
1853// dSlope = dYDiff / dXDiff;
1854//
1855// dIntercept = nY1 - (dSlope * nX1);
1856// }
1857
1858//===============================================================================================================================
1860//===============================================================================================================================
1861CGeom2DIPoint CSimulation::PtiFindClosestCoastPoint(int const nX, int const nY, int& nCoastFound)
1862{
1863 unsigned int nMinSqDist = UINT_MAX;
1864 CGeom2DIPoint PtiCoastPoint;
1865
1866 // Do for every coast
1867 for (int nCoast = 0; nCoast < static_cast<int>(m_VCoast.size()); nCoast++)
1868 {
1869 for (int j = 0; j < m_VCoast[nCoast].nGetCoastlineSize(); j++)
1870 {
1871 // Get the coords of the grid cell marked as coastline for the coastal landform object
1872 int const nXCoast = m_VCoast[nCoast].pPtiGetCellMarkedAsCoastline(j)->nGetX();
1873 int const nYCoast = m_VCoast[nCoast].pPtiGetCellMarkedAsCoastline(j)->nGetY();
1874
1875 // Calculate the squared distance between this point and the given point
1876 int const nXDist = nX - nXCoast;
1877 int const nYDist = nY - nYCoast;
1878
1879 unsigned int const nSqDist = (nXDist * nXDist) + (nYDist * nYDist);
1880
1881 // Is this the closest so dar?
1882 if (nSqDist < nMinSqDist)
1883 {
1884 nMinSqDist = nSqDist;
1885 PtiCoastPoint.SetXY(nXCoast, nYCoast);
1886 nCoastFound = nCoast;
1887 }
1888 }
1889 }
1890
1891 return PtiCoastPoint;
1892}
1893
1894//===============================================================================================================================
1896//===============================================================================================================================
1897int CSimulation::nFindClosestCoastPoint(int const nX, int const nY, int& nCoastFound)
1898{
1899 unsigned int nMinSqDist = UINT_MAX;
1900 int nCoastPoint = INT_NODATA;
1901
1902 // Do for every coast
1903 for (int nCoast = 0; nCoast < static_cast<int>(m_VCoast.size()); nCoast++)
1904 {
1905 for (int j = 0; j < m_VCoast[nCoast].nGetCoastlineSize(); j++)
1906 {
1907 // Get the coords of the grid cell marked as coastline for the coastal landform object
1908 int const nXCoast = m_VCoast[nCoast].pPtiGetCellMarkedAsCoastline(j)->nGetX();
1909 int const nYCoast = m_VCoast[nCoast].pPtiGetCellMarkedAsCoastline(j)->nGetY();
1910
1911 // Calculate the squared distance between this point and the given point
1912 int const nXDist = nX - nXCoast;
1913 int const nYDist = nY - nYCoast;
1914
1915 unsigned int const nSqDist = (nXDist * nXDist) + (nYDist * nYDist);
1916
1917 // Is this the closest so dar?
1918 if (nSqDist < nMinSqDist)
1919 {
1920 nMinSqDist = nSqDist;
1921 nCoastPoint = j;
1922 nCoastFound = nCoast;
1923 }
1924 }
1925 }
1926
1927 return nCoastPoint;
1928}
1929
1930//===============================================================================================================================
1932//===============================================================================================================================
1933int CSimulation::nConvertMetresToNumCells(double const dLen) const
1934{
1935 return nRound(dLen / m_dCellSide);
1936}
1937
1938//===============================================================================================================================
1940//===============================================================================================================================
1941void CSimulation::FindClosestPointOnStraightLine(double const dAx, double const dAy, double const dBx, double const dBy, double const dPx, double const dPy, double& dXRet, double& dYRet)
1942{
1943 double const dAPx = dPx - dAx;
1944 double const dAPy = dPy - dAy;
1945 double const dABx = dBx - dAx;
1946 double const dABy = dBy - dAy;
1947 double const dMagAB2 = (dABx * dABx) + (dABy * dABy);
1948 double const dABdotAP = (dABx * dAPx) + (dABy * dAPy);
1949 double const dT = dABdotAP / dMagAB2;
1950
1951 if ( dT < 0)
1952 {
1953 dXRet = dAx;
1954 dYRet = dAy;
1955 }
1956 else if (dT > 1)
1957 {
1958 dXRet = dBx;
1959 dYRet = dBy;
1960 }
1961 else
1962 {
1963 dXRet = dAx + (dABx * dT);
1964 dYRet = dAy + (dABy * dT);
1965 }
1966}
1967
1968//===============================================================================================================================
1970//===============================================================================================================================
1972{
1973 int const nX1 = pPt1->nGetX();
1974 int const nY1 = pPt1->nGetY();
1975 int const nX2 = pPt2->nGetX();
1976 int const nY2 = pPt2->nGetY();
1977
1978 if (nX1 == nX2)
1979 {
1980 if ((nY1 == nY2 + 1) || (nY1 == nY2 - 1))
1981 return true;
1982 }
1983
1984 if (nY1 == nY2)
1985 {
1986 if ((nX1 == nX2 + 1) || (nX1 == nX2 - 1))
1987 return true;
1988 }
1989
1990 return false;
1991}
Contains CGeom2DPoint definitions.
Contains CGeom2DIPoint definitions.
Geometry class used to represent 2D point objects with integer coordinates.
Definition 2di_point.h:25
void SetY(int const)
The integer parameter sets a value for the CGeom2DIPoint object's Y coordinate.
Definition 2di_point.cpp:75
int nGetY(void) const
Returns the CGeom2DIPoint object's integer Y coordinate.
Definition 2di_point.cpp:51
void SetXY(int const, int const)
The two integer parameters set values for the CGeom2DIPoint object's X and Y coordinates.
Definition 2di_point.cpp:81
void DivXDivY(double const, double const)
Divides the CGeom2DIPoint object's X coordinate by the first double parameter (rounded),...
void SetX(int const)
The integer parameter sets a value for the CGeom2DIPoint object's X coordinate.
Definition 2di_point.cpp:69
int * pnGetY()
Returns a reference to the CGeom2DIPoint object's integer Y coordinate.
Definition 2di_point.cpp:63
void AddXAddY(int const, int const)
The parameter is a pointer to a CGeom2DIPoint object, this is used to set values for the CGeom2DIPoin...
Definition 2di_point.cpp:95
int * pnGetX()
Returns a reference to the CGeom2DIPoint object's integer X coordinate.
Definition 2di_point.cpp:57
int nGetX(void) const
Returns the CGeom2DIPoint object's integer X coordinate.
Definition 2di_point.cpp:45
Geometry class used to represent 2D point objects with floating-point coordinates.
Definition 2d_point.h:25
double dGetY(void) const
Returns the CGeom2DPoint object's double Y coordinate.
Definition 2d_point.cpp:51
double dGetX(void) const
Returns the CGeom2DPoint object's double X coordinate.
Definition 2d_point.cpp:45
bool m_bCliffCollapseSave
Save cliff collapse raster GIS files?
Definition simulation.h:216
bool m_bAvgSeaDepthSave
Save average sea depth raster GIS files?
Definition simulation.h:105
bool m_bDeepWaterWaveAngleSave
Save deep water wave angle raster GIS files?
Definition simulation.h:243
bool bWriteRasterGISFile(int const, string const *, int const =0, double const =0)
Writes GIS raster files using GDAL, using data from the RasterGrid array.
bool m_bSlopeConsSedSave
Save slope of consolidated sediment raster GIS files?
Definition simulation.h:174
string m_strOGRVectorOutputExtension
GDAL-OGR vector output drive file extension.
bool bCheckRasterGISOutputFormat(void)
Checks whether the selected raster GDAL driver supports file creation, 32-bit doubles,...
bool m_bFineUnconsSedSave
Save fine unconsolidated sediment raster GIS files?
Definition simulation.h:189
void GetRasterOutputMinMax(int const, double &, double &, int const, double const)
Finds the max and min values in order to scale raster output if we cannot write doubles.
bool m_bCoastSave
Save coastline as vector GIS file?
Definition simulation.h:261
CGeomRasterGrid * m_pRasterGrid
Pointer to the raster grid object.
int m_nXGridSize
The size of the grid in the x direction.
Definition simulation.h:477
static int nGetOppositeDirection(int const)
Returns the opposite direction.
CGeom2DIPoint PtiExtCRSToGridRound(CGeom2DPoint const *) const
Transforms a pointer to a CGeom2DPoint in the external CRS to the equivalent CGeom2DIPoint in the ras...
bool m_bTopSurfIncSeaSave
Save top surface (sediment, talus, and sea) raster DEMs?
Definition simulation.h:87
bool m_bSedIncTalusTopSurfSave
Save sediment (inc talus) top surface raster DEMs?
Definition simulation.h:84
vector< CRWCoast > m_VCoast
The coastline objects.
bool bSaveAllVectorGISFiles(void)
The bSaveAllvectorGISFiles member function saves the vector GIS files TODO 081 Choose more files to o...
int m_nThisSave
Used in calculations of GIS save intervals.
Definition simulation.h:525
static string strGetComputerName(void)
Returns a string, hopefully giving the name of the computer on which the simulation is running.
Definition utils.cpp:1474
long m_lGDALMaxCanWrite
The maximum integer value which GDAL can write, can be UINT8_MAX, INT16_MAX, UINT16_MAX,...
Definition simulation.h:603
bool m_bAvgWaveAngleAndHeightSave
Save average wave angle and average wave height raster GIS files?
Definition simulation.h:123
bool m_bAvgWaveAngleSave
Save average wave angle raster GIS files?
Definition simulation.h:117
bool m_bCliffEdgeSave
Save cliff edge vector GIS files?
Definition simulation.h:264
bool bWriteVectorGISFile(int const, string const *)
Writes vector GIS files using GDAL/OGR.
bool m_bInvalidNormalsSave
Save invalid coastline-normal vector GIS files?
Definition simulation.h:270
bool m_bShadowBoundarySave
Save wave shadow boundary vector GIS files?
Definition simulation.h:288
bool m_bActualBeachErosionSave
Save actual (supply-limited) beach (unconsolidated sediment) erosion raster GIS files?
Definition simulation.h:156
void FindClosestPointOnStraightLine(double const, double const, double const, double const, double const, double const, double &, double &)
For the straight line between two points A and B, and given another point C, this finds the closest p...
int m_nYGridSize
The size of the grid in the y direction.
Definition simulation.h:480
bool m_bWaveTransectPointsSave
Save wave transect points vector GIS files?
Definition simulation.h:285
GDALDataType m_GDALWriteIntDataType
The data type used by GDAL for integer operations, can be GDT_Byte, GDT_Int16, GDT_UInt16,...
Definition simulation.h:597
bool m_bCliffCollapseDepositionSave
Save cliff collapse deposition raster GIS files?
Definition simulation.h:222
bool m_bTotalActualPlatformErosionSave
Save total actual (supply-limited) shore platform erosion raster GIS files?
Definition simulation.h:150
bool m_bPolygonUnconsSedUpOrDownDriftSave
Save polygon unconsolidated sediment up- or down-drift raster GIS files?
Definition simulation.h:252
double m_dGeoTransform[6]
GDAL geotransformation info (see http://www.gdal.org/classGDALDataset.html)
Definition simulation.h:711
static CGeom2DIPoint PtiGetPerpendicular(CGeom2DIPoint const *, CGeom2DIPoint const *, double const, int const)
Returns a CGeom2DIPoint (grid CRS) which is the 'other' point of a two-point vector passing through P...
string m_strVectorGISOutFormat
Vector GIS output format.
double m_dMissingValue
Used by CoastalME for floating-point missing values.
Definition simulation.h:981
static double dGetDistanceBetween(CGeom2DPoint const *, CGeom2DPoint const *)
Returns the distance (in external CRS) between two points.
bool m_bBasementElevSave
Save basement raster DEMs?
Definition simulation.h:81
static CGeom2DPoint PtAverage(CGeom2DPoint const *, CGeom2DPoint const *)
Returns a point (external CRS) which is the average of (i.e. is midway between) two other external CR...
bool m_bCoarseUnconsSedSave
Save coarse unconsolidated sediment raster GIS files?
Definition simulation.h:195
static double dKeepWithin360(double const)
Constrains the supplied angle to be within 0 and 360 degrees.
bool m_bPotentialPlatformErosionMaskSave
Save potential platform erosion mask raster GIS files?
Definition simulation.h:231
void KeepWithinValidGrid(int &, int &) const
Constrains the supplied point (in the grid CRS) to be a valid cell within the raster grid.
bool m_bWaveHeightSave
Save wave height raster GIS files?
Definition simulation.h:108
bool m_bGDALCanCreate
Is the selected GDAL output file format capable of writing files?
Definition simulation.h:372
bool m_bLandformSave
Save coast landform raster GIS files?
Definition simulation.h:171
static string strTrim(string const *)
Trims whitespace from both sides of a string, does not change the original string.
Definition utils.cpp:2514
string m_strRasterGISOutFormat
Raster GIS output format.
bool m_bTotalBeachDepositionSave
Save total beach (unconsolidated sediment) deposition raster GIS files?
Definition simulation.h:168
bool m_bSandUnconsSedSave
Save sand unconsolidated sediment raster GIS files?
Definition simulation.h:192
bool bIsAdjacentEdgeCell(CGeom2DIPoint const *, CGeom2DIPoint const *)
Given two different edge cells, returns true if they are adjacent.
bool m_bTotCliffCollapseSave
Save total cliff collapse raster GIS files?
Definition simulation.h:219
bool m_bDoShorePlatformErosion
Simulate shore platform erosion?
Definition simulation.h:363
bool m_bSliceSave
Save slices?
Definition simulation.h:99
bool m_bRasterPolygonSave
Save raster polygon raster GIS files?
Definition simulation.h:228
bool m_bBeachDepositionSave
Save beach (unconsolidated sediment) deposition raster GIS files?
Definition simulation.h:165
bool m_bSaveRegular
Save GIS files at regular intervals?
Definition simulation.h:258
int m_nLayers
The number of sediment layers.
Definition simulation.h:483
bool m_bSedimentInput
Do we have sediment input events?
Definition simulation.h:393
bool m_bNormalsSave
Save coastline-normal vector GIS files?
Definition simulation.h:267
bool m_bAvgWaveHeightSave
Save wave height raster GIS files?
Definition simulation.h:111
static string strToLower(string const *)
Returns the lower case version of an string, leaving the original unchanged.
Definition utils.cpp:2539
string m_strGDALBasementDEMDriverCode
GDAL code for the basement DEM raster file type.
bool m_bHaveSandSediment
Does this simulation consider sand-sized sediment?
Definition simulation.h:75
int m_nUSave
If user-defined GIS save intervals, the number of these.
Definition simulation.h:522
bool m_bSeaDepthSave
Save sea depth raster GIS files?
Definition simulation.h:102
bool m_bWaveAngleAndHeightSave
Save wave angle and wave height raster GIS files?
Definition simulation.h:120
bool m_bWorldFile
Write a GIS World file?
Definition simulation.h:384
bool bCheckVectorGISOutputFormat(void)
Checks whether the selected vector GDAL/OGR driver supports file creation etc.
bool m_bSlopeSaveForCliffToe
Save slope raster grids (used for cliff toe location)?
Definition simulation.h:93
bool m_bCliffCollapseTimestepSave
Are we saving the timestep at which each cliff occurred?
Definition simulation.h:447
bool m_bRasterNormalProfileSave
Save rasterized coastline-normal profiles GIS files?
Definition simulation.h:210
bool m_bActiveZoneSave
Save active zone raster GIS files?
Definition simulation.h:213
double dGridCentroidYToExtCRSY(int const) const
Given the integer Y-axis ordinate of a cell in the raster grid CRS, returns the external CRS Y-axis o...
Definition gis_utils.cpp:74
bool m_bDeepWaterWaveHeightSave
Save deep water wave height raster GIS files?
Definition simulation.h:246
bool m_bMeanWaveEnergySave
Save mean wave energy raster GIS files?
Definition simulation.h:132
bool m_bCliffToeSave
Save cliff toe raster grids?
Definition simulation.h:96
double m_dSimElapsed
Time simulated so far, in hours.
Definition simulation.h:693
bool m_bCoastCurvatureSave
Save coastline-curvature vector GIS files?
Definition simulation.h:273
int nConvertMetresToNumCells(double const) const
Given a length in m, this returns the rounded equivalent number of cells.
double dGridYToExtCRSY(double const) const
Given a real-valued Y-axis ordinate in the raster grid CRS (i.e. not the centroid of a cell),...
bool m_bGDALCanWriteInt32
Is the selected GDAL output file format capable of writing 32-bit integers to files?
Definition simulation.h:378
bool m_bBreakingWaveHeightSave
Save breaking wave height raster GIS files?
Definition simulation.h:135
double m_dRegularSaveTime
The time of the next save, in hours from the start of the simulation, if we are saving regularly.
Definition simulation.h:696
bool m_bPolygonNodeSave
Save polygon node vector GIS files?
Definition simulation.h:276
bool m_bTotalPotentialPlatformErosionSave
Save total potential shore platform erosion raster GIS files?
Definition simulation.h:147
void SetRasterFileCreationDefaults(void)
Sets per-driver defaults for raster files created using GDAL.
bool m_bWaveEnergySinceCollapseSave
Save wave energy since cliff collapse raster GIS files?
Definition simulation.h:129
bool m_bPotentialBeachErosionSave
Save potential beach (unconsolidated sediment) erosion raster GIS files?
Definition simulation.h:153
bool m_bSuspSedSave
Save suspended sediment raster GIS files?
Definition simulation.h:183
static double dAngleSubtended(CGeom2DIPoint const *, CGeom2DIPoint const *, CGeom2DIPoint const *)
Returns the signed angle BAC (in radians) subtended between three CGeom2DIPoints B A C....
bool m_bPolygonBoundarySave
Save polygon boundary vector GIS files?
Definition simulation.h:279
double dExtCRSXToGridX(double const) const
Transforms an X-axis ordinate in the external CRS to the equivalent X-axis ordinate in the raster gri...
bool m_bActualPlatformErosionSave
Save actual (supply-limited) shore platform erosion raster GIS files?
Definition simulation.h:144
bool bSaveAllRasterGISFiles(void)
The bSaveAllRasterGISFiles member function saves the raster GIS files using values from the RasterGri...
bool m_bHaveFineSediment
Does this simulation consider fine-sized sediment?
Definition simulation.h:72
static string strGetBuild(void)
Returns the date and time on which the program was compiled.
Definition utils.cpp:1744
bool m_bTotalPotentialBeachErosionSave
Save total potential beach (unconsolidated sediment) erosion raster GIS files?
Definition simulation.h:159
int m_nGISSave
The save number for GIS files (can be sequential, or the iteration number)
Definition simulation.h:519
static CGeom2DIPoint PtiWeightedAverage(CGeom2DIPoint const *, CGeom2DIPoint const *, double const)
Returns an integer point (grid CRS) which is the weighted average of two other grid CRS integer point...
bool m_bSeaMaskSave
Save sea mask raster GIS files?
Definition simulation.h:234
bool m_bInterventionClassSave
Save intervention class raster GIS files?
Definition simulation.h:177
bool m_bTotalActualBeachErosionSave
Save total actual (supply-limited) beach (unconsolidated sediment) erosion raster GIS files?
Definition simulation.h:162
bool m_bRasterCoastlineSave
Save rasterized coastline GIS files?
Definition simulation.h:207
bool m_bInterventionHeightSave
Save intervention height raster GIS files?
Definition simulation.h:180
long m_lGDALMinCanWrite
The minimum integer value which GDAL can write, can be zero, INT16_MIN, INT32_MIN.
Definition simulation.h:606
bool m_bSandConsSedSave
Save sand consolidated sediment raster GIS files?
Definition simulation.h:201
bool m_bSedimentInputEventSave
Save sediment inut data?
Definition simulation.h:405
static double dTriangleAreax2(CGeom2DPoint const *, CGeom2DPoint const *, CGeom2DPoint const *)
Returns twice the signed area of a triangle, defined by three points.
bool m_bHaveCoarseSediment
Does this simulation consider coarse-sized sediment?
Definition simulation.h:78
double m_dRegularSaveInterval
The interval between regular saves, in hours.
Definition simulation.h:699
bool m_bPolygonUnconsSedGainOrLossSave
Save polygon unconsolidated sediment gain or loss raster GIS files?
Definition simulation.h:255
string m_strGDALRasterOutputDriverLongname
GDAL raster output driver long name.
bool m_bDeepWaterWaveAngleAndHeightSave
Save deep water wave angle and wave height raster GIS files?
Definition simulation.h:126
double dExtCRSYToGridY(double const) const
Transforms a Y-axis ordinate in the external CRS to the equivalent Y-axis ordinate in the raster grid...
vector< double > m_VdSliceElev
Elevations for raster slice output.
bool m_bBeachProtectionSave
Save beach protection raster GIS files>
Definition simulation.h:138
bool m_bFineConsSedSave
Save fine consolidated sediment raster GIS files?
Definition simulation.h:198
bool m_bShadowDowndriftBoundarySave
Save wave shadow downdrift boundary vector GIS files?
Definition simulation.h:291
bool bIsWithinValidGrid(int const, int const) const
Checks whether the supplied point (an x-y pair, in the grid CRS) is within the raster grid,...
bool m_bDeepWaterWavePeriodSave
Save deep water wave period raster GIS files?
Definition simulation.h:249
int nFindClosestCoastPoint(int const, int const, int &)
Finds the number of the closest point on any coastline to a given point, or INT_NODATA in case of err...
double dGridXToExtCRSX(double const) const
Given a real-valued X-axis ordinate in the raster grid CRS (i.e. not the centroid of a cell),...
Definition gis_utils.cpp:95
bool m_bCoarseConsSedSave
Save coarse consolidated sediment raster GIS files?
Definition simulation.h:204
CGeom2DPoint PtGridCentroidToExt(CGeom2DIPoint const *) const
Transforms a pointer to a CGeom2DIPoint in the raster grid CRS (assumed to be the centroid of a cell)...
Definition gis_utils.cpp:83
string m_strGDALRasterOutputDriverExtension
GDAL raster output driver file extension.
bool m_bBeachMaskSave
Save beach mask raster GIS files?
Definition simulation.h:237
unsigned long m_ulIter
The number of the current iteration (time step)
Definition simulation.h:609
bool m_bAvgSuspSedSave
Save average suspended sediment raster GIS files?
Definition simulation.h:186
double dGridCentroidXToExtCRSX(int const) const
Given the integer X-axis ordinate of a cell in the raster grid CRS, returns the external CRS X-axis o...
Definition gis_utils.cpp:65
bool m_bTalusSave
Save talus depth?
Definition simulation.h:90
bool m_bCliffNotchAllSave
Are we saving all cliff notches?
Definition simulation.h:444
double m_dCellSide
Length of a cell side (in external CRS units)
Definition simulation.h:672
bool bIsInterventionCell(int const, int const) const
Returns true if the cell is an intervention.
Definition utils.cpp:3066
bool m_bTotCliffCollapseDepositionSave
Save total cliff collapse deposition raster GIS files?
Definition simulation.h:225
double m_dUSaveTime[SAVEMAX]
Save time, in hours from the start of the simukation, if we are not saving regularly.
Definition simulation.h:702
bool m_bDoCliffCollapse
Simulate cliff collapse?
Definition simulation.h:366
CGeom2DIPoint PtiFindClosestCoastPoint(int const, int const, int &)
Finds the closest point on any coastline to a given point.
bool m_bShadowZoneCodesSave
Save wave shadow zones raster GIS files?
Definition simulation.h:240
bool m_bPotentialPlatformErosionSave
Save potential shore platform erosion raster GIS files?
Definition simulation.h:141
static CGeom2DIPoint PtiPolygonCentroid(vector< CGeom2DIPoint > *)
Returns an integer point (grid CRS) which is the centroid of a polygon, given by a vector of grid CRS...
GDALDataType m_GDALWriteFloatDataType
Thw data type used by GDAL for floating point operations, can be GDT_Byte, GDT_Int16,...
Definition simulation.h:600
bool m_bGDALCanWriteFloat
Is the selected GDAL output file format capable of writing floating-point values to files?
Definition simulation.h:375
char ** m_papszGDALRasterOptions
Options for GDAL when handling raster files.
Definition simulation.h:471
bool m_bCliffNotchSave
Save cliff notch incision depth vector GIS files?
Definition simulation.h:282
bool m_bWaveAngleSave
Save wave angle raster GIS files?
Definition simulation.h:114
This file contains global definitions for CoastalME.
int const INT_NODATA
Definition cme.h:380
int const RASTER_PLOT_POLYGON
Definition cme.h:522
string const RASTER_PLOT_POLYGON_UPDRIFT_OR_DOWNDRIFT_TITLE
Definition cme.h:1010
string const RASTER_PLOT_CLIFF_COLLAPSE_EROSION_COARSE_TITLE
Definition cme.h:987
string const VECTOR_PLOT_BREAKING_WAVE_HEIGHT_TITLE
Definition cme.h:1100
int const RASTER_PLOT_CONS_SED_SLOPE
Definition cme.h:539
T tMin(T a, T b)
Definition cme.h:1175
string const VECTOR_PLOT_INVALID_NORMALS_TITLE
Definition cme.h:1113
double const TOLERANCE
Definition cme.h:725
string const VECTOR_PLOT_NORMALS_TITLE
Definition cme.h:1115
string const RASTER_PLOT_CLIFF_COLLAPSE_EROSION_SAND_TITLE
Definition cme.h:989
string const RASTER_PLOT_COAST_TITLE
Definition cme.h:997
string const RASTER_PLOT_POLYGON_TITLE
Definition cme.h:1009
int const VECTOR_PLOT_BREAKING_WAVE_HEIGHT
Definition cme.h:559
int const VECTOR_PLOT_POLYGON_NODES
Definition cme.h:573
int const RASTER_PLOT_TOTAL_ACTUAL_BEACH_EROSION
Definition cme.h:542
string const RASTER_PLOT_BEACH_DEPOSITION_TITLE
Definition cme.h:982
int const RASTER_PLOT_CLIFF_COLLAPSE_DEPOSITION_SAND
Definition cme.h:499
int const RASTER_PLOT_BEACH_DEPOSITION
Definition cme.h:495
int const RASTER_PLOT_SUSPENDED_SEDIMENT
Definition cme.h:540
int const VECTOR_PLOT_NORMALS
Definition cme.h:571
string const RASTER_PLOT_DEEP_WATER_WAVE_ORIENTATION_TITLE
Definition cme.h:999
string const RASTER_PLOT_FINE_CONSOLIDATED_SEDIMENT_TITLE
Definition cme.h:1001
int const RASTER_PLOT_FINE_UNCONSOLIDATED_SEDIMENT
Definition cme.h:515
string const RASTER_PLOT_AVG_SEA_DEPTH_TITLE
Definition cme.h:977
string const ERR
Definition cme.h:805
string const RASTER_PLOT_BEACH_MASK_TITLE
Definition cme.h:983
int const RASTER_PLOT_INUNDATION_MASK
Definition cme.h:518
string const RASTER_PLOT_AVG_WAVE_ORIENTATION_TITLE
Definition cme.h:980
string const RASTER_PLOT_ACTUAL_BEACH_EROSION_TITLE
Definition cme.h:975
string const RASTER_PLOT_TOTAL_CLIFF_COLLAPSE_EROSION_FINE_TITLE
Definition cme.h:1035
int const RASTER_PLOT_ACTUAL_BEACH_EROSION
Definition cme.h:488
string const RASTER_PLOT_POLYGON_GAIN_OR_LOSS_TITLE
Definition cme.h:1008
string const VECTOR_PLOT_MEAN_WAVE_ENERGY_TITLE
Definition cme.h:1114
string const RASTER_PLOT_DEEP_WATER_WAVE_HEIGHT_TITLE
Definition cme.h:998
string const RASTER_PLOT_WAVE_ORIENTATION_TITLE
Definition cme.h:1041
string const VECTOR_PLOT_SHADOW_ZONE_BOUNDARY_TITLE
Definition cme.h:1119
string const RASTER_PLOT_BASEMENT_ELEVATION_TITLE
Definition cme.h:981
string const RASTER_PLOT_INTERVENTION_CLASS_TITLE
Definition cme.h:1003
string const RASTER_PLOT_COARSE_UNCONSOLIDATED_SEDIMENT_TITLE
Definition cme.h:996
string const RASTER_PLOT_SLOPE_FOR_CLIFF_TOE_TITLE
Definition cme.h:1024
string const VECTOR_PLOT_POLYGON_NODES_TITLE
Definition cme.h:1117
int const RASTER_PLOT_POTENTIAL_PLATFORM_EROSION_MASK
Definition cme.h:527
string const RASTER_PLOT_AVG_WAVE_HEIGHT_TITLE
Definition cme.h:979
string const RASTER_PLOT_TOTAL_POTENTIAL_PLATFORM_EROSION_TITLE
Definition cme.h:1038
string const RASTER_PLOT_SED_TOP_INC_TALUS_ELEV_TITLE
Definition cme.h:1018
int const RASTER_PLOT_SAND_CONSOLIDATED_SEDIMENT
Definition cme.h:528
int const RASTER_PLOT_AVG_WAVE_HEIGHT
Definition cme.h:492
string const VECTOR_PLOT_CLIFF_EDGE_TITLE
Definition cme.h:1101
int const SOUTH_EAST
Definition cme.h:402
int const RASTER_PLOT_FINE_CONSOLIDATED_SEDIMENT
Definition cme.h:514
int const RASTER_PLOT_TOP_ELEV_INC_SEA
Definition cme.h:521
int const VECTOR_PLOT_COAST_SWL_HIGHEST
Definition cme.h:564
int const RASTER_PLOT_ACTIVE_ZONE
Definition cme.h:487
int const VECTOR_PLOT_COAST_CURVATURE
Definition cme.h:563
string const RASTER_PLOT_CLIFF_COLLAPSE_DEPOSITION_COARSE_TITLE
Definition cme.h:985
string const RASTER_PLOT_TOTAL_CLIFF_COLLAPSE_EROSION_SAND_TITLE
Definition cme.h:1036
string const RASTER_PLOT_LANDFORM_TITLE
Definition cme.h:1006
int const RASTER_PLOT_DEEP_WATER_WAVE_PERIOD
Definition cme.h:513
int const RASTER_PLOT_TALUS
Definition cme.h:541
int const RASTER_PLOT_TOTAL_CLIFF_COLLAPSE_DEPOSITION_SAND
Definition cme.h:546
string const RASTER_PLOT_WAVE_HEIGHT_TITLE
Definition cme.h:1040
int const RASTER_PLOT_COAST
Definition cme.h:510
string const VECTOR_PLOT_AVG_WAVE_ANGLE_AND_HEIGHT_TITLE
Definition cme.h:1099
string const RASTER_PLOT_TOTAL_POTENTIAL_BEACH_EROSION_TITLE
Definition cme.h:1037
string const VECTOR_PLOT_COAST_CURVATURE_TITLE
Definition cme.h:1103
string const RASTER_PLOT_POTENTIAL_BEACH_EROSION_TITLE
Definition cme.h:1011
bool bFPIsEqual(const T d1, const T d2, const T dEpsilon)
Definition cme.h:1213
int const RASTER_PLOT_POLYGON_UPDRIFT_OR_DOWNDRIFT
Definition cme.h:524
int const VECTOR_PLOT_SHADOW_ZONE_BOUNDARY
Definition cme.h:575
int const VECTOR_PLOT_INVALID_NORMALS
Definition cme.h:569
string const RASTER_PLOT_CLIFF_NOTCH_ALL_TITLE
Definition cme.h:993
int const RASTER_PLOT_POLYGON_GAIN_OR_LOSS
Definition cme.h:523
string const RASTER_PLOT_SAND_CONSOLIDATED_SEDIMENT_TITLE
Definition cme.h:1014
int const RASTER_PLOT_DEEP_WATER_WAVE_ORIENTATION
Definition cme.h:512
int const RASTER_PLOT_TOTAL_CLIFF_COLLAPSE_EROSION_COARSE
Definition cme.h:547
string const RASTER_PLOT_ACTUAL_PLATFORM_EROSION_TITLE
Definition cme.h:976
int const RASTER_PLOT_COARSE_UNCONSOLIDATED_SEDIMENT
Definition cme.h:509
int const RASTER_PLOT_AVG_WAVE_ORIENTATION
Definition cme.h:493
string const RASTER_PLOT_TOTAL_BEACH_DEPOSITION_TITLE
Definition cme.h:1031
int const RASTER_PLOT_WAVE_ORIENTATION
Definition cme.h:554
int const SOUTH_WEST
Definition cme.h:404
int const RASTER_PLOT_BEACH_MASK
Definition cme.h:496
int const RASTER_PLOT_WAVE_FLOOD_LINE
Definition cme.h:552
T tMax(T a, T b)
Definition cme.h:1162
int const RASTER_PLOT_SEDIMENT_INPUT
Definition cme.h:531
int const RASTER_PLOT_POTENTIAL_BEACH_EROSION
Definition cme.h:525
int const NORTH_WEST
Definition cme.h:406
int const VECTOR_PLOT_DOWNDRIFT_ZONE_BOUNDARY
Definition cme.h:567
int const RASTER_PLOT_AVG_SUSPENDED_SEDIMENT
Definition cme.h:491
int const VECTOR_PLOT_COAST
Definition cme.h:562
string const RASTER_PLOT_BEACH_PROTECTION_TITLE
Definition cme.h:984
int const RASTER_PLOT_TOTAL_POTENTIAL_PLATFORM_EROSION
Definition cme.h:551
int const RASTER_PLOT_SLOPE_FOR_CLIFF_TOE
Definition cme.h:538
string const RASTER_PLOT_POTENTIAL_PLATFORM_EROSION_MASK_TITLE
Definition cme.h:1012
int const VECTOR_PLOT_COAST_SWL_LOWEST
Definition cme.h:565
string const VECTOR_PLOT_COAST_TITLE
Definition cme.h:1106
string const PROGRAM_NAME
Definition cme.h:738
int const RASTER_PLOT_INTERVENTION_CLASS
Definition cme.h:516
string const VECTOR_PLOT_WAVE_TRANSECT_POINTS_TITLE
Definition cme.h:1124
int const VECTOR_PLOT_WAVE_ENERGY_SINCE_COLLAPSE
Definition cme.h:578
int const RASTER_PLOT_BEACH_PROTECTION
Definition cme.h:497
int const RASTER_PLOT_AVG_SEA_DEPTH
Definition cme.h:490
int const VECTOR_PLOT_WAVE_ANGLE_AND_HEIGHT
Definition cme.h:577
int const VECTOR_PLOT_CLIFF_NOTCH_ACTIVE
Definition cme.h:561
string const RASTER_PLOT_SHADOW_ZONE_TITLE
Definition cme.h:1022
int const RASTER_PLOT_TOTAL_BEACH_DEPOSITION
Definition cme.h:544
int const RASTER_PLOT_TOTAL_CLIFF_COLLAPSE_EROSION_FINE
Definition cme.h:548
int const RASTER_PLOT_TOTAL_CLIFF_COLLAPSE_EROSION_SAND
Definition cme.h:549
int const RASTER_PLOT_CLIFF_TOE
Definition cme.h:507
int const RASTER_PLOT_CLIFF_COLLAPSE_EROSION_FINE
Definition cme.h:501
int const SOUTH
Definition cme.h:403
int const NO_DIRECTION
Definition cme.h:398
int const VECTOR_PLOT_CLIFF_EDGE
Definition cme.h:560
string const RASTER_PLOT_AVG_SUSPENDED_SEDIMENT_TITLE
Definition cme.h:978
int const RASTER_PLOT_ACTUAL_PLATFORM_EROSION
Definition cme.h:489
int const VECTOR_PLOT_POLYGON_BOUNDARY
Definition cme.h:572
string const RASTER_PLOT_INTERVENTION_HEIGHT_TITLE
Definition cme.h:1004
int const VECTOR_PLOT_MEAN_WAVE_ENERGY
Definition cme.h:570
int const EAST
Definition cme.h:401
int const RASTER_PLOT_SED_TOP_INC_TALUS_ELEV
Definition cme.h:532
string const VECTOR_PLOT_POLYGON_BOUNDARY_TITLE
Definition cme.h:1116
string const RASTER_PLOT_SEA_DEPTH_TITLE
Definition cme.h:1016
string const RASTER_PLOT_TOTAL_ACTUAL_PLATFORM_EROSION_TITLE
Definition cme.h:1030
string const RASTER_PLOT_SEDIMENT_INPUT_EVENT_TITLE
Definition cme.h:1017
int const RASTER_PLOT_INTERVENTION_HEIGHT
Definition cme.h:517
int const RASTER_PLOT_TOTAL_POTENTIAL_BEACH_EROSION
Definition cme.h:550
int const NORTH_EAST
Definition cme.h:400
int const RASTER_PLOT_SAND_UNCONSOLIDATED_SEDIMENT
Definition cme.h:529
string const RASTER_PLOT_TOP_ELEV_INC_SEA_TITLE
Definition cme.h:1028
string const RASTER_PLOT_FINE_UNCONSOLIDATED_SEDIMENT_TITLE
Definition cme.h:1002
string const RASTER_PLOT_NORMAL_PROFILE_TITLE
Definition cme.h:1007
string const RASTER_PLOT_DEEP_WATER_WAVE_PERIOD_TITLE
Definition cme.h:1000
int const RASTER_PLOT_POTENTIAL_PLATFORM_EROSION
Definition cme.h:526
string const RASTER_PLOT_INUNDATION_MASK_TITLE
Definition cme.h:1005
string const RASTER_PLOT_POTENTIAL_PLATFORM_EROSION_TITLE
Definition cme.h:1013
string const RASTER_PLOT_TALUS_TITLE
Definition cme.h:1027
int const RASTER_PLOT_NORMAL_PROFILE
Definition cme.h:520
double const DBL_NODATA
Definition cme.h:736
int const RASTER_PLOT_CLIFF_NOTCH_ALL
Definition cme.h:506
int const RASTER_PLOT_CLIFF_COLLAPSE_DEPOSITION_COARSE
Definition cme.h:498
int const RASTER_PLOT_WAVE_HEIGHT
Definition cme.h:553
int const RASTER_PLOT_SHADOW_ZONE
Definition cme.h:536
int const VECTOR_PLOT_AVG_WAVE_ANGLE_AND_HEIGHT
Definition cme.h:558
int const RASTER_PLOT_DEEP_WATER_WAVE_HEIGHT
Definition cme.h:511
string const RASTER_PLOT_CLIFF_COLLAPSE_DEPOSITION_SAND_TITLE
Definition cme.h:986
string const RASTER_PLOT_CLIFF_TOE_TITLE
Definition cme.h:994
string const VECTOR_PLOT_DOWNDRIFT_ZONE_BOUNDARY_TITLE
Definition cme.h:1108
int const VECTOR_PLOT_DEEP_WATER_WAVE_ANGLE_AND_HEIGHT
Definition cme.h:566
string const RASTER_PLOT_TOTAL_CLIFF_COLLAPSE_EROSION_COARSE_TITLE
Definition cme.h:1034
string const RASTER_PLOT_SHADOW_DOWNDRIFT_ZONE_TITLE
Definition cme.h:1021
int const RASTER_PLOT_COARSE_CONSOLIDATED_SEDIMENT
Definition cme.h:508
string const RASTER_PLOT_TOTAL_CLIFF_COLLAPSE_DEPOSITION_SAND_TITLE
Definition cme.h:1033
int const RASTER_PLOT_TOTAL_ACTUAL_PLATFORM_EROSION
Definition cme.h:543
int const RASTER_PLOT_SETUP_SURGE_RUNUP_FLOOD_MASK
Definition cme.h:534
int const RASTER_PLOT_BASEMENT_ELEVATION
Definition cme.h:494
string const VECTOR_PLOT_WAVE_ENERGY_SINCE_COLLAPSE_TITLE
Definition cme.h:1122
int const VECTOR_PLOT_WAVE_TRANSECT_POINTS
Definition cme.h:580
string const VECTOR_PLOT_CLIFF_NOTCH_ACTIVE_TITLE
Definition cme.h:1102
int const RASTER_PLOT_LANDFORM
Definition cme.h:519
string const RASTER_PLOT_SUSPENDED_SEDIMENT_TITLE
Definition cme.h:1026
int const RIGHT_HANDED
Definition cme.h:413
string const RASTER_PLOT_CLIFF_COLLAPSE_EROSION_FINE_TITLE
Definition cme.h:988
int const RASTER_PLOT_CLIFF_COLLAPSE_EROSION_COARSE
Definition cme.h:500
string const RASTER_PLOT_TOTAL_CLIFF_COLLAPSE_DEPOSITION_COARSE_TITLE
Definition cme.h:1032
string const VECTOR_PLOT_COAST_SWL_LOWEST_TITLE
Definition cme.h:1105
int const RASTER_PLOT_SHADOW_DOWNDRIFT_ZONE
Definition cme.h:535
string const RASTER_PLOT_SAND_UNCONSOLIDATED_SEDIMENT_TITLE
Definition cme.h:1015
int const RASTER_PLOT_SLICE
Definition cme.h:537
string const RASTER_PLOT_ACTIVE_ZONE_TITLE
Definition cme.h:974
string const VECTOR_PLOT_WAVE_ANGLE_AND_HEIGHT_TITLE
Definition cme.h:1121
string const RASTER_PLOT_CONS_SED_SLOPE_TITLE
Definition cme.h:1025
int const RASTER_PLOT_TOTAL_CLIFF_COLLAPSE_DEPOSITION_COARSE
Definition cme.h:545
string const RASTER_PLOT_TOTAL_ACTUAL_BEACH_EROSION_TITLE
Definition cme.h:1029
string const VECTOR_PLOT_COAST_SWL_HIGHEST_TITLE
Definition cme.h:1104
int const RASTER_PLOT_CLIFF_COLLAPSE_EROSION_SAND
Definition cme.h:502
string const RASTER_PLOT_SLICE_TITLE
Definition cme.h:1023
int const RASTER_PLOT_SEA_DEPTH
Definition cme.h:530
int const NORTH
Definition cme.h:399
string const RASTER_PLOT_COARSE_CONSOLIDATED_SEDIMENT_TITLE
Definition cme.h:995
string const VECTOR_PLOT_DEEP_WATER_WAVE_ANGLE_AND_HEIGHT_TITLE
Definition cme.h:1107
int const RASTER_PLOT_SETUP_SURGE_FLOOD_MASK
Definition cme.h:533
int const WEST
Definition cme.h:405
char const SPACE
Definition cme.h:357
Contains CRWCoast definitions.
Contains CGeomRasterGrid definitions.
int nRound(double const d)
Correctly rounds doubles, returns an int.