CoastalME (Coastal Modelling Environment)
Simulates the long-term behaviour of complex coastlines
Loading...
Searching...
No Matches
simulation.h
Go to the documentation of this file.
1
14
15#ifndef SIMULATION_H
16#define SIMULATION_H
17/*===============================================================================================================================
18
19This file is part of CoastalME, the Coastal Modelling Environment.
20
21CoastalME 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.
22
23This 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.
24
25You 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.
26
27===============================================================================================================================*/
28#include <ctime>
29using std::localtime;
30using std::time;
31using std::time_t;
32
33#include <fstream>
34using std::ofstream;
35
36#include <string>
37using std::string;
38
39#include <utility>
40using std::pair;
41
42#include <stack>
43using std::stack;
44
45#include <random>
46using std::default_random_engine;
47using std::normal_distribution;
48
49#include <gdal_priv.h>
50
51#include "line.h"
52#include "i_line.h"
53
54#include "inc/cshore.h"
55
56int const NRNG = 2;
57int const SAVEMAX = 100000;
58
59class CGeomRasterGrid; // Forward declarations
60class CRWCoast;
61class CGeomProfile;
63class CRWCliff;
64class CSedInputEvent;
65class CRWCellLandform;
66
68{
69private:
72
75
78
81
84
87
90
93
96
99
102
105
108
111
114
117
120
123
126
129
132
135
138
141
144
147
150
153
156
159
162
165
168
171
174
177
180
183
186
189
192
195
198
201
204
207
210
213
216
219
222
225
228
231
234
237
240
243
246
249
252
255
258
261
264
267
270
273
276
279
282
285
288
291
294
297
300
303
306
309
312
315
318
321
324
327
330
333
336
339
342
345
348
351
354
357
360
363
366
369
372
375
378
381
384
387
390
393
396
399
402
405
408
411
414
417
420
423
426
429
432
435
438
441
444
447
450
453
456
459
462
465
468
471
474
475 // NOT USED
476 // int m_nNThisIterCliffCollapse;
477
478 // NOT USED
479 //int m_nNTotCliffCollapse;
480
483
486
489
492
495
498
501
504
507
510
513
516
519
522
525
528
531
534
537
539
542
545
548
551
553 unsigned long m_ulIter;
554
556 unsigned long m_ulTotTimestep;
557
559 unsigned long m_ulRandSeed[NRNG];
560
562 unsigned long m_ulNumCells;
563
566
569
572
575
578
581
584
587
590
593
596
599
602
605
608
611
614
617
620
623
626
629
632
635
638
641
644
647
650
653
656
659
662
665
668
671
674
676 double m_dMinSWL;
677
679 double m_dMaxSWL;
680
683
686
689
692
695
697 double m_dC_0;
698
700 double m_dL_0;
701
704
707
710
713
716
719
722
724 double m_dR;
725
728
731
734
737
740
743
746
749
752
755
758
760 double m_dKLS;
761
764
766 double m_dG;
767
770
773
776
779
782
785
788
791
794
797
800
803
806
809
812
815
818
821
824
827
830
833
836
839
842
845
848
851
854
857
860
863
866
869
872
875
878
881
884
887
890
893
896
899
902
905
908
911
914
917
920
923
926
929
932
935
938
941
944
947
950
953
956
959
962
965
968
971
974
977
978 // These grand totals are all long doubles. The aim is to minimize rounding errors when many very small numbers are added to a single much larger number, see e.g. http://www.ddj.com/cpp/184403224
979
982
985
988
991
994
997
1000
1003
1006
1009
1012
1015
1018
1021
1024
1027
1030
1033
1036
1039
1042
1045
1048
1051
1054
1057
1060
1063
1066
1069
1072
1075
1078
1081
1084
1087
1089
1092
1095
1098
1101
1104
1107
1110
1113
1116
1119
1122
1125
1128
1131
1134
1137
1140
1143
1146
1149
1152
1155
1158
1161
1164
1167
1170
1173
1176
1179
1182
1185
1188
1191
1194
1197
1198 // GDAL description for the deep water wave stations vector file
1200
1203
1206
1209
1210 // GDAL description for the sediment input event locations vector file
1212
1215
1218
1221
1222 // GDAL description for the flood input locations point or vector file
1224
1227
1230
1233
1236
1239
1242
1245
1248
1251
1254
1257
1260
1263
1266
1268 ofstream OutStream;
1269
1272
1275
1278
1281
1284
1287
1290
1293
1296
1299
1302
1305
1308
1311
1314
1317
1320
1323
1326
1328 vector<unsigned long> m_VulProfileTimestep;
1329
1332
1334 vector<double> m_VdSliceElev;
1335
1337 vector<double> m_VdErosionPotential;
1338
1340 vector<double> m_VdDepthOverDB;
1341
1343 vector<double> m_VdSavGolFCRWCoast;
1344
1347
1349 vector<double> m_VdTideData;
1350
1353
1356
1359
1362
1365
1368
1371
1374
1377
1380
1382 vector<double> m_VdFloodLocationX;
1383
1385 vector<double> m_VdFloodLocationY;
1386
1389
1392
1395
1398
1401
1404
1407
1410
1413
1416
1419
1422
1425
1428
1431
1434
1437
1440
1443
1446
1449
1452
1455
1458
1461
1464
1467
1470
1473
1476
1479
1481 vector<CRWCoast> m_VCoast;
1482
1484 vector<CRWCoast> m_VFloodWaveSetupSurge;
1485
1488
1490 vector<CGeomCoastPolygon*> m_pVCoastPolygon;
1491
1493 vector<CGeom2DIPoint> m_VEdgeCell;
1494
1496 vector<int> m_VEdgeCellEdge;
1497
1500
1502 vector<CSedInputEvent*> m_pVSedInputEvent;
1503
1505 default_random_engine m_Rand[NRNG];
1506
1508 normal_distribution<double> m_dUnitNormalDist{0.0, 1.0};
1509
1510private:
1511 // Input and output routines
1512 int nHandleCommandLineParams(int, char const* []);
1513 bool bReadIniFile(void);
1514 bool bReadRunDataFile(void);
1515 bool bOpenLogFile(void);
1516 bool bSetUpTSFiles(void);
1517 void WriteStartRunDetails(void);
1518 bool bWritePerTimestepResults(void);
1519 bool bWriteTSFiles(void);
1520 int nWriteEndRunDetails(void);
1521 int nReadShapeFunctionFile(void);
1522 int nReadWaveStationInputFile(int const);
1524 int nReadTideDataFile(void);
1525 int nSaveProfile(int const, CGeomProfile const*, int const, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<CGeom2DIPoint>* const, vector<double> const*) const;
1526 bool bWriteProfileData(int const, CGeomProfile const*, int const, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<CGeom2DIPoint>* const, vector<double> const*) const;
1527 int nSaveParProfile(int const, CGeomProfile const*, int const, int const, int const, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<CGeom2DIPoint>*const, vector<double> const*) const;
1528 bool bWriteParProfileData(int const, int const, int const, int const, int const, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*, vector<CGeom2DIPoint>*const, vector<double> const*) const;
1529 void WriteLookUpData(void) const;
1530
1531 // GIS input and output stuff
1532 int nReadRasterBasementDEM(void);
1533 int nReadRasterGISFile(int const, int const);
1534 int nReadVectorGISFile(int const);
1535 bool bWriteRasterGISFile(int const, string const*, int const = 0, double const = 0);
1536 bool bWriteVectorGISFile(int const, string const*);
1537 void GetRasterOutputMinMax(int const, double&, double&, int const, double const);
1539 int nInterpolateWavesToPolygonCells(vector<double> const*, vector<double> const*, vector<double> const*, vector<double> const*);
1540
1541 // Initialization
1542 bool bCreateErosionPotentialLookUp(vector<double>*, vector<double>*, vector<double>*);
1543
1544 // Top-level simulation routines
1545 static int nUpdateIntervention(void);
1547 int nCalcExternalForcing(void);
1549 int nLocateSeaAndCoasts(int&);
1550 int nLocateFloodAndCoasts(void);
1553 int nDoAllPropagateWaves(void);
1556 int nDoCliffCollapse(int const, CRWCliff*, double&, double&, double&, double&, double&);
1557 int nDoCliffCollapseDeposition(int const, CRWCliff const*, double const, double const, double const, double const);
1558 int nUpdateGrid(void);
1559
1560 // Lower-level simulation routines
1561 void FindAllSeaCells(void);
1562 int FindAllInundatedCells(void);
1563 void FloodFillSea(int const, int const);
1564 void FloodFillLand(int const, int const);
1565 int nTraceCoastLine(unsigned int const, int const, int const, vector<bool>*, vector<CGeom2DIPoint> const*);
1566 int nTraceAllCoasts(int&);
1567 int nTraceFloodCoastLine(unsigned int const, int const, int const, vector<bool>*, vector<CGeom2DIPoint> const*);
1568 int nTraceAllFloodCoasts(void);
1569 void DoCoastCurvature(int const, int const);
1570 int nCheckAllProfiles(void);
1571 int nCreateAllProfiles(void);
1572 void LocateAndCreateProfiles(int const, int&, int&, int const, vector<bool>*, vector<pair<int, double>> const*);
1573 int nCreateProfile(int const, int const, int const, int const, int&, bool const, CGeom2DIPoint const*);
1574 int nLocateAndCreateGridEdgeProfile(bool const, int const, int&, int&);
1575 void MarkProfilesOnGrid(int const, int&);
1577 static bool bCheckForIntersection(CGeomProfile* const, CGeomProfile* const, int&, int&, double&, double&, double&, double&);
1578 void MergeProfilesAtFinalLineSegments(int const, CGeomProfile*, CGeomProfile*, int const, int const, double const, double const, double const, double const);
1579 void TruncateOneProfileRetainOtherProfile(int const, CGeomProfile*, CGeomProfile*, double const, double const, int const, int const, bool const);
1580 int nInsertPointIntoProfilesIfNeededThenUpdate(int const, CGeomProfile*, double const, double const, int const, CGeomProfile*, int const, bool const);
1581 void TruncateProfileAndAppendNew(int const, CGeomProfile*, int const, vector<CGeom2DPoint> const*, vector<vector<pair<int, int>>> const*);
1582 void CreateRasterizedProfile(int const, CGeomProfile*, CGeom2DIPoint const*, CGeom2DIPoint const*, vector<CGeom2DIPoint>*, vector<bool>*, bool&, bool&, bool&, bool& /*, bool&*/); // TODO 044
1583 static void CalcDeanProfile(vector<double>*, double const, double const, double const, bool const, int const, double const);
1584 static double dSubtractProfiles(vector<double> const*, vector<double> const*, vector<bool> const*);
1585 void RasterizeCliffCollapseProfile(vector<CGeom2DPoint> const*, vector<CGeom2DIPoint>*) const;
1588 void ConstructParallelProfile(int const, int const, int const, int const, int const, vector<CGeom2DIPoint>* const, vector<CGeom2DIPoint>*, vector<CGeom2DPoint>*);
1589 double dCalcBeachProtectionFactor(int const, int const, double const);
1590 void FillInBeachProtectionHoles(void);
1592 void DoActualPlatformErosionOnCell(int const, int const);
1593 double dLookUpErosionPotential(double const) const;
1594 static CGeom2DPoint PtChooseEndPoint(int const, CGeom2DPoint const*, CGeom2DPoint const*, double const, double const, double const, double const);
1595 int nGetCoastNormalEndPoint(int const, int const, int const, CGeom2DPoint const*, double const, CGeom2DPoint*, CGeom2DIPoint*, bool const);
1596 int nLandformToGrid(int const, int const);
1597 int nCalcWavePropertiesOnProfile(int const, int const, CGeomProfile*, vector<double>*, vector<double>*, vector<double>*, vector<double>*, vector<bool>*);
1598 int nGetThisProfileElevationsForCShore(int const, CGeomProfile*, int const, vector<double>*, vector<double>*, vector<double>*);
1599 int nCreateCShoreInfile(int const, int const, int const, int const, int const, int const, int const, int const, int const, int const, int const, int const, int const, double const, double const, double const, double const, double const, double const, double const, double const, vector<double> const*, vector<double> const*, vector<double> const*);
1600 int nReadCShoreOutput(int const, string const*, int const, int const, vector<double> const*, vector<double>*);
1601/* static */ void InterpolateCShoreOutput(vector<double> const*, int const, int const, vector<double>*, vector<double>*, vector<double>*, vector<double>*, vector<double>*, vector<double>*, vector<double>*, vector<double>*, vector<double>*);
1602 static double dCalcWaveAngleToCoastNormal(double const, double const, int const);
1603 void CalcCoastTangents(int const);
1604 void InterpolateWavePropertiesBetweenProfiles(int const, int const);
1605 void InterpolateWaveHeightToCoastPoints(int const);
1606 // void InterpolateWavePropertiesToCells(int const, int const, int const);
1608 static double dCalcCurvature(int const, CGeom2DPoint const*, CGeom2DPoint const*, CGeom2DPoint const*);
1609 void CalcD50AndFillWaveCalcHoles(void);
1610 int nDoAllShadowZones(void);
1611 static bool bOnOrOffShoreAndUpOrDownCoast(double const, double const, int const, bool&);
1612 static CGeom2DIPoint PtiFollowWaveAngle(CGeom2DIPoint const*, double const, double&);
1613 // int nFindAllShadowZones(void);
1614 int nFloodFillShadowZone(int const, CGeom2DIPoint const*, CGeom2DIPoint const*, CGeom2DIPoint const*);
1615 void DoShadowZoneAndDownDriftZone(int const, int const, int const, int const);
1616 void ProcessDownDriftCell(int const, int const, int const, double const, int const);
1617 void ProcessShadowZoneCell(int const, int const, int const, CGeom2DIPoint const*, int const, int const, int const);
1618 int nCreateAllPolygons(void);
1619 void RasterizePolygonJoiningLine(CGeom2DIPoint const*, CGeom2DIPoint const*, int const);
1620 static bool bIsWithinPolygon(CGeom2DPoint const*, vector<CGeom2DPoint> const*);
1621 static CGeom2DPoint PtFindPointInPolygon(vector<CGeom2DPoint> const*, int const);
1622 void MarkPolygonCells(void);
1624 void DoAllPotentialBeachErosion(void);
1626 int nDoParallelProfileUnconsErosion(CGeomCoastPolygon*, int const, int const, int const, int const, int const, int const, vector<CGeom2DIPoint> const*, vector<double> const*, double&, double&, double&);
1627 void ErodeCellBeachSedimentSupplyLimited(int const, int const, int const, int const, double const, double&);
1628 int nDoUnconsErosionOnPolygon(int const, CGeomCoastPolygon*, int const, double const, double&);
1629 int nDoUnconsDepositionOnPolygon(int const, CGeomCoastPolygon*, int const, double, double&);
1630 void CalcDepthOfClosure(void);
1633 int nDoSedimentInputEvent(int const);
1634 void AllPolygonsUpdateStoredUncons(int const);
1635 bool bIsInterventionCell(int const, int const) const;
1636 bool bSurroundedByDriftCells(int const, int const);
1637 bool bElevAboveDeanElev(int const, int const, double const, CRWCellLandform const*);
1638 // void CreatePolygonIndexIDSeq(int const);
1639
1640 // GIS utility routines
1641 int nMarkBoundingBoxEdgeCells(void);
1642 bool bCheckRasterGISOutputFormat(void);
1643 bool bCheckVectorGISOutputFormat(void);
1644 bool bSaveAllRasterGISFiles(void);
1645 bool bSaveAllVectorGISFiles(void);
1646 bool bIsWithinValidGrid(int const, int const) const;
1647 bool bIsWithinValidGrid(CGeom2DIPoint const*) const;
1648 double dGridCentroidXToExtCRSX(int const) const;
1649 double dGridCentroidYToExtCRSY(int const) const;
1650 double dGridXToExtCRSX(double const) const;
1651 double dGridYToExtCRSY(double const) const;
1652 // double dExtCRSXToGridCentroidX(double const) const;
1653 // double dExtCRSYToGridCentroidY(double const) const;
1656 double dExtCRSXToGridX(double const) const;
1657 double dExtCRSYToGridY(double const) const;
1658 static double dGetDistanceBetween(CGeom2DPoint const*, CGeom2DPoint const*);
1659 static double dGetDistanceBetween(CGeom2DIPoint const*, CGeom2DIPoint const*);
1660 static double dTriangleAreax2(CGeom2DPoint const*, CGeom2DPoint const*, CGeom2DPoint const*);
1661 void KeepWithinValidGrid(int, int, int&, int&) const;
1662 void KeepWithinValidGrid(CGeom2DIPoint const*, CGeom2DIPoint*) const;
1663 static double dKeepWithin360(double const);
1664 // vector<CGeom2DPoint> VGetPerpendicular(CGeom2DPoint const*, CGeom2DPoint const*, double const, int const);
1665 static CGeom2DPoint PtGetPerpendicular(CGeom2DPoint const*, CGeom2DPoint const*, double const, int const);
1666 static CGeom2DIPoint PtiGetPerpendicular(CGeom2DIPoint const*, CGeom2DIPoint const*, double const, int const);
1667 static CGeom2DIPoint PtiGetPerpendicular(int const, int const, int const, int const, double const, int const);
1668 static CGeom2DPoint PtAverage(CGeom2DPoint const*, CGeom2DPoint const*);
1669 static CGeom2DPoint PtAverage(vector<CGeom2DPoint>*);
1670 // static CGeom2DIPoint PtiAverage(CGeom2DIPoint const*, CGeom2DIPoint const*);
1671 // static CGeom2DIPoint PtiAverage(vector<CGeom2DIPoint>*);
1672 static CGeom2DIPoint PtiWeightedAverage(CGeom2DIPoint const*, CGeom2DIPoint const*, double const);
1673 static CGeom2DIPoint PtiPolygonCentroid(vector<CGeom2DIPoint>*);
1674 static double dAngleSubtended(CGeom2DIPoint const*, CGeom2DIPoint const*, CGeom2DIPoint const*);
1675 static int nGetOppositeDirection(int const);
1676 // static void GetSlopeAndInterceptFromPoints(CGeom2DIPoint const*, CGeom2DIPoint const*, double&, double&);
1677 CGeom2DIPoint PtiFindClosestCoastPoint(int const, int const);
1678 int nConvertMetresToNumCells(double const) const;
1679
1680 // Utility routines
1681 static void AnnounceStart(void);
1682 void AnnounceLicence(void);
1683 void AnnounceReadBasementDEM(void) const;
1684 static void AnnounceAddLayers(void);
1685 static void AnnounceReadRasterFiles(void);
1686 static void AnnounceReadVectorFiles(void);
1687 void AnnounceReadLGIS(void) const;
1688 void AnnounceReadICGIS(void) const;
1689 void AnnounceReadIHGIS(void) const;
1690 static void AnnounceInitializing(void);
1691 void AnnounceReadInitialSuspSedGIS(void) const;
1692 void AnnounceReadInitialFineUnconsSedGIS(int const) const;
1693 void AnnounceReadInitialSandUnconsSedGIS(int const) const;
1694 void AnnounceReadInitialCoarseUnconsSedGIS(int const) const;
1695 void AnnounceReadInitialFineConsSedGIS(int const) const;
1696 void AnnounceReadInitialSandConsSedGIS(int const) const;
1697 void AnnounceReadInitialCoarseConsSedGIS(int const) const;
1698 void AnnounceReadDeepWaterWaveValuesGIS(void) const;
1700 void AnnounceReadFloodLocationGIS(void) const;
1701 void AnnounceReadTideData(void) const;
1702 static void AnnounceReadSCAPEShapeFunctionFile(void);
1703 static void AnnounceAllocateMemory(void);
1704 static void AnnounceIsRunning(void);
1705 static void AnnounceSimEnd(void);
1706 void StartClock(void);
1707 bool bFindExeDir(char const*);
1708 bool bTimeToQuit(void);
1709 static int nDoTimeUnits(string const*);
1710 int nDoSimulationTimeMultiplier(string const*);
1711 static double dGetTimeMultiplier(string const*);
1712 static bool bParseDate(string const*, int&, int&, int&);
1713 static bool bParseTime(string const*, int&, int&, int&);
1714 void DoEndOfTimestepTotals(void);
1715 static string strGetBuild(void);
1716 static string strGetComputerName(void);
1717 void DoCPUClockReset(void);
1718 void CalcTime(double const);
1719 static string strDispTime(double const, bool const, bool const);
1720 static string strDispSimTime(double const);
1721 void AnnounceProgress(void);
1722 static string strGetErrorText(int const);
1723 string strListRasterFiles(void) const;
1724 string strListVectorFiles(void) const;
1725 string strListTSFiles(void) const;
1726 void CalcProcessStats(void);
1727 void CalcSavitzkyGolayCoeffs(void);
1728 CGeomLine LSmoothCoastSavitzkyGolay(CGeomLine*, int const, int const) const;
1730 vector<double> dVSmoothProfileSlope(vector<double>*) const;
1731 // vector<double> dVCalCGeomProfileSlope(vector<CGeom2DPoint>*, vector<double>*); // TODO 007 Why was this removed?
1732 // vector<double> dVSmoothProfileSavitzkyGolay(vector<double>*, vector<double>*); // TODO 007 was this removed?
1733 // vector<double> dVSmoothProfileRunningMean(vector<double>*); // TODO 007 was this removed?
1734 static void CalcSavitzkyGolay(double[], int const, int const, int const, int const, int const);
1735 static string pstrChangeToBackslash(string const*);
1736 static string pstrChangeToForwardSlash(string const*);
1737 static string strTrim(string const*);
1738 static string strTrimLeft(string const*);
1739 static string strTrimRight(string const*);
1740 static string strToLower(string const*);
1741 // static string strToUpper(string const*);
1742 static string strRemoveSubstr(string *, string const*);
1743 static vector<string>* VstrSplit(string const*, char const, vector<string>*);
1744 static vector<string> VstrSplit(string const*, char const);
1745 // static double dCrossProduct(double const, double const, double const, double const, double const, double const);
1746 // static double dGetMean(vector<double> const*);
1747 // static double dGetStdDev(vector<double> const*);
1748 static void AppendEnsureNoGap(vector<CGeom2DIPoint>*, CGeom2DIPoint const*);
1749 // static bool bIsNumeric(string const*);
1750 unsigned long ulConvertToTimestep(string const*) const;
1751 void WritePolygonShareTable(int const);
1754 void WritePolygonShorePlatformErosion(int const);
1755 void WritePolygonCliffCollapseErosion(int const);
1756 void WritePolygonSedimentBeforeMovement(int const);
1757 void WritePolygonPotentialErosion(int const);
1758 // void WritePolygonUnconsErosion(int const);
1759 void WritePolygonUnsortedSequence(int const, vector<vector<int> >&);
1760 void WritePolygonSortedSequence(int const, vector<vector<int> >&);
1761 void WritePolygonEstimatedMovement(int const, vector<vector<int> >&);
1762 void WritePolygonActualMovement(int const, vector<vector<int> > const&);
1763 void DoEndOfRunDeletes(void);
1764
1765protected:
1766
1767public:
1768 ofstream LogStream;
1769
1770 CSimulation(void);
1771 ~CSimulation(void);
1772
1774 int nGetCoastPolygonSize(void) const;
1775
1777 CGeomCoastPolygon* pGetPolygon(int const) const;
1778
1781
1783 double dGetMissingValue(void) const;
1784
1786 double dGetThisIterSWL(void) const;
1787
1789 double dGetThisIterTotWaterLevel(void) const;
1790
1791 // //! Returns the vertical tolerance for beach cells to be included in smoothing
1792 // double dGetMaxBeachElevAboveSWL(void) const;
1793
1795 // double dGetCellSide(void) const;
1796
1798 int nGetGridXMax(void) const;
1799
1801 int nGetGridYMax(void) const;
1802
1804 double dGetD50Fine(void) const;
1805
1807 double dGetD50Sand(void) const;
1808
1810 double dGetD50Coarse(void) const;
1811
1813 int nDoSimulation(int, char const* []);
1814
1816 void DoSimulationEnd(int const);
1817};
1818#endif // SIMULATION_H
Geometry class used to represent 2D point objects with integer coordinates.
Definition 2di_point.h:29
Geometry class used to represent 2D point objects with floating-point coordinates.
Definition 2d_point.h:27
Geometry class used for coast polygon objects.
Geometry class used to represent 2D vector line objects.
Definition line.h:31
Geometry class used to represent coast profile objects.
Definition profile.h:35
Geometry cass used to represent the raster grid of cell objects.
Definition raster_grid.h:35
Real-world class used to represent the landform of a cell.
Real-world class used to represent the 'cliff' category of coastal landform object.
Definition cliff.h:33
Real-world class used to represent coastline objects.
Definition coast.h:39
Class used to represent a sediment input event.
double m_dThisIterPotentialBeachErosion
Total potential beach erosion (all size classes of unconsolidated sediment) for this iteration (depth...
Definition simulation.h:796
double m_dCliffDepositionPlanviewWidth
Planview width of cliff collapse talus (m)
Definition simulation.h:871
void CalcDepthOfClosure(void)
Calculate the depth of closure.
Definition utils.cpp:2437
bool m_bCliffCollapseSave
Save cliff collapse raster GIS files?
Definition simulation.h:206
int m_nLogFileDetail
The level of detail in the log file output. Can be LOG_FILE_LOW_DETAIL, LOG_FILE_MIDDLE_DETAIL,...
Definition simulation.h:530
bool m_bAvgSeaDepthSave
Save average sea depth raster GIS files?
Definition simulation.h:95
void WritePolygonSedimentInputEventTable(int const)
Writes to the log file a table showing per-polygon sediment input event totals.
string m_strGDALISSDriverCode
GDAL code for the initial suspended sediment raster file.
vector< string > m_VstrInitialSandConsSedimentFile
The name of the initial sand-sized consolidated sediment GIS file.
double m_dAllCellsDeepWaterWaveHeight
Deep water wave height (m) for all sea cells.
Definition simulation.h:709
bool m_bDeepWaterWaveAngleSave
Save deep water wave angle raster GIS files?
Definition simulation.h:233
bool bWriteRasterGISFile(int const, string const *, int const =0, double const =0)
Writes GIS raster files using GDAL, using data from the RasterGrid array.
static void AnnounceInitializing(void)
Tells the user that we are now initializing.
Definition utils.cpp:618
int m_nGISMaxSaveDigits
The maximum number of digits in GIS filenames. These can be sequential, or the iteration number.
Definition simulation.h:458
bool m_bTopSurfSave
Save fop surface (sediment and sea) raster DEMs?
Definition simulation.h:86
int m_nYMinBoundingBox
The minimum y value of the bounding box.
Definition simulation.h:500
static void AnnounceReadSCAPEShapeFunctionFile(void)
Now reading the SCAPE shape function file.
Definition utils.cpp:610
string m_strGDALISSProjection
GDAL projection string for the initial suspended sediment raster file.
void CalcSavitzkyGolayCoeffs(void)
Calculates the Savitzky-Golay smoothing coefficients for a given size of smoothing window....
string m_strInitialSuspSedimentFile
Name of initial suspended sediment file.
void FindAllSeaCells(void)
Finds and flags all sea areas which have at least one cell at a grid edge (i.e. does not flag 'inland...
vector< double > m_VdTSDeepWaterWaveStationPeriod
Time series of wave period at deep water wave station.
double dLookUpErosionPotential(double const) const
The erosion potential lookup: it returns a value for erosion potential given a value of Depth Over DB...
static void CalcSavitzkyGolay(double[], int const, int const, int const, int const, int const)
CalcSavitzkyGolay uses LULinearSolve and LUDecomp. It returns dFilterCoeffsArray[nWindowSize],...
string m_strOGRVectorOutputExtension
GDAL-OGR vector output drive file extension.
CGeomLine LSmoothCoastRunningMean(CGeomLine *) const
Does running-mean smoothing of a CGeomLine coastline vector (is in external CRS coordinates)
bool bCheckRasterGISOutputFormat(void)
Checks whether the selected raster GDAL driver supports file creation, 32-bit doubles,...
void InterpolateCShoreOutput(vector< double > const *, int const, int const, vector< double > *, vector< double > *, vector< double > *, vector< double > *, vector< double > *, vector< double > *, vector< double > *, vector< double > *, vector< double > *)
static void AnnounceIsRunning(void)
Tell the user that the simulation is now running.
Definition utils.cpp:627
void DoCPUClockReset(void)
Resets the CPU clock timer to prevent it 'rolling over', as can happen during long runs....
Definition utils.cpp:1343
bool m_bSedimentTopSurfSave
Save sediment top surface raster DEMs?
Definition simulation.h:83
bool m_bFineUnconsSedSave
Save fine unconsolidated sediment raster GIS files?
Definition simulation.h:179
string m_strSedimentInputEventFile
The name of the sediment input events time series file.
int nLocateSeaAndCoasts(int &)
First find all connected sea areas, then locate the vector coastline(s), then put these onto the rast...
double m_dNotchDepthAtCollapse
Notch overhang (i.e. length of horizontal incision) to initiate collapse (m)
Definition simulation.h:862
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_bFloodSWLSetupSurgeLine
Are we saving the flood still water level setup surge line? TODO 007.
Definition simulation.h:416
time_t m_tSysEndTime
System finish-simulation time.
string m_strCMEIni
Folder for the CME .ini file.
double m_dTotalCoarseUnconsInPolygons
Total coarse unconsolidated sediment in all polygons, before polygon-to-polygon movement (only cells ...
Definition simulation.h:961
bool m_bSedimentInputAtPoint
Do we have sediment inputat a point?
Definition simulation.h:371
double m_dMinSWL
Minimum still water level.
Definition simulation.h:676
double m_dG
Gravitational acceleration (m**2/sec)
Definition simulation.h:766
vector< string > m_VstrInitialCoarseUnconsSedimentFile
The name of the initial coarse-sized unconsolidated sediment GIS file.
void AnnounceReadDeepWaterWaveValuesGIS(void) const
Tells the user that we are now reading the deep water wave values GIS file.
Definition utils.cpp:465
double m_dThisIterSWL
The still water level for this timestep (this includes tidal changes and any long-term SWL change)
Definition simulation.h:667
string m_strGDALICDataType
GDAL data type of the initial intervention class raster file.
int m_nBeachErosionDepositionEquation
Which beach erosion-deposition equation is used. Possible values are UNCONS_SEDIMENT_EQUATION_CERC an...
Definition simulation.h:488
ofstream CliffCollapseNetChangeTSStream
Cliff collapse net change (erosion - deposition) time series file output stream.
bool m_bCoastSave
Save.
Definition simulation.h:251
bool m_bBeachDepositionTSSave
Save the beach (unconsolidated sediment) deposition time series file?
Definition simulation.h:299
bool bWritePerTimestepResults(void)
Write the results for this timestep to the .out file.
static string strRemoveSubstr(string *, string const *)
Returns a string with a substring removed, and with whitespace trimmed.
Definition utils.cpp:2239
CGeomRasterGrid * m_pRasterGrid
Pointer to the raster grid object.
vector< string > m_VstrGDALICCDataType
GDAL data type for the initial consolidated coarse sediment GIS data.
int nCalcWavePropertiesOnProfile(int const, int const, CGeomProfile *, vector< double > *, vector< double > *, vector< double > *, vector< double > *, vector< bool > *)
Calculates wave properties along a coastline-normal profile using either the COVE linear wave theory ...
static bool bParseTime(string const *, int &, int &, int &)
Parses a time string into hours, minutes, and seconds, and checks each of them.
Definition utils.cpp:2541
double dGetThisIterTotWaterLevel(void) const
Returns this timestep's total water level.
void AnnounceReadInitialSandConsSedGIS(int const) const
Tells the user that we are now reading the initial sand consolidated sediment depth GIS file.
Definition utils.cpp:572
int m_nXGridSize
The size of the grid in the x direction.
Definition simulation.h:431
string strListRasterFiles(void) const
Return a space-separated string containing the names of the raster GIS output files.
Definition utils.cpp:635
vector< double > m_VdDeepWaterWaveStationY
Y coordinate (grid CRS) for deep water wave station.
int nAssignLandformsForAllCoasts(void)
Each timestep, classify coastal landforms and assign a coastal landform object to every point on ever...
void RasterizeCliffCollapseProfile(vector< CGeom2DPoint > const *, vector< CGeom2DIPoint > *) const
Given the start and end points of a cliff-collapse normal profile, returns an output vector of cells ...
int nFloodFillShadowZone(int const, CGeom2DIPoint const *, CGeom2DIPoint const *, CGeom2DIPoint const *)
Flood fills a shadow zone from the centroid.
void AnnounceReadInitialSandUnconsSedGIS(int const) const
Tells the user that we are now reading the initial sand unconsolidated sediment depth GIS file.
Definition utils.cpp:533
double m_dL_0
Deep water wave length (m)
Definition simulation.h:700
double m_dWaveDataWrapHours
Number of hours after which deep water wave data wraps.
Definition simulation.h:913
static int nGetOppositeDirection(int const)
Returns the opposite direction.
string strListVectorFiles(void) const
Return a space-separated string containing the names of the vector GIS output files.
Definition utils.cpp:888
int nDoCliffCollapseDeposition(int const, CRWCliff const *, double const, double const, double const, double const)
Redistributes the sand-sized and coarse-sized sediment from a cliff collapse onto the foreshore,...
CGeom2DIPoint PtiExtCRSToGridRound(CGeom2DPoint const *) const
Transforms a pointer to a CGeom2DPoint in the external CRS to the equivalent CGeom2DIPoint in the ras...
vector< string > m_VstrGDALICFDataType
GDAL data type for the initial consolidated fine sediment GIS data.
vector< int > m_VCellFloodLocation
The location to compute the total water level for flooding.
double m_dMaxUserInputWavePeriod
Used to constrain depth of closure.
Definition simulation.h:721
bool m_bFloodSetupSurgeRunupTSSave
Save the flood setup surge runup time series file? TODO 007 Does this work correctly?
Definition simulation.h:311
bool bSetUpTSFiles(void)
The bSetUpTSFiles member function sets up the time series files.
Definition utils.cpp:1078
ofstream LogStream
int nWriteEndRunDetails(void)
Writes end-of-run information to Out, Log and time-series files.
void AnnounceReadInitialFineUnconsSedGIS(int const) const
Tells the user that we are now reading the initial fine unconsolidated sediment depth GIS file.
Definition utils.cpp:520
vector< string > m_VstrGDALICCDriverCode
GDAL driver code for the initial consolidated coarse sediment GIS data.
double m_dThisIterBeachErosionCoarse
Total actual beach erosion (coarse unconsolidated sediment) for this iteration (depth in m)
Definition simulation.h:805
char ** m_papszGDALVectorOptions
Options for GDAL when handling vector files.
Definition simulation.h:428
vector< CRWCoast > m_VFloodWaveSetupSurgeRunup
TODO 007 Info needed.
double m_dStartIterUnconsCoarseAllCells
Depth (m) of coarse unconsolidated sediment at the start of the simulation, all cells (both inside an...
Definition simulation.h:943
int nReadRasterGISFile(int const, int const)
Reads all other raster GIS datafiles into the RasterGrid array.
vector< CRWCoast > m_VCoast
The coastline objects.
double m_dCoastNormalLength
Length of the cost-normal profiles, in m.
Definition simulation.h:778
void WritePolygonCliffCollapseErosion(int const)
Writes to the log file a table showing per-polygon per-polygon cliff collapse.
int nLocateFloodAndCoasts(void)
First find all connected sea areas, then locate the vector coastline(s), then put these onto the rast...
static string pstrChangeToForwardSlash(string const *)
Swaps all backslashes in the input string to forward slashes, leaving the original unchanged.
Definition utils.cpp:2153
vector< string > m_VstrGDALIUFDataType
GDAL data type for the initial unconsolidated fine sediment GIS data.
bool m_bSaveGISThisIter
Save GIS files this iteration?
Definition simulation.h:314
long double m_ldGTotCliffCollapseCoarseErodedDuringDeposition
All-simulation total of coarse sediment eroded during talus deposition following cliff collapse (m)
int nTraceCoastLine(unsigned int const, int const, int const, vector< bool > *, vector< CGeom2DIPoint > const *)
Traces a coastline (which is defined to be just above still water level) on the grid using the 'wall ...
default_random_engine m_Rand[NRNG]
The c++11 random number generators.
double m_dUnconsCoarseNotDepositedLastIter
Depth of unconsolidated coarse sediment that could not be deposited during the last iteration,...
Definition simulation.h:976
string m_strOGRSedInputDataType
GDAL data type for the sediment input event locations vector file.
void DoSimulationEnd(int const)
Carries out end-of-simulation tidying (error messages etc.)
Definition utils.cpp:2030
double m_dCPUClock
Total elapsed CPU time.
Definition simulation.h:649
bool bSaveAllVectorGISFiles(void)
The bSaveAllvectorGISFiles member function saves the vector GIS files TODO 081 Choose more files to o...
bool m_bSingleDeepWaterWaveValues
Do we have just a point source for (i.e. only a single measurement of) deep water wave values.
Definition simulation.h:362
string m_strRunName
The name of this simulation.
int m_nThisSave
Used in calculations of GIS save intervals.
Definition simulation.h:467
static string strGetComputerName(void)
Returns a string, hopefully giving the name of the computer on which the simulation is running.
Definition utils.cpp:1320
long m_lGDALMaxCanWrite
The maximum integer value which GDAL can write, can be UINT8_MAX, INT16_MAX, UINT16_MAX,...
Definition simulation.h:547
static string strDispSimTime(double const)
strDispSimTime returns a string formatted as year Julian_day hour, given a parameter in hours
Definition utils.cpp:1473
bool m_bAvgWaveAngleAndHeightSave
Save average wave angle and average wave height raster GIS files?
Definition simulation.h:113
int nTraceAllFloodCoasts(void)
Locates all the potential coastline start points on the edges of the raster grid, then traces vector ...
void FloodFillLand(int const, int const)
Use the sealevel, wave set-up and run-up to evaluate flood hydraulically connected TODO 007 Not clear...
void MarkProfilesOnGrid(int const, int &)
For this coastline, marks all coastline-normal profiles (apart from the two 'special' ones at the sta...
double m_dSouthEastXExtCRS
The south-east x coordinate, in the external coordinate reference system (CRS)
Definition simulation.h:604
string m_strSCAPEShapeFunctionFile
Name of SCAPE shape function file.
int m_nMissingValue
The value used for integer missing values.
Definition simulation.h:491
bool m_bAvgWaveAngleSave
Save average wave angle raster GIS files?
Definition simulation.h:107
long double m_ldGTotActualCoarseBeachErosion
All-simulation total of coarse sediment eroded during beach (unconsolidated sediment) movement (m)
int nCalcExternalForcing(void)
Calculate external forcings: change in still water level, tide level and deep water waves height,...
int nReadWaveStationInputFile(int const)
Reads the deep water wave station input data. Each point in m_strDeepWaterWavesInputFile is a triad o...
double m_dFineErodibilityNormalized
Relative erodibility of fine unconsolidated beach sediment, normalized.
Definition simulation.h:751
void WriteLookUpData(void) const
Output the erosion potential look-up values, for checking purposes.
bool bWriteVectorGISFile(int const, string const *)
Writes vector GIS files using OGR.
double m_dThisIterCliffCollapseFineErodedDuringDeposition
Total fine sediment eroded during Dean profile deposition of talus following cliff collapse (depth in...
Definition simulation.h:829
void WritePolygonSedimentBeforeMovement(int const)
Writes to the log file a table showing per-polygon totals of stored unconsolidated beach sediment pri...
bool m_bInvalidNormalsSave
Save invalid coastline-normal vector GIS files?
Definition simulation.h:257
bool m_bShadowBoundarySave
Save wave shadow boundary vector GIS files?
Definition simulation.h:272
double m_dThisIterDiffWaveSetupSurgeWaterLevel
TODO 007 Info needed.
Definition simulation.h:688
static CGeom2DPoint PtChooseEndPoint(int const, CGeom2DPoint const *, CGeom2DPoint const *, double const, double const, double const, double const)
Choose which end point to use for the coastline-normal profile.
int nDoSimulationTimeMultiplier(string const *)
Given a string containing time units, this sets up the appropriate multiplier and display units for t...
Definition utils.cpp:298
int nReadSedimentInputEventFile(void)
Reads the sediment input event file.
vector< int > m_VnSedimentInputLocationID
ID for sediment input location, this corresponds with the ID in the sediment input time series file.
double m_dThisIterActualPlatformErosionCoarseCons
Total actual platform erosion (coarse consolidated sediment) for this iteration (depth in m)
Definition simulation.h:793
bool m_bActualBeachErosionSave
Save actual (supply-limited) beach (unconsolidated sediment) erosion raster GIS files?
Definition simulation.h:146
int nDoAllShorePlatFormErosion(void)
Does platform erosion on all coastlines by first calculating platform erosion on coastline-normal pro...
int m_nYGridSize
The size of the grid in the y direction.
Definition simulation.h:434
vector< double > m_VdErosionPotential
For erosion potential lookup.
double m_dThisIterTopElevMin
This-iteration lowest elevation of DEM.
Definition simulation.h:919
vector< double > m_VdTideData
Tide data: one record per timestep, is the change (m) from still water level for that timestep.
vector< string > m_VstrGDALIUFProjection
GDAL projection for the initial unconsolidated fine sediment GIS data.
bool m_bStillWaterLevelTSSave
Save the still water level time series file?
Definition simulation.h:281
bool m_bRunUpSave
Are we saving runup? TODO 007.
Definition simulation.h:395
void ConstructParallelProfile(int const, int const, int const, int const, int const, vector< CGeom2DIPoint > *const, vector< CGeom2DIPoint > *, vector< CGeom2DPoint > *)
Constructs a parallel coastline-normal profile.
GDALDataType m_GDALWriteIntDataType
The data type used by GDAL for integer operations, can be GDT_Byte, GDT_Int16, GDT_UInt16,...
Definition simulation.h:541
long double m_ldGTotCliffTalusFineToSuspension
All-simulation total of fine sediment moved to suspension, due to cliff collapse (m)
void ErodeCellBeachSedimentSupplyLimited(int const, int const, int const, int const, double const, double &)
Erodes the unconsolidated beach sediment on a single cell, for a single size class,...
bool m_bCliffCollapseDepositionSave
Save cliff collapse deposition raster GIS files?
Definition simulation.h:212
int nSaveParProfile(int const, CGeomProfile const *, int const, int const, int const, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< CGeom2DIPoint > *const, vector< double > const *) const
Save a coastline-normal parallel profile.
string m_strSedimentInputEventShapefile
The name of the sediment input events shape file.
bool m_bTotalActualPlatformErosionSave
Save total actual (supply-limited) shore platform erosion raster GIS files?
Definition simulation.h:140
static string strTrimLeft(string const *)
Trims whitespace from the left side of a string, does not change the original string.
Definition utils.cpp:2163
bool m_bPolygonUnconsSedUpOrDownDriftSave
Save polygon unconsolidated sediment up- or down-drift raster GIS files?
Definition simulation.h:242
double m_dCoarseErodibility
The relative erodibility (0- 1) of coarse unconsolidated beach sediment.
Definition simulation.h:748
double m_dGeoTransform[6]
GDAL geotransformation info (see http://www.gdal.org/classGDALDataset.html)
Definition simulation.h:652
double m_dCliffTalusMinDepositionLength
Planview length of cliff deposition talus (m)
Definition simulation.h:874
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...
void WritePolygonSortedSequence(int const, vector< vector< int > > &)
Writes to the log file a table showing the sorted sequence of polygon processing, and any circulariti...
string m_strVectorGISOutFormat
Base name for CME vector GIS output files.
void ProcessShadowZoneCell(int const, int const, int const, CGeom2DIPoint const *, int const, int const, int const)
Process a single cell which is in the shadow zone, changing its wave height and orientation.
vector< string > m_VstrGDALICFDriverCode
GDAL driver code for the initial consolidated fine sediment GIS data.
vector< string > m_VstrInitialFineConsSedimentFile
The name of the initial fine-sized consolidated sediment GIS file.
double m_dMissingValue
Missing value.
Definition simulation.h:910
int m_nUnconsSedimentHandlingAtGridEdges
How sediment which moves off an edge of the grid is handled. Possible values are GRID_EDGE_CLOSED,...
Definition simulation.h:485
static double dGetDistanceBetween(CGeom2DPoint const *, CGeom2DPoint const *)
Returns the distance (in external CRS) between two points.
double m_dSandErodibility
The relative erodibility (0- 1) of sand unconsolidated beach sediment.
Definition simulation.h:745
bool m_bBasementElevSave
Save basement raster DEMs?
Definition simulation.h:80
double m_dInvCellDiagonal
Inverse of m_dCellDiagonal.
Definition simulation.h:625
int m_nSimStartHour
Start time of the simulation (hours)
Definition simulation.h:515
int nDoAllShadowZones(void)
Finds wave shadow zones and modifies waves in and near them. Note that where up-coast and down-coast ...
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...
double m_dCoarseErodibilityNormalized
Relative erodibility of coarse unconsolidated beach sediment, normalized.
Definition simulation.h:757
bool m_bCoarseUnconsSedSave
Save coarse unconsolidated sediment raster GIS files?
Definition simulation.h:185
string m_strGDALLDriverCode
GDAL code for the for the initial landform class raster file.
bool m_bSuspSedTSSave
Save the suspended sediment time series file?
Definition simulation.h:305
static double dKeepWithin360(double const)
Constrains the supplied angle to be within 0 and 360 degrees.
string m_strDeepWaterWaveStationsShapefile
The name of the deep water wave stations shape file.
vector< double > m_VdDepthOverDB
For erosion potential lookup.
void CalcCoastTangents(int const)
Calculates tangents to a coastline: the tangent is assumed to be the orientation of energy/sediment f...
bool m_bCliffCollapseNetTSSave
Save the cliff collapse net change time series file?
Definition simulation.h:293
int nLocateAndCreateGridEdgeProfile(bool const, int const, int &, int &)
Creates a 'special' profile at each end of a coastline, at the edge of the raster grid....
double m_dStartIterSuspFineInPolygons
Depth (m) of fine suspended sediment at the start of the simulation (only cells in polygons)
Definition simulation.h:934
bool m_bPotentialPlatformErosionMaskSave
Save potential platform erosion mask raster GIS files?
Definition simulation.h:221
unsigned long m_ulThisIterNumCoastCells
The number of grid cells which are marked as coast, for this iteration.
Definition simulation.h:568
void FillInBeachProtectionHoles(void)
Fills in 'holes' in the beach protection i.e. orphan cells which get omitted because of rounding prob...
static void AnnounceSimEnd(void)
Announce the end of the simulation.
Definition utils.cpp:1392
vector< string > m_VstrGDALICSDriverCode
GDAL driver code for the initial consolidated sand sediment GIS data.
unsigned long m_ulNumCells
The number of cells in the grid.
Definition simulation.h:562
double m_dTotPotentialPlatformErosionBetweenProfiles
Total potential platform erosion between profiles.
Definition simulation.h:850
int nGetGridYMax(void) const
Returns the size of the grid in the Y direction.
void AnnounceReadICGIS(void) const
Tells the user that we are now reading the Intervention class GIS file.
Definition utils.cpp:437
int nCheckForSedimentInputEvent(void)
Check to see if we have any sediment input events this timestep, if so then do the event(s)
double m_dWaveDepthRatioForWaveCalcs
Start depth for wave calculations.
Definition simulation.h:703
bool m_bWaveHeightSave
Save wave height raster GIS files?
Definition simulation.h:98
bool m_bFloodLocation
Are we saving the flood location? TODO 007.
Definition simulation.h:410
string m_strGDALISSDataType
GDAL data type for the initial suspended sediment raster file.
void WriteStartRunDetails(void)
Writes beginning-of-run information to Out and Log files.
double m_dThisIterCliffCollapseErosionCoarseUncons
This-iteration total of coarse unconsolidated sediment produced by cliff collapse (m^3)
Definition simulation.h:886
vector< string > m_VstrGDALIUCProjection
GDAL projection for the initial unconsolidated coarse sediment GIS data.
bool m_bGDALCanCreate
Is the selected GDAL output file format capable of writing files?
Definition simulation.h:347
bool m_bLandformSave
Save coast landform raster GIS files?
Definition simulation.h:161
ofstream CliffCollapseDepositionTSStream
Cliff collapse deposition time series file output stream.
string m_strGDALLDataType
GDAL data type for the initial landform class raster file.
vector< CRWCoast > m_VFloodWaveSetupSurge
TODO 007 Info needed.
static string strTrim(string const *)
Trims whitespace from both sides of a string, does not change the original string.
Definition utils.cpp:2194
string m_strRasterGISOutFormat
Base name for CME raster GIS output files.
string m_strGDALIWDriverCode
GDAL code for the initial water depth raster file.
long double m_ldGTotCoarseBeachDeposition
All-simulation total of coarse sediment deposited during beach (unconsolidated sediment) movement (m)
double m_dDepositionCoarseDiff
Error term: if we are unable to deposit enough unconslidated coarse on polygon(s),...
Definition simulation.h:841
bool m_bTotalBeachDepositionSave
Save total beach (unconsolidated sediment) deposition raster GIS files?
Definition simulation.h:158
int m_nCoastMax
Maximum valid coast length when searching for coasts, actually is COAST_LENGTH_MAX * tMax(m_nXGridSiz...
Definition simulation.h:470
int nCreateCShoreInfile(int const, int const, int const, int const, int const, int const, int const, int const, int const, int const, int const, int const, int const, double const, double const, double const, double const, double const, double const, double const, double const, vector< double > const *, vector< double > const *, vector< double > const *)
double m_dBeachSedimentPorosity
The porosity of unconsolidated beach sediment (0 - 1)
Definition simulation.h:739
double dGetMissingValue(void) const
Returns the NODATA value.
int m_nSimStartSec
Start time of the simulation (seconds)
Definition simulation.h:509
int m_nSimStartDay
Start date of the simulation (day)
Definition simulation.h:518
unsigned long m_ulThisIterNumActualPlatformErosionCells
The number of grid cells on which actual platform erosion occurs, for this iteration.
Definition simulation.h:574
bool m_bSandUnconsSedSave
Save sand unconsolidated sediment raster GIS files?
Definition simulation.h:182
bool m_bActualPlatformErosionTSSave
Save the actual (supply-limited) shore platform erosion time series file?
Definition simulation.h:284
unsigned long ulConvertToTimestep(string const *) const
For sediment input events, parses a string that may be relative (a number of hours or days after the ...
Definition utils.cpp:2602
int m_nXMaxBoundingBox
The maximum x value of the bounding box.
Definition simulation.h:497
string m_strOGRFloodGeometry
GDAL geometry for the flood input locations point or vector file.
long double m_ldGTotPotentialSedLostBeachErosion
All-simulation total of potential sediment lost via beach (unconsolidated) sediment movement (m),...
Definition simulation.h:993
double m_dThisIterCliffCollapseErosionFineUncons
This-iteration total of fine unconsolidated sediment produced by cliff collapse (m^3)
Definition simulation.h:880
double dGetD50Coarse(void) const
Returns the global d50 value for coarse sediment.
string m_strOGRSedInputDriverDesc
bool m_bTotCliffCollapseSave
Save total cliff collapse raster GIS files?
Definition simulation.h:209
int nReadShapeFunctionFile(void)
Reads the shape of the erosion potential distribution (see shape function in Walkden & Hall,...
void MergeProfilesAtFinalLineSegments(int const, CGeomProfile *, CGeomProfile *, int const, int const, double const, double const, double const, double const)
Merges two profiles which intersect at their final (most seaward) line segments, seaward of their poi...
static void AnnounceAddLayers(void)
Tells the user that we are now adding layers.
Definition utils.cpp:398
vector< unsigned long > m_VlDeepWaterWaveValuesAtTimestep
Calculate deep water wave values at these timesteps.
int nDoAllActualBeachErosionAndDeposition(void)
Does between-polygon and within-polygon actual (supply-limited) redistribution of transported beach s...
bool m_bDoShorePlatformErosion
Simulate shore platform erosion?
Definition simulation.h:338
bool m_bSliceSave
Save slices?
Definition simulation.h:89
void WritePolygonShorePlatformErosion(int const)
Writes to the log file a table showing per-polygon unconsolidated sand/coarse sediment derived from e...
bool m_bRasterPolygonSave
Save raster polygon raster GIS files?
Definition simulation.h:218
string m_strOGRDWWVDriverCode
GDAL code for the deep water wave stations vector file.
bool m_bBeachDepositionSave
Save beach (unconsolidated sediment) deposition raster GIS files?
Definition simulation.h:155
vector< double > m_VdFloodLocationX
X coordinate (grid CRS) for total water level flooding.
double m_dNorthWestYExtCRS
The north-west y coordinate, in the external coordinate reference system (CRS)
Definition simulation.h:601
vector< string > m_VstrGDALICSDriverDesc
GDAL driver description for the initial consolidated sand sediment GIS data.
bool m_bSaveRegular
Save GIS files at regular intervals?
Definition simulation.h:248
int m_nLayers
The number of sediment layers.
Definition simulation.h:437
bool bFindExeDir(char const *)
Finds the folder (directory) in which the CoastalME executable is located.
Definition utils.cpp:208
double m_dStartIterConsCoarseAllCells
Depth (m) of coarse consolidated sediment at the start of the simulation, all cells (both inside and ...
Definition simulation.h:952
double m_dStartIterConsSandAllCells
Depth (m) of sand consolidated sediment at the start of the simulation, all cells (both inside and ou...
Definition simulation.h:949
void InterpolateWavePropertiesBetweenProfiles(int const, int const)
Interpolates wave properties from profiles to the coastline points between two profiles....
static double dSubtractProfiles(vector< double > const *, vector< double > const *, vector< bool > const *)
Calculate the total elevation difference between every point in two elevation profiles (first profile...
Definition utils.cpp:2394
bool m_bSedimentInputAlongLine
Do we have sediment input along a line?
Definition simulation.h:377
bool m_bSedimentInput
Do we have sediment input events?
Definition simulation.h:368
unsigned long m_ulThisIterNumBeachDepositionCells
The number of grid cells on which beach (unconsolidated sediment) deposition occurs,...
Definition simulation.h:583
double m_dInmersedToBulkVolumetric
For beach erosion/deposition, conversion from immersed weight to bulk volumetric (sand and voids) tra...
Definition simulation.h:769
bool m_bNormalsSave
Save coastline-normal vector GIS files?
Definition simulation.h:254
bool m_bAvgWaveHeightSave
Save wave height raster GIS files?
Definition simulation.h:101
static string strToLower(string const *)
Returns the lower case version of an string, leaving the original unchanged.
Definition utils.cpp:2219
int nInterpolateWavesToPolygonCells(vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *)
Interpolates wave properties from all profiles to all within-polygon sea cells, using GDALGridCreate(...
vector< double > m_VdThisIterDeepWaterWaveStationHeight
This-iteration wave height at deep water wave station.
int nCalcPotentialPlatformErosionBetweenProfiles(int const, CGeomProfile *, int const)
Calculates potential platform erosion on cells to one side of a given coastline-normal profile,...
CSimulation(void)
The CSimulation constructor.
double m_dDepositionSandDiff
Error term: if we are unable to deposit enough unconslidated sand on polygon(s), this is held over to...
Definition simulation.h:838
string m_strGDALBasementDEMDriverCode
GDAL code for the basement DEM raster file type.
vector< string > m_VstrGDALIUFDriverCode
GDAL driver code for the initial unconsolidated fine sediment GIS data.
bool m_bHaveSandSediment
Does this simulation consider sand-sized sediment?
Definition simulation.h:74
unsigned long m_ulMissingValueBasementCells
The number of basement cells marked with as missing value.
Definition simulation.h:592
void WritePolygonEstimatedMovement(int const, vector< vector< int > > &)
int nSaveProfile(int const, CGeomProfile const *, int const, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< CGeom2DIPoint > *const, vector< double > const *) const
Save a coastline-normal profile.
double m_dThisiterUnconsCoarseInput
Depth (m) of coarse unconsolidated sediment added, at this iteration.
Definition simulation.h:928
int nReadVectorGISFile(int const)
Reads vector GIS datafiles using OGR.
bool bReadRunDataFile(void)
Reads the run details input file and does some initialization.
int m_nUSave
If user-defined GIS save intervals, the number of these.
Definition simulation.h:464
double m_dThisIterPotentialPlatformErosion
Total potential platform erosion (all size classes of consolidated sediment) for this iteration (dept...
Definition simulation.h:784
bool bWriteTSFiles(void)
Write the results for this timestep to the time series CSV files.
double m_dDeanProfileStartAboveSWL
Berm height i.e. height above SWL of start of depositional Dean profile.
Definition simulation.h:907
int m_nSavGolCoastPoly
The order of the coastline profile smoothing polynomial if Savitsky-Golay smoothing is used (usually ...
Definition simulation.h:446
void WritePolygonPotentialErosion(int const)
Writes to the log file a table showing per-polygon potential erosion of all size classes of unconsoli...
bool m_bSeaDepthSave
Save sea depth raster GIS files?
Definition simulation.h:92
bool m_bWaveAngleAndHeightSave
Save wave angle and wave height raster GIS files?
Definition simulation.h:110
ofstream BeachDepositionTSStream
Beach sediment deposition time series file output stream.
double m_dNotchBaseBelowSWL
Notch base below SWL (m)
Definition simulation.h:865
vector< double > dVSmoothProfileSlope(vector< double > *) const
Does running-mean smoothing of the slope of a coastline-normal profile.
string m_strInitialBasementDEMFile
Name of initial basement DEM file.
int m_nRunUpEquation
The run-up equation used TODO 007.
Definition simulation.h:533
long double m_ldGTotCliffCollapseFine
All-simulation total of fine sediment from cliff collapse (m)
void WritePolygonActualMovement(int const, vector< vector< int > > const &)
Writes to the log file a table showing per-polygon actual movement of unconsolidated beach sediment.
bool m_bWorldFile
Write a GIS World file?
Definition simulation.h:359
long double m_ldGTotCliffTalusCoarseDeposition
All-simulation total of coarse sediment deposited as talus following cliff collapse (m)
int nDoSedimentInputEvent(int const)
Do a sediment input event.
static string strTrimRight(string const *)
Trims whitespace from the right side of a string, does not change the original string.
Definition utils.cpp:2176
vector< double > m_VdDeepWaterWaveStationX
X coordinate (grid CRS) for deep water wave station.
void AnnounceReadSedimentEventInputValuesGIS(void) const
Tells the user that we are now reading the sediment input events GIS file.
Definition utils.cpp:479
double dGetD50Fine(void) const
Returns the global d50 value for fine sediment.
string m_strLogFile
Name of output log file.
double m_dDepthOverDBMax
Maximum value of deoth over DB, is used in erosion potential look-up function.
Definition simulation.h:844
long double m_ldGTotCliffCollapseCoarse
All-simulation total of coarse sediment from cliff collapse (m)
double m_dFineErodibility
The relative erodibility (0- 1) of fine unconsolidated beach sediment.
Definition simulation.h:742
double m_dThisiterUnconsFineInput
Depth (m) of fine unconsolidated sediment added, at this iteration.
Definition simulation.h:922
time_t m_tSysStartTime
System start-simulation time.
double m_dBreakingWaveHeightDepthRatio
Breaking wave height-to-depth ratio.
Definition simulation.h:706
void AnnounceReadFloodLocationGIS(void) const
Tells the user that we are now reading the flood location GIS file.
Definition utils.cpp:493
int nDoUnconsDepositionOnPolygon(int const, CGeomCoastPolygon *, int const, double, double &)
Deposits unconsolidated beach sediment (sand or coarse) on the cells within a polygon....
bool bCheckVectorGISOutputFormat(void)
Checks whether the selected vector OGR driver supports file creation etc.
bool m_bHaveWaveStationData
Do we have wave station data?
Definition simulation.h:365
double m_dSeaWaterDensity
Density of sea water in kg/m**3.
Definition simulation.h:655
int nDoSimulation(int, char const *[])
Runs the simulation.
int nReadCShoreOutput(int const, string const *, int const, int const, vector< double > const *, vector< double > *)
void StartClock(void)
Starts the clock ticking.
Definition utils.cpp:185
long double m_ldGTotSuspendedSediment
All-simulation total of suspended sediment (m)
vector< double > m_VdSedimentInputLocationX
X coordinate (grid CRS) for sediment input event.
int nCreateAllPolygons(void)
Create polygons, and mark the polygon boundaries on the raster grid.
double m_dSandErodibilityNormalized
Relative erodibility of sand unconsolidated beach sediment, normalized.
Definition simulation.h:754
double m_dR
Coast platform resistance to erosion R, see Walkden & Hall, 2011.
Definition simulation.h:724
string m_strGDALLProjection
GDAL projection string for the initial landform class raster file.
bool m_bActiveZoneSave
Save active zone raster GIS files?
Definition simulation.h:203
void AnnounceReadInitialSuspSedGIS(void) const
Tells the user that we are now reading the initial suspended sediment depth GIS file.
Definition utils.cpp:507
ofstream StillWaterLevelTSStream
SWL time series file output stream.
void DoEndOfRunDeletes(void)
Do end-of-run memory clearance.
Definition utils.cpp:2772
vector< string > m_VstrInitialSandUnconsSedimentFile
The name of the initial sand-sized unconsolidated sediment GIS file.
double m_dThisIterActualPlatformErosionSandCons
Total actual platform erosion (sand consolidated sediment) for this iteration (depth in m)
Definition simulation.h:790
bool m_bOmitSearchWestEdge
Omit the west edge of the grid from coast-end searches?
Definition simulation.h:332
string m_strGDALIHProjection
GDAL projection string for the initial intervention height raster file.
void CheckForIntersectingProfiles(void)
Checks all coastline-normal profiles for intersection, and modifies those that intersect.
int nTraceFloodCoastLine(unsigned int const, int const, int const, vector< bool > *, vector< CGeom2DIPoint > const *)
Traces a coastline (which is defined to be just above still water level) on the grid using the 'wall ...
static int nDoTimeUnits(string const *)
This finds time units in a string.
Definition utils.cpp:337
string m_strGDALICDriverDesc
GDAL description of the initial intervention class raster file.
bool m_bSedimentInputAtCoast
Do we have sediment input at the coast?
Definition simulation.h:374
static void AnnounceReadRasterFiles(void)
Now reading raster GIS files.
Definition utils.cpp:407
double m_dCliffErosionResistance
Resistance of cliff to notch erosion.
Definition simulation.h:859
vector< string > m_VstrGDALIUSProjection
GDAL projection for the initial unconsolidated sand sediment GIS data.
double m_dThisIterBeachDepositionCoarse
Total beach deposition (coarse unconsolidated sediment) for this iteration (depth in m)
Definition simulation.h:811
vector< double > m_VdThisIterDeepWaterWaveStationPeriod
This-iteration wave period at deep water wave station.
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:71
double m_dD50Sand
The D50 for sand sediment.
Definition simulation.h:730
void CalcD50AndFillWaveCalcHoles(void)
Calculates an average d50 for each polygon. Also fills in 'holes' in active zone and wave calcs i....
long double m_ldGTotCliffTalusSandDeposition
All-simulation total of sand sediment deposited as talus following cliff collapse (m)
double m_dCellDiagonal
Length of a cell's diagonal (in external CRS units)
Definition simulation.h:619
string m_strCMEDir
The CME folder.
vector< int > m_VnProfileToSave
The numbers of the profiles which are to be saved.
bool m_bDeepWaterWaveHeightSave
Save deep water wave height raster GIS files?
Definition simulation.h:236
bool m_bMeanWaveEnergySave
Save mean wave energy raster GIS files?
Definition simulation.h:122
double m_dKLS
Transport parameter KLS in the CERC equation.
Definition simulation.h:760
ofstream BeachSedimentNetChangeTSStream
Beach sediment net change (erosion - deposition) time series file output stream.
double m_dInvCellSide
Inverse of m_dCellSide.
Definition simulation.h:622
string m_strOGRDWWVDriverDesc
int m_nWavePropagationModel
The wave propagation model used. Possible values are WAVE_MODEL_CSHORE and WAVE_MODEL_COVE.
Definition simulation.h:506
long double m_ldGTotCliffCollapseSandErodedDuringDeposition
All-simulation total of sand sediment eroded during talus deposition following cliff collapse (m)
double m_dAllCellsDeepWaterWaveAngle
Deep water wave angle for all sea cells.
Definition simulation.h:712
static bool bOnOrOffShoreAndUpOrDownCoast(double const, double const, int const, bool &)
Determines whether the wave orientation at this point on a coast is onshore or offshore,...
vector< string > m_VstrGDALICFProjection
GDAL projection for the initial consolidated fine sediment GIS data.
string m_strOGRFloodDriverCode
GDAL code for the flood input locations point or vector file.
double m_dSimElapsed
Time simulated so far, in hours.
Definition simulation.h:634
vector< string > m_VstrGDALICCProjection
GDAL projection for the initial consolidated coarse sediment GIS data.
static double dGetTimeMultiplier(string const *)
Given a string containing time units, this returns the appropriate multiplier.
Definition utils.cpp:263
bool m_bCoastCurvatureSave
Save coastline-curvature vector GIS files?
Definition simulation.h:260
int nGetGridXMax(void) const
Returns the cell size.
int m_nDeepWaterWaveDataNumTimeSteps
The duration of data for deep water waves, expressed as a number of time steps.
Definition simulation.h:527
bool bWriteProfileData(int const, CGeomProfile const *, int const, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< CGeom2DIPoint > *const, vector< double > const *) const
Writes values for a single profile, for checking purposes.
void AnnounceReadBasementDEM(void) const
Tells the user that we are now reading the DEM file.
Definition utils.cpp:377
int nConvertMetresToNumCells(double const) const
Given a length in m, this returns the rounded equivalent number of cells.
long double m_ldGTotCoarseSedimentInput
All-simulation total of coarse sediment input (m)
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 bWriteParProfileData(int const, int const, int const, int const, int const, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< double > const *, vector< CGeom2DIPoint > *const, vector< double > const *) const
Writes values for a single parallel profile, for checking purposes.
double m_dThisIterUnconsCoarseCliffDeposition
This-iteration total of coarse unconsolidated sediment deposited due to cliff collapse (m^3)
Definition simulation.h:901
bool m_bGDALCanWriteInt32
Is the selected GDAL output file format capable of writing 32-bit integers to files?
Definition simulation.h:353
bool m_bRasterWaveFloodLineSave
Are we saving the raster wave flood line? TODO 007.
Definition simulation.h:404
bool m_bBreakingWaveHeightSave
Save breaking wave height raster GIS files?
Definition simulation.h:125
double m_dThisIterTotSeaDepth
Total sea depth (m) for this iteration.
Definition simulation.h:781
double m_dTotPotentialPlatformErosionOnProfiles
Total potential platform erosion on profiles.
Definition simulation.h:847
int m_nCoastMin
Minimum valid coast legth when searching for coass, actualli is tMin(m_nXGridSize,...
Definition simulation.h:473
string m_strMailAddress
An email addresx to which to send end-of-simulation messages.
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:637
bool m_bPolygonNodeSave
Save polygon node vector GIS files?
Definition simulation.h:263
bool m_bOmitSearchNorthEdge
Omit the north edge of the grid from coast-end searches?
Definition simulation.h:326
long double m_ldGTotCoarseActualPlatformErosion
All-simulation total of coarse sediment actual platform erosion (m)
Definition simulation.h:990
long double m_ldGTotFineActualPlatformErosion
All-simulation total of fine sediment actual platform erosion (m)
Definition simulation.h:984
double m_dThisIterUnconsSandCliffDeposition
This-iteration total of sand unconsolidated sediment deposited due to cliff collapse (m^3)
Definition simulation.h:898
vector< int > m_VnDeepWaterWaveStationID
ID for deep water wave station, this corresponds with the ID in the wave time series file.
long double m_ldGTotSandSedimentInput
All-simulation total of sand sediment input (m)
string strListTSFiles(void) const
Return a space-separated string containing the names of the time series output files.
Definition utils.cpp:992
bool m_bFloodSetupSurgeTSSave
Save the flood setup surge time series file? TODO 007 Does this work correctly?
Definition simulation.h:308
bool m_bTotalPotentialPlatformErosionSave
Save total potential shore platform erosion raster GIS files?
Definition simulation.h:137
string m_strGDALLDriverDesc
GDAL description of the initial landform class raster file.
void SetRasterFileCreationDefaults(void)
Sets per-driver defaults for raster files created using GDAL.
bool m_bSetupSurgeFloodMaskSave
Are we saving the setup surge flood mask? TODO 007.
Definition simulation.h:398
bool m_bWaveEnergySinceCollapseSave
Save wave energy since cliff collapse raster GIS files?
Definition simulation.h:119
double m_dNorthWestXExtCRS
The north-west x coordinate, in the external coordinate reference system (CRS)
Definition simulation.h:598
double m_dD50Coarse
The D50 for coarse sediment.
Definition simulation.h:733
vector< CGeom2DIPoint > m_VEdgeCell
Edge cells.
static bool bCheckForIntersection(CGeomProfile *const, CGeomProfile *const, int &, int &, double &, double &, double &, double &)
Checks all line segments of a pair of coastline-normal profiles for intersection. If the lines inters...
int nInterpolateAllDeepWaterWaveValues(void)
If the user supplies multiple deep water wave height and angle values, this routine interplates these...
static void AnnounceStart(void)
Tells the user that we have started the simulation.
Definition utils.cpp:177
static int nUpdateIntervention(void)
Check to see if we have a new intervention in place (not yet implemented)
vector< int > m_VnFloodLocationID
ID for flood location.
int m_nCoastCurvatureInterval
Coast curvature interval is a length, measured in coastline points.
Definition simulation.h:455
vector< unsigned long > m_VulProfileTimestep
Timesteps at which to save profiles.
int nCalcPotentialPlatformErosionOnProfile(int const, CGeomProfile *)
Calculates potential (i.e. unconstrained by available sediment) erosional lowering of the shore platf...
ofstream CliffCollapseErosionTSStream
Cliff collapse erosion time series file output stream.
double m_dThisIterPotentialSedLostBeachErosion
Total unconsolidated sediment from beach erosion (all size classes) lost from the grid this iteration...
Definition simulation.h:817
double m_dSouthEastYExtCRS
The south-east y coordinate, in the external coordinate reference system (CRS)
Definition simulation.h:607
bool m_bPotentialBeachErosionSave
Save potential beach (unconsolidated sediment) erosion raster GIS files?
Definition simulation.h:143
normal_distribution< double > m_dUnitNormalDist
c++11 unit normal distribution (mean = 0, stdev = 1)
static bool bParseDate(string const *, int &, int &, int &)
Parses a date string into days, months, and years, and checks each of them.
Definition utils.cpp:2480
string m_strFloodLocationShapefile
The name of the flood loction events shape file.
long double m_ldGTotCliffCollapseFineErodedDuringDeposition
All-simulation total of fine sediment eroded during talus deposition following cliff collapse (m)
~CSimulation(void)
The CSimulation destructor.
long double m_ldGTotCoarseDepositionDiff
All-simulation total of shortfall in unconsolidated coarse sediment deposition (m,...
int nSetAllCoastpointDeepWaterWaveValues(void)
Give every coast point a value for deep water wave height and direction TODO 005 This may not be real...
double m_dCoastNormalAvgSpacing
Average spacing of the cost-normal profiles, in m.
Definition simulation.h:775
double m_dStartIterConsFineAllCells
Depth (m) of fine consolidated sediment at the start of the simulation, all cells (both inside and ou...
Definition simulation.h:946
long double m_ldGTotActualCoarseLostBeachErosion
All-simulation total of coarse sediment lost via beach (unconsolidated) sediment movement (m)
bool m_bGISSaveDigitsSequential
Are the GIS save digits (which are part of each GIS file name) sequential, or are they the iteration ...
Definition simulation.h:422
static CGeom2DIPoint PtiFollowWaveAngle(CGeom2DIPoint const *, double const, double &)
Given a cell and a wave orientation, finds the 'upwave' cell.
vector< double > m_VdSavGolFCRWCoast
Savitzky-Golay filter coefficients for the coastline vector(s)
vector< double > m_VdSavGolFCGeomProfile
Savitzky-Golay filter coefficients for the profile vectors.
int nUpdateGrid(void)
Update all cells in the raster raster grid and do some per-timestep accounting.
int m_nSimStartMonth
Start date of the simulation (month)
Definition simulation.h:521
int nGetCoastNormalEndPoint(int const, int const, int const, CGeom2DPoint const *, double const, CGeom2DPoint *, CGeom2DIPoint *, bool const)
Finds the end point of a coastline-normal line, given the start point on the vector coastline....
int nAssignLandformsForAllCells(void)
Each timestep, classify landforms for cells that are not on the coastline.
string m_strGDALIWProjection
GDAL projection string for the initial water depth raster file.
bool m_bSuspSedSave
Save suspended sediment raster GIS files?
Definition simulation.h:173
double m_dMinCliffTalusHeightFrac
Minimum height of the landward end of cliff collapse talus, as a fraction of cliff elevation.
Definition simulation.h:877
double m_dThisIterBeachErosionSand
Total actual beach erosion (sand unconsolidated sediment) for this iteration (depth in m)
Definition simulation.h:802
static double dAngleSubtended(CGeom2DIPoint const *, CGeom2DIPoint const *, CGeom2DIPoint const *)
Returns the signed angle BAC (in radians) subtended between three CGeom2DIPoints B A C....
vector< string > m_VstrGDALIUSDriverDesc
GDAL driver description for the initial unconsolidated sand sediment GIS data.
string m_strDurationUnits
The duration units for this simulation.
void TruncateProfileAndAppendNew(int const, CGeomProfile *, int const, vector< CGeom2DPoint > const *, vector< vector< pair< int, int > > > const *)
Truncate a profile at the point of intersection, and do the same for all its co-incident profiles.
double m_dAccumulatedSeaLevelChange
If long-term SWL changes, the total change so far since the start of simulation.
Definition simulation.h:673
bool m_bPolygonBoundarySave
Save polygon boundary vector GIS files?
Definition simulation.h:266
bool bElevAboveDeanElev(int const, int const, double const, CRWCellLandform const *)
Return true if the given elevation is higher than the Dean elevation (and other conditions are met),...
void WritePolygonPreExistingSedimentTable(int const)
Writes to the log file a table showing per-polygon pre-existing unconsolidated sediment.
bool m_bStormSurgeSave
Are we saving the storm surge? TODO 007.
Definition simulation.h:392
void DoActualPlatformErosionOnCell(int const, int const)
Calculates actual (constrained by available sediment) erosion of the consolidated shore platform on a...
vector< string > m_VstrInitialFineUnconsSedimentFile
The name of the initial fine-sized unconsolidated sediment GIS file.
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:134
double m_dThisIterFineSedimentToSuspension
Total fine unconsolidated sediment in suspension for this iteration (depth in m)
Definition simulation.h:814
int nDoAllWaveEnergyToCoastLandforms(void)
Update accumulated wave energy in coastal landform objects.
vector< CSedInputEvent * > m_pVSedInputEvent
Sediment input events.
int nInsertPointIntoProfilesIfNeededThenUpdate(int const, CGeomProfile *, double const, double const, int const, CGeomProfile *, int const, bool const)
Inserts an intersection point into the profile that is to be retained, if that point is not already p...
double m_dThisIterBeachErosionFine
Total actual beach erosion (fine unconsolidated sediment) for this iteration (depth in m)
Definition simulation.h:799
void KeepWithinValidGrid(int, int, int &, int &) const
Given two points in the grid CRS (the points assumed not to be coincident), this routine modifies the...
void DoAllPotentialBeachErosion(void)
Uses either the CERC equation or the Kamphuis (1990) equation to calculate potential (unconstrained) ...
void FillPotentialPlatformErosionHoles(void)
Fills in 'holes' in the potential platform erosion i.e. orphan cells which get omitted because of rou...
bool bSurroundedByDriftCells(int const, int const)
Returns true if this cell has four drift cells surrounding it.
double m_dThisIterCliffCollapseSandErodedDuringDeposition
Total sand sediment eroded during Dean profile deposition of talus following cliff collapse (depth in...
Definition simulation.h:832
int m_nSimStartMin
Start time of the simulation (minutes)
Definition simulation.h:512
bool bSaveAllRasterGISFiles(void)
The bSaveAllRasterGISFiles member function saves the raster GIS files using values from the RasterGri...
string m_strOGRDWWVGeometry
GDAL geometry for the deep water wave stations vector file.
unsigned long m_ulThisIterNumActualBeachErosionCells
The number of grid cells on which actual beach (unconsolidated sediment) erosion occurs,...
Definition simulation.h:580
bool bReadIniFile(void)
The bReadIniFile member function reads the initialization file.
bool m_bHaveFineSediment
Does this simulation consider fine-sized sediment?
Definition simulation.h:71
double m_dCliffDepositionA
Scale parameter A for cliff deposition (m^(1/3)), may be zero for auto-calculation.
Definition simulation.h:868
string m_strDataPathName
Folder in which the CME data file is found.
string m_strOutPath
Path for all output files.
bool m_bBeachErosionTSSave
Save the beach (unconsolidated sediment) erosion time series file?
Definition simulation.h:296
int m_nYMaxBoundingBox
The maximum y value of the bounding box.
Definition simulation.h:503
int nDoCliffCollapse(int const, CRWCliff *, double &, double &, double &, double &, double &)
Simulates cliff collapse on a single cliff object, which happens when when a notch (incised into a co...
vector< string > m_VstrGDALICFDriverDesc
GDAL driver description for the initial consolidated fine sediment GIS data.
vector< double > m_VdSedimentInputLocationY
X coordinate (grid CRS) for sediment input event.
static string strGetBuild(void)
Returns the date and time on which the program was compiled.
Definition utils.cpp:1584
string m_strTideDataFile
Name of tide data file.
vector< string > m_VstrGDALIUFDriverDesc
GDAL driver description for the initial unconsolidated fine sediment GIS data.
vector< bool > m_bUnconsChangedThisIter
One element per layer: has the consolidated sediment of this layer been changed during this iteration...
bool m_bTotalPotentialBeachErosionSave
Save total potential beach (unconsolidated sediment) erosion raster GIS files?
Definition simulation.h:149
unsigned long m_ulThisIterNumPotentialPlatformErosionCells
The number of grid cells on which potential platform erosion occurs, for this iteration.
Definition simulation.h:571
string m_strGDALIWDataType
GDAL data type for the initial water depth raster file.
vector< string > m_VstrGDALIUSDataType
GDAL data type for the initial unconsolidated sand sediment GIS data.
double m_dStartIterUnconsSandAllCells
Depth (m) of sand unconsolidated sediment at the start of the simulation, all cells (both inside and ...
Definition simulation.h:940
double m_dMaxSWL
Maximum still water level.
Definition simulation.h:679
long double m_ldGTotActualSandLostBeachErosion
All-simulation total of sand sediment lost via beach (unconsolidated) sediment movement (m)
Definition simulation.h:999
double m_dDeltaSWLPerTimestep
If long-term SWL changes, the increment per timestep.
Definition simulation.h:664
void AllPolygonsUpdateStoredUncons(int const)
Before simulating beach erosion, update the per-polygon values of pre-existing unconsolidated sedimen...
int m_nGISSave
The save number for GIS files (can be sequential, or the iteration number)
Definition simulation.h:461
string m_strInterventionClassFile
Name of intervention class file.
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...
long double m_ldGTotActualFineLostBeachErosion
All-simulation total of fine sediment lost via beach (unconsolidated) sediment movement (m)
Definition simulation.h:996
string m_strGDALBasementDEMDataType
GDAL data type for the basement DEM raster file.
bool m_bFloodSWLSetupLine
Are we saving the flood still water level setup line? TODO 007.
Definition simulation.h:413
bool m_bSeaMaskSave
Save sea mask raster GIS files?
Definition simulation.h:224
void CalcTime(double const)
Calculates and displays time elapsed in terms of CPU time and real time, also calculates time per tim...
Definition utils.cpp:1401
double m_dThisIterCliffCollapseErosionFineCons
This-iteration total of fine consolidated sediment produced by cliff collapse (m^3)
Definition simulation.h:889
bool m_bInterventionClassSave
Save intervention class raster GIS files?
Definition simulation.h:167
long double m_ldGTotCoarseSedLostCliffCollapse
All-simulation total of coarse sediment lost via cliff collapse (m)
long double m_ldGTotSandActualPlatformErosion
All-simulation total of sand sediment actual platform erosion (m)
Definition simulation.h:987
int m_nSimStartYear
Start date of the simulation (year)
Definition simulation.h:524
CGeom2DIPoint PtiFindClosestCoastPoint(int const, int const)
Finds the closest point on any coastline to a given point.
bool m_bTotalActualBeachErosionSave
Save total actual (supply-limited) beach (unconsolidated sediment) erosion raster GIS files?
Definition simulation.h:152
void WritePolygonShareTable(int const)
Writes to the log file a table showing polygon to polygon shares of unconsolidated sediment transport...
bool m_bRasterCoastlineSave
Save raster coastline GIS files?
Definition simulation.h:197
string m_strDeepWaterWavesInputFile
The name of the deep water wave stations time series file.
static string pstrChangeToBackslash(string const *)
Changes all forward slashes in the input string to backslashes, leaving the original unchanged.
Definition utils.cpp:2143
static CGeom2DPoint PtFindPointInPolygon(vector< CGeom2DPoint > const *, int const)
Finds a point in a polygon: is guaranteed to succeed, as every strictly closed polygon has at least o...
long double m_ldGTotSandSedLostCliffCollapse
All-simulation total of sand sediment lost via cliff collapse (m)
double m_dUnconsSandNotDepositedLastIter
Depth of unconsolidated sand sediment that could not be deposited during the last iteration,...
Definition simulation.h:973
double m_dThisIterDiffWaveSetupWaterLevel
TODO 007 Info needed.
Definition simulation.h:685
static CGeom2DPoint PtGetPerpendicular(CGeom2DPoint const *, CGeom2DPoint const *, double const, int const)
Returns a CGeom2DPoint which is the 'other' point of a two-point vector passing through PtStart,...
string m_strGDALICProjection
GDAL projection string for the initial intervention class raster file.
bool m_bInterventionHeightSave
Save intervention height raster GIS files?
Definition simulation.h:170
double m_dTotalCoarseConsInPolygons
Total coarse consolidated sediment in all polygons, before polygon-to-polygon movement (only cells in...
Definition simulation.h:970
int m_nCoastCurvatureMovingWindowSize
Definition simulation.h:538
ofstream SeaAreaTSStream
Sea area time series file output stream.
unsigned long m_ulThisIterNumSeaCells
The number of grid cells which are marked as sea, for this iteration.
Definition simulation.h:565
int nDoUnconsErosionOnPolygon(int const, CGeomCoastPolygon *, int const, double const, double &)
Erodes unconsolidated beach sediment of one texture class on the cells within a polygon....
void AnnounceReadIHGIS(void) const
Tells the user that we are now reading the Intervention height GIS file.
Definition utils.cpp:451
double m_dTotalSandUnconsInPolygons
Total sand unconsolidated sediment in all polygons, before polygon-to-polygon movement (only cells in...
Definition simulation.h:958
void AnnounceReadInitialCoarseConsSedGIS(int const) const
Tells the user that we are now reading the initial coarse consolidated sediment depth GIS file.
Definition utils.cpp:585
bool m_bRiverineFlooding
Are we doing flooding? TODO 007.
Definition simulation.h:386
void RasterizePolygonJoiningLine(CGeom2DIPoint const *, CGeom2DIPoint const *, int const)
Puts a polygon 'joining line' (the line which is the seaward boundary of the polygon,...
long m_lGDALMinCanWrite
The minimum integer value which GDAL can write, can be zero, INT16_MIN, INT32_MIN.
Definition simulation.h:550
int nCreateAllProfiles(void)
Create coastline-normal profiles for all coastlines. The first profiles are created 'around' the most...
double m_dStartIterUnconsFineAllCells
Depth (m) of fine unconsolidated sediment at the start of the simulation, all cells (both inside and ...
Definition simulation.h:937
static string strGetErrorText(int const)
Returns an error message given an error code.
Definition utils.cpp:1818
bool m_bSandConsSedSave
Save sand consolidated sediment raster GIS files?
Definition simulation.h:191
CGeomCoastPolygon * pGetPolygon(int const) const
Returns a pointer to a coast polygon, in down-coast sequence.
Definition utils.cpp:2755
bool m_bSedimentInputEventSave
Save sediment inut data?
Definition simulation.h:380
static double dTriangleAreax2(CGeom2DPoint const *, CGeom2DPoint const *, CGeom2DPoint const *)
Returns twice the signed area of a triangle.
bool m_bHaveCoarseSediment
Does this simulation consider coarse-sized sediment?
Definition simulation.h:77
double m_dRegularSaveInterval
The interval between regular saves, in hours.
Definition simulation.h:640
long double m_ldGTotActualFineBeachErosion
All-simulation total of fine sediment eroded during beach (unconsolidated sediment) movement (m)
CGeomLine LSmoothCoastSavitzkyGolay(CGeomLine *, int const, int const) const
Does smoothing of a CGeomLine coastline vector (is in external CRS coordinates) using a Savitzky-Gola...
bool m_bPolygonUnconsSedGainOrLossSave
Save polygon unconsolidated sediment gain or loss raster GIS files?
Definition simulation.h:245
void MarkPolygonCells(void)
Marks cells of the raster grid that are within each coastal polygon. The flood fill code used here is...
string m_strGDALISSDriverDesc
GDAL description for the initial suspended sediment raster file.
double m_dTimeStep
The length of an iteration (a time step) in hours.
Definition simulation.h:631
void TruncateOneProfileRetainOtherProfile(int const, CGeomProfile *, CGeomProfile *, double const, double const, int const, int const, bool const)
Truncates one intersecting profile at the point of intersection, and retains the other profile.
static void AnnounceAllocateMemory(void)
Tells the user that we are now allocating memory.
Definition utils.cpp:390
unsigned long m_ulTotPotentialPlatformErosionOnProfiles
The number of cells on which on-profile average potential shore platform erosion occurs.
Definition simulation.h:586
bool m_bCliffCollapseDepositionTSSave
Save the cliff collapse deposition time series file?
Definition simulation.h:290
int m_nXMinBoundingBox
The minimum x value of the bounding box.
Definition simulation.h:494
double dGetD50Sand(void) const
Returns the global d50 value for sand sediment.
string m_strGDALRasterOutputDriverLongname
GDAL raster output driver long name.
vector< string > m_VstrGDALIUCDriverCode
GDAL driver code for the initial unconsolidated coarse sediment GIS data.
ofstream PlatformErosionTSStream
Shore platform erosion time series file output stream.
void DoEndOfTimestepTotals(void)
Update and print totals at the end of each timestep.
double m_dCellArea
Area of a cell (in external CRS units)
Definition simulation.h:616
bool m_bDeepWaterWaveAngleAndHeightSave
Save deep water wave angle and wave height raster GIS files?
Definition simulation.h:116
unsigned long m_ulTotPotentialPlatformErosionBetweenProfiles
The number of cells on which between-profile average potential shore platform erosion occurs.
Definition simulation.h:589
double dExtCRSYToGridY(double const) const
Transforms a Y-axis ordinate in the external CRS to the equivalent Y-axis ordinate in the raster grid...
double m_dCoastNormalRandSpacingFactor
Random factor for spacing of along-coast normals.
Definition simulation.h:904
double m_dMaxUserInputWaveHeight
Maximum deep water wave height.
Definition simulation.h:718
vector< double > m_VdSliceElev
Elevations for raster slice output.
void AnnounceReadInitialFineConsSedGIS(int const) const
Tells the user that we are now reading the initial fine consolidated sediment depth GIS file.
Definition utils.cpp:559
bool m_bBeachProtectionSave
Save beach protection raster GIS files>
Definition simulation.h:128
static string strDispTime(double const, bool const, bool const)
strDispTime returns a string formatted as h:mm:ss, given a parameter in seconds, with rounding and fr...
Definition utils.cpp:1522
string m_strOGRSedInputGeometry
GDAL geometry for the sediment input event locations vector file.
string m_strOGRSedInputDriverCode
GDAL code for the sediment input event locations vector file.
double m_dFinalSWL
The end-of-simulation still water (m), is same as m_dOrigSWL unless SWL changes.
Definition simulation.h:661
bool m_bDoBeachSedimentTransport
Simulate unconsolidated sediment (beach) transport?
Definition simulation.h:344
int nGetThisProfileElevationsForCShore(int const, CGeomProfile *, int const, vector< double > *, vector< double > *, vector< double > *)
Get profile horizontal distance and bottom elevation vectors in CShore units.
unsigned long m_ulThisIterNumPotentialBeachErosionCells
The number of grid cells on which potential beach (unconsolidated sediment) erosion occurs,...
Definition simulation.h:577
bool m_bFineConsSedSave
Save fine consolidated sediment raster GIS files?
Definition simulation.h:188
string m_strGDALIWDriverDesc
GDAL description for the initial water depth raster file.
bool m_bShadowDowndriftBoundarySave
Save wave shadow downdrift boundary vector GIS files?
Definition simulation.h:275
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,...
int nDoParallelProfileUnconsErosion(CGeomCoastPolygon *, int const, int const, int const, int const, int const, int const, vector< CGeom2DIPoint > const *, vector< double > const *, double &, double &, double &)
This routine erodes unconsolidated beach sediment (either fine, sand, or coarse) on a parallel profil...
int m_nCoastSmooth
Which method to use for coast smoothing.
Definition simulation.h:440
bool m_bDeepWaterWavePeriodSave
Save deep water wave period raster GIS files?
Definition simulation.h:239
string m_strGDALBasementDEMDriverDesc
GDAL description of the basement DEM raster file type.
void ModifyBreakingWavePropertiesWithinShadowZoneToCoastline(int const, int const)
Modifies the wave breaking properties at coastline points of profiles within the shadow zone.
string m_strInitialLandformFile
Name of initial landform file.
double m_dC_0
Deep water wave speed (m/s)
Definition simulation.h:697
string m_strInterventionHeightFile
Name of intervention height file.
double m_dThisIterCliffCollapseErosionSandCons
This-iteration total of sand consolidated sediment produced by cliff collapse (m^3)
Definition simulation.h:892
bool m_bBeachSedimentChangeNetTSSave
Save the beach (unconsolidated sediment) net change time series file?
Definition simulation.h:302
int nCreateProfile(int const, int const, int const, int const, int &, bool const, CGeom2DIPoint const *)
Creates a single coastline-normal profile (which may be an intervention profile)
vector< string > m_VstrGDALICSDataType
GDAL data type for the initial consolidated sand sediment GIS data.
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:92
double m_dThisIterBeachDepositionSand
Total beach deposition (sand unconsolidated sediment) for this iteration (depth in m)
Definition simulation.h:808
int nReadTideDataFile(void)
Reads the tide time series data.
bool m_bCoarseConsSedSave
Save coarse consolidated sediment raster GIS files?
Definition simulation.h:194
string m_strOGRFloodDriverDesc
bool m_bSeaAreaTSSave
Save the sea area time series file?
Definition simulation.h:278
bool m_bScaleRasterOutput
Scale raster output?
Definition simulation.h:356
double dGetThisIterSWL(void) const
Returns this timestep's still water level.
vector< string > m_VstrGDALICSProjection
GDAL dprojection for the initial consolidated sand sediment GIS data.
void CreateRasterizedProfile(int const, CGeomProfile *, CGeom2DIPoint const *, CGeom2DIPoint const *, vector< CGeom2DIPoint > *, vector< bool > *, bool &, bool &, bool &, bool &)
Given a pointer to a coastline-normal profile, returns an output vector of cells which are 'under' ev...
void CalcProcessStats(void)
This calculates and displays process statistics.
Definition utils.cpp:1636
void DoShadowZoneAndDownDriftZone(int const, int const, int const, int const)
Traverse the shadow zone, changing wave orientation and height, and the down-drift zone,...
static double dCalcCurvature(int const, CGeom2DPoint const *, CGeom2DPoint const *, CGeom2DPoint const *)
Calculates signed Menger curvature (https://en.wikipedia.org/wiki/Menger_curvature#Definition) from t...
int nGetCoastPolygonSize(void) const
Returns the size of the coast polygon vector.
Definition utils.cpp:2747
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:80
double m_dThisIterLeftGridUnconsCoarse
Total coarse unconsolidated sediment lost from the grid this iteration (depth in m)
Definition simulation.h:826
void DoCoastCurvature(int const, int const)
Calculates both detailed and smoothed curvature for every point on a coastline.
double m_dThisIterCliffCollapseErosionSandUncons
This-iteration total of sand unconsolidated sediment produced by cliff collapse (m^3)
Definition simulation.h:883
string m_strOGRFloodDataType
GDAL data type for the flood input locations point or vector file.
double m_dD50Fine
The D50 for fine sediment.
Definition simulation.h:727
bool m_bOmitSearchSouthEdge
Omit the south edge of the grid from coast-end searches?
Definition simulation.h:329
string m_strGDALRasterOutputDriverExtension
GDAL raster output driver file extension.
bool m_bBeachMaskSave
Save beach mask raster GIS files?
Definition simulation.h:227
string m_strGDALBasementDEMProjection
GDAL projection string for the basement DEM raster file.
bool bOpenLogFile(void)
Opens the log file.
Definition utils.cpp:354
double m_dTotalSandConsInPolygons
Total sand consolidated sediment in all polygons, before polygon-to-polygon movement (only cells in p...
Definition simulation.h:967
static void AppendEnsureNoGap(vector< CGeom2DIPoint > *, CGeom2DIPoint const *)
Appends a CGeom2DIPoint to a vector<CGeom2DIPoint>, making sure that the new end point touches the pr...
Definition utils.cpp:2317
unsigned long m_ulIter
The number of the current iteration (time step)
Definition simulation.h:553
void AnnounceLicence(void)
Tells the user about the licence.
Definition utils.cpp:242
bool m_bAvgSuspSedSave
Save average suspended sediment raster GIS files?
Definition simulation.h:176
void WritePolygonUnsortedSequence(int const, vector< vector< int > > &)
Writes to the log file a table showing the unsorted sequence of polygon processing.
double m_dBeachSedimentDensity
The density of unconsolidated beach sediment (kg/m**3)
Definition simulation.h:736
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:62
vector< string > m_VstrInitialCoarseConsSedimentFile
The name of the initial coarse-sized consolidated sediment GIS file.
string m_strOGRDWWVDataType
GDAL data type for the deep water wave stations vector file.
int m_nNumPolygonGlobal
Number of global (all coasts) polygons.
Definition simulation.h:482
int nReadRasterBasementDEM(void)
Reads a raster DEM of basement elevation data to the Cell array.
double m_dSimDuration
Duration of simulation, in hours.
Definition simulation.h:628
double m_dTotalFineUnconsInPolygons
Total fine unconsolidated sediment in all polygons, before polygon-to-polygon movement (only cells in...
Definition simulation.h:955
int nTraceAllCoasts(int &)
Locates all the potential coastline start points on the edges of the raster grid, then traces vector ...
string m_strGDALIHDriverCode
GDAL code for the initial intervention height raster file.
long double m_ldGTotSandDepositionDiff
All-simulation total of shortfall in unconsolidated sand sediment deposition (m, not currently used)
string m_strGDALICDriverCode
GDAL code for the initial intervention class raster file.
ofstream FloodSetupSurgeTSStream
Flood setup surge time series file output stream.
double dCalcBeachProtectionFactor(int const, int const, double const)
Calculates the (inverse) beach protection factor as in SCAPE: 0 is fully protected,...
bool bTimeToQuit(void)
Checks to see if the simulation has gone on too long, amongst other things.
Definition utils.cpp:1292
double m_dTotalFineConsInPolygons
Total fine consolidated sediment in all polygons, before polygon-to-polygon movement (only cells in p...
Definition simulation.h:964
double m_dThisiterUnconsSandInput
Depth (m) of sand unconsolidated sediment added, at this iteration.
Definition simulation.h:925
double m_dExtCRSGridArea
The area of the grid (in external CRS units)
Definition simulation.h:610
void ProcessDownDriftCell(int const, int const, int const, double const, int const)
Process a single cell which is in the downdrift zone, changing its wave height.
int nDoAllPropagateWaves(void)
Simulates wave propagation along all coastline-normal profiles.
vector< double > m_VdThisIterDeepWaterWaveStationAngle
This-iteration wave orientation at deep water wave station.
bool m_bOutputProfileData
Output profile data?
Definition simulation.h:317
double m_dThisIterDiffWaveSetupSurgeRunupWaterLevel
TODO 007 Info needed.
Definition simulation.h:691
double m_dCellSide
Length of a cell side (in external CRS units)
Definition simulation.h:613
vector< string > m_VstrGDALIUCDriverDesc
GDAL driver description for the initial unconsolidated coarse sediment GIS data.
bool bIsInterventionCell(int const, int const) const
Returns true if the cell is an intervention.
Definition utils.cpp:2736
double m_dMaxBeachElevAboveSWL
Maximum elevation of beach above SWL (m)
Definition simulation.h:856
long double m_ldGTotActualSandBeachErosion
All-simulation total of sand sediment eroded during beach (unconsolidated sediment) movement (m)
bool m_bRasterNormalSave
Save raster coastline-normal GIS files?
Definition simulation.h:200
double m_dBreakingWaveHeight
The height of breaking waves (m)
Definition simulation.h:694
void AnnounceProgress(void)
Displays information regarding the progress of the simulation.
Definition utils.cpp:1601
bool m_bTotCliffCollapseDepositionSave
Save total cliff collapse deposition raster GIS files?
Definition simulation.h:215
double m_dAllCellsDeepWaterWavePeriod
Deep water wave period for all sea cells.
Definition simulation.h:715
double m_dThisIterCliffCollapseCoarseErodedDuringDeposition
Total coarse sediment eroded during Dean profile deposition of talus following cliff collapse (depth ...
Definition simulation.h:835
int nMarkBoundingBoxEdgeCells(void)
Mark cells which are at the edge of a bounding box which represents the valid part of the grid,...
bool bCreateErosionPotentialLookUp(vector< double > *, vector< double > *, vector< double > *)
Creates a look-up table for erosion potential, given depth over DB.
long double m_ldGTotSandBeachDeposition
All-simulation total of sand sediment deposited during beach (unconsolidated sediment) movement (m)
double m_dUSaveTime[SAVEMAX]
Save time, in hours from the start of the simukation, if we are not saving regularly.
Definition simulation.h:643
vector< int > m_VnSavGolIndexCoast
Savitzky-Golay shift index for the coastline vector(s)
void AnnounceReadLGIS(void) const
Tells the user that we are now reading the Landscape category GIS file.
Definition utils.cpp:423
int m_nCoastNormalAvgSpacing
Average spacing between cost normals, measured in cells.
Definition simulation.h:452
double m_dThisIterTopElevMax
This-iteration highest elevation of DEM.
Definition simulation.h:916
double m_dThisIterActualPlatformErosionFineCons
Total actual platform erosion (fine consolidated sediment) for this iteration (depth in m)
Definition simulation.h:787
vector< double > m_VdFloodLocationY
X coordinate (grid CRS) for total water level flooding.
long double m_ldGTotPotentialBeachErosion
All-simulation total of potential beach erosion (m), all size classes.
double m_dProfileMaxSlope
Maximum slope on costline-normal profiles.
Definition simulation.h:853
int m_nProfileSmoothWindow
The size of the window used for running-mean coast-normal profile smoothing (must be odd)
Definition simulation.h:449
string m_strGDALIHDriverDesc
GDAL description for the initial intervention height raster file.
bool m_bDoCliffCollapse
Simulate cliff collapse?
Definition simulation.h:341
static void AnnounceReadVectorFiles(void)
Now reading vector GIS files.
Definition utils.cpp:415
double m_dThisIterCliffCollapseErosionCoarseCons
This-iteration total of coarse consolidated sediment produced by cliff collapse (m^3)
Definition simulation.h:895
double m_dThisIterLeftGridUnconsSand
Total sand unconsolidated sediment lost from the grid this iteration (depth in m)
Definition simulation.h:823
unsigned long m_ulRandSeed[NRNG]
A seed for each of the NRNG random number generators.
Definition simulation.h:559
double m_dDepthOfClosure
Depth of closure (in m) TODO 007 can be calculated using Hallermeier, R.J. (1978) or Birkemeier (1985...
Definition simulation.h:772
int nCheckAllProfiles(void)
Check all coastline-normal profiles and modify the profiles if they intersect, then mark valid profil...
int nInitGridAndCalcStillWaterLevel(void)
At the beginning of each timestep: clear vector coasts, profiles, and polygons, initialize the raster...
Definition init_grid.cpp:46
vector< string > m_VstrGDALICCDriverDesc
GDAL driver decription for the initial consolidated coarse sediment GIS data.
double m_dOrigSWL
The start-of-simulation still water level (m)
Definition simulation.h:658
double m_dStartIterSuspFineAllCells
Depth (m) of fine suspended sediment at the start of the simulation, all cells (both inside and outsi...
Definition simulation.h:931
string m_strOutFile
Name of main output file.
bool m_bSetupSurgeRunupFloodMaskSave
Are we saving the setup surge runup flood mask? TODO 007.
Definition simulation.h:401
void FloodFillSea(int const, int const)
Flood-fills all sea cells starting from a given cell. The flood fill code used here is adapted from a...
static void CalcDeanProfile(vector< double > *, double const, double const, double const, bool const, int const, double const)
Calculates a Dean equilibrium profile h(y) = A * y^(2/3) where h(y) is the distance below the highest...
Definition utils.cpp:2356
double m_dThisIterDiffTotWaterLevel
TODO 007 Info needed.
Definition simulation.h:682
vector< bool > m_bConsChangedThisIter
One element per layer: has the consolidated sediment of this layer been changed during this iteration...
string m_strGDALIHDataType
GDAL data type for the initial intervention height raster file.
ofstream FineSedSuspensionTSStream
Fine sediment in suspension time series file output stream.
ofstream FloodSetupSurgeRunupTSStream
Flood setup surge runup time series file output stream.
bool m_bWaveSetupSave
Are we saving the wave setup? TODO 007.
Definition simulation.h:389
int nDoPolygonSharedBoundaries(void)
For between-polygon potential sediment routing: find which are the adjacent polygons,...
void AnnounceReadTideData(void) const
Now reading tide data file.
Definition utils.cpp:598
static bool bIsWithinPolygon(CGeom2DPoint const *, vector< CGeom2DPoint > const *)
Determines whether a point is within a polygon: however if the point is exactly on the edge of the po...
vector< double > m_VdTSDeepWaterWaveStationHeight
Time series of wave heights at deep water wave station.
double m_dClkLast
Last value returned by clock()
Definition simulation.h:646
bool m_bShadowZoneCodesSave
Save wave shadow zones raster GIS files?
Definition simulation.h:230
long double m_ldGTotFineSedimentInput
All-simulation total of fine sediment input (m)
double m_dThisIterMeanSWL
The mean still water level for this timestep (does not include tidal changes, but includes any long-t...
Definition simulation.h:670
vector< string > m_VstrGDALIUCDataType
GDAL data type for the initial unconsolidated coarse sediment GIS data.
bool m_bCliffCollapseErosionTSSave
Save the cliff collapse erosion time series file?
Definition simulation.h:287
bool m_bPotentialPlatformErosionSave
Save potential shore platform erosion raster GIS files?
Definition simulation.h:131
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...
long double m_ldGTotPotentialPlatformErosion
All-simulation total of potential platform erosion (m), all size classes.
Definition simulation.h:981
GDALDataType m_GDALWriteFloatDataType
Thw data type used by GDAL for floating point operations, can be GDT_Byte, GDT_Int16,...
Definition simulation.h:544
void InterpolateWaveHeightToCoastPoints(int const)
Linearly interpolates wave properties from profiles to the coastline cells between two profiles.
ofstream BeachErosionTSStream
Beach sediment erosion time series file output stream.
bool m_bGDALCanWriteFloat
Is the selected GDAL output file format capable of writing floating-point values to files?
Definition simulation.h:350
void AnnounceReadInitialCoarseUnconsSedGIS(int const) const
Tells the user that we are now reading the initial coarse unconsolidated sediment depth GIS file.
Definition utils.cpp:546
int nHandleCommandLineParams(int, char const *[])
Handles command-line parameters.
Definition utils.cpp:86
long double m_ldGTotCliffCollapseSand
All-simulation total of sand sediment from cliff collapse (m)
int FindAllInundatedCells(void)
Finds and flags all sea areas which have at least one cell at a grid edge (i.e. does not flag 'inland...
char ** m_papszGDALRasterOptions
Options for GDAL when handling raster files.
Definition simulation.h:425
ofstream OutStream
The main output file stream.
double m_dDurationUnitsMult
Multiplier for duration units, to convert to hours.
Definition simulation.h:595
bool m_bFloodSWLSetupSurgeRunupLine
Are we saving the flood still water level setup surge runup line? TODO 007.
Definition simulation.h:419
void LocateAndCreateProfiles(int const, int &, int &, int const, vector< bool > *, vector< pair< int, double > > const *)
For a single coastline, locate the start points for all coastline-normal profiles (except the grid-ed...
unsigned long m_ulTotTimestep
The target number of iterations.
Definition simulation.h:556
bool m_bOutputParallelProfileData
Output parallel profile data?
Definition simulation.h:320
int nLandformToGrid(int const, int const)
At the end of each timestep, this routine stores the attributes from a single coastal landform object...
double m_dKamphuis
Transport parameter for the Kamphuis equation.
Definition simulation.h:763
vector< double > m_VdTSDeepWaterWaveStationAngle
Time series of wave orientation at deep water wave station.
static vector< string > * VstrSplit(string const *, char const, vector< string > *)
From http://stackoverflow.com/questions/236129/split-a-string-in-c They implement (approximately) Pyt...
Definition utils.cpp:2259
bool m_bOutputErosionPotentialData
Output erosion potential data?
Definition simulation.h:323
int m_nLevel
TODO 007 Used in WAVESETUP + SURGE + RUNUP.
Definition simulation.h:536
bool m_bSedimentInputThisIter
Do we have a sediment input event this iteration?
Definition simulation.h:383
bool m_bCliffNotchSave
Save cliff notch incision depth vector GIS files?
Definition simulation.h:269
bool m_bVectorWaveFloodLineSave
Are we saving the vector wave flood line? TODO 007.
Definition simulation.h:407
vector< string > m_VstrGDALIUSDriverCode
GDAL driver code for the initial unconsolidated sand sediment GIS data.
bool m_bWaveAngleSave
Save wave angle raster GIS files?
Definition simulation.h:104
void AppendPolygon(CGeomCoastPolygon *)
Appends a pointer to a coast polygon to the coast polygon vector.
Definition utils.cpp:2764
static double dCalcWaveAngleToCoastNormal(double const, double const, int const)
Calculates the angle between the wave direction and a normal to the coastline tangent....
bool m_bOmitSearchEastEdge
Omit the east edge of the grid from coast-end searches?
Definition simulation.h:335
vector< int > m_VEdgeCellEdge
The grid edge that each edge cell belongs to.
vector< CGeomCoastPolygon * > m_pVCoastPolygon
Pointers to coast polygon objects, in down-coast sequence TODO 044 Will need to use global polygon ID...
int m_nCoastSmoothWindow
The size of the window used for coast smoothing. Must be an odd number.
Definition simulation.h:443
double m_dThisIterLeftGridUnconsFine
Total fine unconsolidated sediment lost from the grid this iteration (depth in m)
Definition simulation.h:820
bool m_bLocalSlopeSave
Save local slope raster GIS files?
Definition simulation.h:164
Contains CGeomILine definitions.
Contains CGeomLine definitions.
int const SAVEMAX
Definition simulation.h:57
int const NRNG
Definition simulation.h:56