CoastalME (Coastal Modelling Environment)
Simulates the long-term behaviour of complex coastlines
Loading...
Searching...
No Matches
calc_waves.cpp
Go to the documentation of this file.
1
12
13/* ==============================================================================================================================
14
15 This file is part of CoastalME, the Coastal Modelling Environment.
16
17 CoastalME is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.
18
19 This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
20
21 You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
22
23==============================================================================================================================*/
24#include <assert.h>
25
26#include <cstdio>
27
28#include <cmath>
29using std::sin;
30using std::cos;
31using std::pow;
32using std::fmod;
33// using std::isfinite;
34
35// #include <unistd.h>
36
37#include <string>
38// using std::stoi;
39using std::to_string;
40
41#include <iostream>
42// using std::cerr;
43// using std::cout;
44using std::endl;
45using std::ios;
46
47#include <fstream>
48using std::ifstream;
49
50// #include <iomanip>
51// using std::resetiosflags;
52// using std::setiosflags;
53// using std::setprecision;
54// using std::setw;
55
56#include <algorithm>
57// using std::remove;
58using std::reverse_copy;
59// using std::sort;
60
61#include "cme.h"
62#include "coast.h"
63#include "simulation.h"
64#include "cshore/cshore.h"
65#include "2d_point.h"
66
67//===============================================================================================================================
69//===============================================================================================================================
71{
72 // For each coastline, put a value for deep water wave height and direction at each coastline point
73 for (int nCoast = 0; nCoast < static_cast<int>(m_VCoast.size()); nCoast++)
74 {
75 int nDistFromPrevProfile = 0;
76 int nDistToNextProfile = 0;
77
78 double dPrevProfileDeepWaterWaveHeight = 0;
79 double dPrevProfileDeepWaterWaveAngle = 0;
80 double dPrevProfileDeepWaterWavePeriod = 0;
81 double dNextProfileDeepWaterWaveHeight = 0;
82 double dNextProfileDeepWaterWaveAngle = 0;
83 double dNextProfileDeepWaterWavePeriod = 0;
84 double dDist = 0;
85
86 for (int nPoint = 0; nPoint < m_VCoast[nCoast].nGetCoastlineSize(); nPoint++)
87 {
88 // We are going down-coast (i.e. along the coast in the direction of increasing coastline point numbers)
89 if (m_VCoast[nCoast].bIsProfileAtCoastPoint(nPoint))
90 {
91 // OK, a coastline-normal profile begins at this coastline point, so set the deep water wave values at this coastline point to be the values at the seaward end of the coastline normal
92 CGeomProfile const *pProfile = m_VCoast[nCoast].pGetProfileAtCoastPoint(nPoint);
93 // int nProfile = pProfile->nGetCoastID();
94
95 double const dThisDeepWaterWaveHeight = pProfile->dGetProfileDeepWaterWaveHeight();
96 double const dThisDeepWaterWaveAngle = pProfile->dGetProfileDeepWaterWaveAngle();
97 double const dThisDeepWaterWavePeriod = pProfile->dGetProfileDeepWaterWavePeriod();
98
99 m_VCoast[nCoast].SetCoastDeepWaterWaveHeight(nPoint, dThisDeepWaterWaveHeight);
100 m_VCoast[nCoast].SetCoastDeepWaterWaveAngle(nPoint, dThisDeepWaterWaveAngle);
101 m_VCoast[nCoast].SetCoastDeepWaterWavePeriod(nPoint, dThisDeepWaterWavePeriod);
102
103 // Reset for next time
104 nDistFromPrevProfile = 0;
105 dPrevProfileDeepWaterWaveHeight = dThisDeepWaterWaveHeight;
106 dPrevProfileDeepWaterWaveAngle = dThisDeepWaterWaveAngle;
107 dPrevProfileDeepWaterWavePeriod = dThisDeepWaterWavePeriod;
108
109 // Find the next profile
110 CGeomProfile const *pNextProfile = pProfile->pGetDownCoastAdjacentProfile();
111
112 if (pNextProfile == NULL)
113 {
114 // We are at the end of the coast
115 break;
116 }
117
118 // And the distance (in along-coast points) to the next profile
119 nDistToNextProfile = pNextProfile->nGetCoastPoint() - nPoint;
120 dDist = nDistToNextProfile;
121
122 // And the next profile's deep water wave values
123 dNextProfileDeepWaterWaveHeight = pNextProfile->dGetProfileDeepWaterWaveHeight();
124 dNextProfileDeepWaterWaveAngle = pNextProfile->dGetProfileDeepWaterWaveAngle();
125 dNextProfileDeepWaterWavePeriod = pNextProfile->dGetProfileDeepWaterWavePeriod();
126
127 // LogStream << m_ulIter << ": coast point = " << nPoint << " is start of profile " << pProfile->nGetCoastID() << ", next profile is " << pNextProfile->nGetCoastID() << ", which starts at coast piint " << pNextProfile->nGetCoastPoint() << ", dThisDeepWaterWaveHeight = " << dThisDeepWaterWaveHeight << ", dThisDeepWaterWaveAngle = " << dThisDeepWaterWaveAngle << " nDistToNextProfile = " << nDistToNextProfile << endl;
128 }
129
130 else
131 {
132 // This coast point is not the start of a coastline normal, so set the deep water wave values at this coastline point to be a weighted average of those from the up-coast and down-coast profiles
133 nDistFromPrevProfile++;
134 nDistToNextProfile--;
135
136 double const dPrevWeight = (dDist - nDistFromPrevProfile) / dDist;
137 double const dNextWeight = (dDist - nDistToNextProfile) / dDist;
138 double const dThisDeepWaterWaveHeight = (dPrevWeight * dPrevProfileDeepWaterWaveHeight) + (dNextWeight * dNextProfileDeepWaterWaveHeight);
139 double const dThisDeepWaterWaveAngle = dKeepWithin360((dPrevWeight * dPrevProfileDeepWaterWaveAngle) + (dNextWeight * dNextProfileDeepWaterWaveAngle));
140 double const dThisDeepWaterWavePeriod = (dPrevWeight * dPrevProfileDeepWaterWavePeriod) + (dNextWeight * dNextProfileDeepWaterWavePeriod);
141
142 m_VCoast[nCoast].SetCoastDeepWaterWaveHeight(nPoint, dThisDeepWaterWaveHeight);
143 m_VCoast[nCoast].SetCoastDeepWaterWaveAngle(nPoint, dThisDeepWaterWaveAngle);
144 m_VCoast[nCoast].SetCoastDeepWaterWavePeriod(nPoint, dThisDeepWaterWavePeriod);
145
146 // LogStream << m_ulIter << ": coast point = " << nPoint << " dThisDeepWaterWaveHeight = " << dThisDeepWaterWaveHeight << " dThisDeepWaterWaveAngle = " << dThisDeepWaterWaveAngle << " nDistFromPrevProfile = " << nDistFromPrevProfile << " nDistToNextProfile = " << nDistToNextProfile << endl;
147 }
148 }
149 }
150
151 return RTN_OK;
152}
153
154//===============================================================================================================================
156//===============================================================================================================================
158{
159 // Set up all-profile vectors to hold the wave attribute data at every profile point on all profiles
160 vector<bool> VbBreakingAll;
161
162 vector<double> VdXAll;
163 vector<double> VdYAll;
164 vector<double> VdHeightXAll;
165 vector<double> VdHeightYAll;
166
167 // Calculate wave properties for every coast
168 bool bSomeNonStartOrEndOfCoastProfiles = false;
169
170 for (int nCoast = 0; nCoast < static_cast<int>(m_VCoast.size()); nCoast++)
171 {
172 int const nCoastSize = m_VCoast[nCoast].nGetCoastlineSize();
173 int const nNumProfiles = m_VCoast[nCoast].nGetNumProfiles();
174
175 static bool bDownCoast = true;
176
177 // Calculate wave properties at every point along each valid profile, and for the cells under the profiles. Do this alternately in up-coast and down-coast sequence
178 for (int nn = 0; nn < nNumProfiles; nn++)
179 {
180 vector<bool> VbBreaking;
181 vector<double> VdX;
182 vector<double> VdY;
183 vector<double> VdHeightX;
184 vector<double> VdHeightY;
185
186 CGeomProfile *pProfile;
187
188 if (bDownCoast)
189 pProfile = m_VCoast[nCoast].pGetProfileWithDownCoastSeq(nn);
190
191 else
192 pProfile = m_VCoast[nCoast].pGetProfileWithUpCoastSeq(nn);
193
194 int const nRet = nCalcWavePropertiesOnProfile(nCoast, nCoastSize, pProfile, &VdX, &VdY, &VdHeightX, &VdHeightY, &VbBreaking);
195
196 if (nRet != RTN_OK)
197 {
198 if (nRet == RTN_ERR_CSHORE_ERROR)
199 {
200 // Abandon calculations on this profile, and flag the profile
201 pProfile->SetCShoreProblem(true);
202
203 // Move on to next profile
204 continue;
205 }
206
207 else
208 {
209 // A serious CShore error, so abort the run
210 return nRet;
211 }
212 }
213
214 // Are the waves off-shore? If so, do nothing more with this profile. The wave values for cells have already been given the off-shore value
215 if (VbBreaking.empty())
216 continue;
217
218 // Is this a start of coast or end of coast profile?
219 if ((!pProfile->bStartOfCoast()) && (!pProfile->bEndOfCoast()))
220 {
221 // It is neither a start of coast or an end of coast profile, so set switch
222 bSomeNonStartOrEndOfCoastProfiles = true;
223 }
224
225 // // DEBUG CODE ===============================================================================================
226 // for (int nn = 0; nn < VdX.size(); nn++)
227 // {
228 // LogStream << "nProfile = " << nProfile << " nn = " << nn << " VdX[nn] = " << VdX[nn] << " VdY[nn] = " << VdY[nn] << " VdHeightX[nn] = " << VdHeightX[nn] << " VdHeightY[nn] = " << VdHeightY[nn] << " VbBreaking[nn] = " << VbBreaking[nn] << endl;
229 // }
230 // LogStream << endl;
231 // // DEBUG CODE ===============================================================================================
232
233 // Append to the all-profile vectors
234 VdXAll.insert(VdXAll.end(), VdX.begin(), VdX.end());
235 VdYAll.insert(VdYAll.end(), VdY.begin(), VdY.end());
236 VdHeightXAll.insert(VdHeightXAll.end(), VdHeightX.begin(), VdHeightX.end());
237 VdHeightYAll.insert(VdHeightYAll.end(), VdHeightY.begin(), VdHeightY.end());
238 VbBreakingAll.insert(VbBreakingAll.end(), VbBreaking.begin(), VbBreaking.end());
239 }
240
241 bDownCoast = !bDownCoast;
242 }
243
244 // OK, do we have some profiles other than start of coast or end of coast profiles in the all-profile vectors? We need to check this, because GDALGridCreate() in nInterpolateWavePropertiesToWithinPolygonCells() does not work if we give it only a start-of-coast or an end-of-coast profile to work with TODO 006 Is this still true?
245 if (!bSomeNonStartOrEndOfCoastProfiles)
246 {
247 LogStream << m_ulIter << ": waves are on-shore only, for start and/or end of coast profiles" << endl;
248
249 return RTN_OK;
250 }
251
252 // We need to also send the deepwater points from the edge of the grid to nInterpolateWavePropertiesToWithinPolygonCells(), this is necessary to prevent GDALGridCreate() leaving holes in the interpolated grid when the polygons are far from regular
253 double dDeepWaterWaveX;
254 double dDeepWaterWaveY;
255
257 {
258 // Just using the same value of deep water wave hehght and angle for all cells
259 dDeepWaterWaveX = m_dAllCellsDeepWaterWaveHeight * sin(m_dAllCellsDeepWaterWaveAngle * PI / 180),
260 dDeepWaterWaveY = m_dAllCellsDeepWaterWaveHeight * cos(m_dAllCellsDeepWaterWaveAngle * PI / 180);
261 }
262
263 for (int nX = 0; nX < m_nXGridSize; nX++)
264 {
265 if (m_pRasterGrid->m_Cell[nX][0].bIsInContiguousSea())
266 {
267 int const nPolyID = m_pRasterGrid->m_Cell[nX][0].nGetPolygonID();
268
269 if (nPolyID == INT_NODATA)
270 {
271 // Not in a polygon
272 VdXAll.push_back(nX);
273 VdYAll.push_back(0);
274
276 {
277 // Not using the same value of deep water height and angle for all cells, so get this cell's deep water height and angle values
278 dDeepWaterWaveX = m_pRasterGrid->m_Cell[nX][0].dGetCellDeepWaterWaveHeight() * sin(m_pRasterGrid->m_Cell[nX][0].dGetCellDeepWaterWaveAngle() * PI / 180);
279 dDeepWaterWaveY = m_pRasterGrid->m_Cell[nX][0].dGetCellDeepWaterWaveHeight() * cos(m_pRasterGrid->m_Cell[nX][0].dGetCellDeepWaterWaveAngle() * PI / 180);
280 }
281
282 VdHeightXAll.push_back(dDeepWaterWaveX);
283 VdHeightYAll.push_back(dDeepWaterWaveY);
284 }
285 }
286
287 if (m_pRasterGrid->m_Cell[nX][m_nYGridSize - 1].bIsInContiguousSea())
288 {
289 int const nPolyID = m_pRasterGrid->m_Cell[nX][m_nYGridSize - 1].nGetPolygonID();
290
291 if (nPolyID == INT_NODATA)
292 {
293 // Not in a polygon
294 VdXAll.push_back(nX);
295 VdYAll.push_back(m_nYGridSize - 1);
296
298 {
299 // Not using the same value of deep water height and angle for all cells, so get this cell's deep water height and angle values
300 dDeepWaterWaveX = m_pRasterGrid->m_Cell[nX][m_nYGridSize - 1].dGetCellDeepWaterWaveHeight() * sin(m_pRasterGrid->m_Cell[nX][m_nYGridSize - 1].dGetCellDeepWaterWaveAngle() * PI / 180);
301 dDeepWaterWaveY = m_pRasterGrid->m_Cell[nX][m_nYGridSize - 1].dGetCellDeepWaterWaveHeight() * cos(m_pRasterGrid->m_Cell[nX][m_nYGridSize - 1].dGetCellDeepWaterWaveAngle() * PI / 180);
302 }
303
304 VdHeightXAll.push_back(dDeepWaterWaveX);
305 VdHeightYAll.push_back(dDeepWaterWaveY);
306 }
307 }
308 }
309
310 for (int nY = 0; nY < m_nYGridSize; nY++)
311 {
312 if (m_pRasterGrid->m_Cell[0][nY].bIsInContiguousSea())
313 {
314 int const nPolyID = m_pRasterGrid->m_Cell[0][nY].nGetPolygonID();
315
316 if (nPolyID == INT_NODATA)
317 {
318 // Not in a polygon
319 VdXAll.push_back(0);
320 VdYAll.push_back(nY);
321
323 {
324 // Not using the same value of deep water height and angle for all cells, so get this cell's deep water height and angle values
325 dDeepWaterWaveX = m_pRasterGrid->m_Cell[0][nY].dGetCellDeepWaterWaveHeight() * sin(m_pRasterGrid->m_Cell[0][nY].dGetCellDeepWaterWaveAngle() * PI / 180);
326 dDeepWaterWaveY = m_pRasterGrid->m_Cell[0][nY].dGetCellDeepWaterWaveHeight() * cos(m_pRasterGrid->m_Cell[0][nY].dGetCellDeepWaterWaveAngle() * PI / 180);
327 }
328
329 VdHeightXAll.push_back(dDeepWaterWaveX);
330 VdHeightYAll.push_back(dDeepWaterWaveY);
331 }
332 }
333
334 if (m_pRasterGrid->m_Cell[m_nXGridSize - 1][nY].bIsInContiguousSea())
335 {
336 int const nPolyID = m_pRasterGrid->m_Cell[m_nXGridSize - 1][nY].nGetPolygonID();
337
338 if (nPolyID == INT_NODATA)
339 {
340 // Not in a polygon
341 VdXAll.push_back(m_nXGridSize - 1);
342 VdYAll.push_back(nY);
343
345 {
346 // Not using the same value of deep water height and angle for all cells, so get this cell's deep water height and angle values
347 dDeepWaterWaveX = m_pRasterGrid->m_Cell[m_nXGridSize - 1][nY].dGetCellDeepWaterWaveHeight() * sin(m_pRasterGrid->m_Cell[m_nXGridSize - 1][nY].dGetCellDeepWaterWaveAngle() * PI / 180);
348 dDeepWaterWaveY = m_pRasterGrid->m_Cell[m_nXGridSize - 1][nY].dGetCellDeepWaterWaveHeight() * cos(m_pRasterGrid->m_Cell[m_nXGridSize - 1][nY].dGetCellDeepWaterWaveAngle() * PI / 180);
349 }
350
351 VdHeightXAll.push_back(dDeepWaterWaveX);
352 VdHeightYAll.push_back(dDeepWaterWaveY);
353 }
354 }
355 }
356
357 // // DEBUG CODE ============================================================================================================
358 // LogStream << "Out of loop" << endl;
359 // for (int nn = 0; nn < VdXAll.size(); nn++)
360 // {
361 // LogStream << "nn = " << nn << " VdXAll[nn] = " << VdXAll[nn] << " VdYAll[nn] = " << VdYAll[nn] << " VdHeightXAll[nn] = " << VdHeightXAll[nn] << " VdHeightYAll[nn] = " << VdHeightYAll[nn] << " VbBreakingAll[nn] = " << VbBreakingAll[nn] << endl;
362 // }
363 // LogStream << endl;
364 // // DEBUG CODE ============================================================================================================
365
366 // Are the waves off-shore for every profile? If so, do nothing more
367 if (VbBreakingAll.empty())
368 {
369 LogStream << m_ulIter << ": waves off-shore for all profiles" << endl;
370
371 return RTN_OK;
372 }
373
374 // Some waves are on-shore, so interpolate the wave attributes from all profile points to all within-polygon sea cells, also update the active zone status for each cell
375 int nRet = nInterpolateWavesToPolygonCells(&VdXAll, &VdYAll, &VdHeightXAll, &VdHeightYAll);
376
377 if (nRet != RTN_OK)
378 return nRet;
379
380 // // DEBUG CODE ===========================================================================================================
381 // string strOutFile = m_strOutPath;
382 // strOutFile += "sea_wave_height_CHECKPOINT_1_";
383 // strOutFile += to_string(m_ulIter);
384 // strOutFile += ".tif";
385 //
386 // GDALDriver* pDriver = GetGDALDriverManager()->GetDriverByName("gtiff");
387 // GDALDataset* pDataSet = pDriver->Create(strOutFile.c_str(), m_nXGridSize, m_nYGridSize, 1, GDT_Float64, m_papszGDALRasterOptions);
388 // pDataSet->SetProjection(m_strGDALBasementDEMProjection.c_str());
389 // pDataSet->SetGeoTransform(m_dGeoTransform);
390 //
391 // int nn = 0;
392 // double* pdRaster = new double[m_nXGridSize * m_nYGridSize];
393 // for (int nY = 0; nY < m_nYGridSize; nY++)
394 // {
395 // for (int nX = 0; nX < m_nXGridSize; nX++)
396 // {
397 // pdRaster[nn++] = m_pRasterGrid->m_Cell[nX][nY].dGetWaveHeight();
398 // }
399 // }
400 //
401 // GDALRasterBand* pBand = pDataSet->GetRasterBand(1);
402 // pBand->SetNoDataValue(m_dMissingValue);
403 // nRet = pBand->RasterIO(GF_Write, 0, 0, m_nXGridSize, m_nYGridSize, pdRaster, m_nXGridSize, m_nYGridSize, GDT_Float64, 0, 0, NULL);
404 //
405 // if (nRet == CE_Failure)
406 // return RTN_ERR_GRIDCREATE;
407 //
408 // GDALClose(pDataSet);
409 // delete[] pdRaster;
410 // // DEBUG CODE ===========================================================================================================
411
412 // // DEBUG CODE ===========================================================================================================
413 // strOutFile = m_strOutPath;
414 // strOutFile += "sea_wave_angle_CHECKPOINT_1_";
415 // strOutFile += to_string(m_ulIter);
416 // strOutFile += ".tif";
417 //
418 // pDriver = GetGDALDriverManager()->GetDriverByName("gtiff");
419 // pDataSet = pDriver->Create(strOutFile.c_str(), m_nXGridSize, m_nYGridSize, 1, GDT_Float64, m_papszGDALRasterOptions);
420 // pDataSet->SetProjection(m_strGDALBasementDEMProjection.c_str());
421 // pDataSet->SetGeoTransform(m_dGeoTransform);
422 //
423 // nn = 0;
424 // pdRaster = new double[m_nXGridSize * m_nYGridSize];
425 // for (int nY = 0; nY < m_nYGridSize; nY++)
426 // {
427 // for (int nX = 0; nX < m_nXGridSize; nX++)
428 // {
429 // pdRaster[nn++] = m_pRasterGrid->m_Cell[nX][nY].dGetWaveAngle();
430 // }
431 // }
432 //
433 // pBand = pDataSet->GetRasterBand(1);
434 // pBand->SetNoDataValue(m_dMissingValue);
435 // nRet = pBand->RasterIO(GF_Write, 0, 0, m_nXGridSize, m_nYGridSize, pdRaster, m_nXGridSize, m_nYGridSize, GDT_Float64, 0, 0, NULL);
436 //
437 // if (nRet == CE_Failure)
438 // return RTN_ERR_GRIDCREATE;
439 //
440 // GDALClose(pDataSet);
441 // delete[] pdRaster;
442 // // DEBUG CODE ===========================================================================================================
443
444 // Find any shadow zones and then modify waves in and adjacent to them
445 nRet = nDoAllShadowZones();
446
447 if (nRet != RTN_OK)
448 return nRet;
449
450 // Calculate the D50 for each polygon. Also fill in any artefactual 'holes' in active zone and wave property patterns
452
453 // // DEBUG CODE ===========================================================================================================
454 // strOutFile = m_strOutPath;
455 // strOutFile += "sea_wave_height_CHECKPOINT_2_";
456 // strOutFile += to_string(m_ulIter);
457 // strOutFile += ".tif";
458 //
459 // pDriver = GetGDALDriverManager()->GetDriverByName("gtiff");
460 // pDataSet = pDriver->Create(strOutFile.c_str(), m_nXGridSize, m_nYGridSize, 1, GDT_Float64, m_papszGDALRasterOptions);
461 // pDataSet->SetProjection(m_strGDALBasementDEMProjection.c_str());
462 // pDataSet->SetGeoTransform(m_dGeoTransform);
463 //
464 // nn = 0;
465 // pdRaster = new double[m_nXGridSize * m_nYGridSize];
466 // for (int nY = 0; nY < m_nYGridSize; nY++)
467 // {
468 // for (int nX = 0; nX < m_nXGridSize; nX++)
469 // {
470 // pdRaster[nn++] = m_pRasterGrid->m_Cell[nX][nY].dGetWaveHeight();
471 // }
472 // }
473 //
474 // pBand = pDataSet->GetRasterBand(1);
475 // pBand->SetNoDataValue(m_dMissingValue);
476 // nRet = pBand->RasterIO(GF_Write, 0, 0, m_nXGridSize, m_nYGridSize, pdRaster, m_nXGridSize, m_nYGridSize, GDT_Float64, 0, 0, NULL);
477 //
478 // if (nRet == CE_Failure)
479 // return RTN_ERR_GRIDCREATE;
480 //
481 // GDALClose(pDataSet);
482 // delete[] pdRaster;
483 // // DEBUG CODE ===========================================================================================================
484
485 // // DEBUG CODE ===========================================================================================================
486 // strOutFile = m_strOutPath;
487 // strOutFile += "sea_wave_angle_CHECKPOINT_2_";
488 // strOutFile += to_string(m_ulIter);
489 // strOutFile += ".tif";
490 //
491 // pDriver = GetGDALDriverManager()->GetDriverByName("gtiff");
492 // pDataSet = pDriver->Create(strOutFile.c_str(), m_nXGridSize, m_nYGridSize, 1, GDT_Float64, m_papszGDALRasterOptions);
493 // pDataSet->SetProjection(m_strGDALBasementDEMProjection.c_str());
494 // pDataSet->SetGeoTransform(m_dGeoTransform);
495 //
496 // nn = 0;
497 // pdRaster = new double[m_nXGridSize * m_nYGridSize];
498 // for (int nY = 0; nY < m_nYGridSize; nY++)
499 // {
500 // for (int nX = 0; nX < m_nXGridSize; nX++)
501 // {
502 // pdRaster[nn++] = m_pRasterGrid->m_Cell[nX][nY].dGetWaveAngle();
503 // }
504 // }
505 //
506 // pBand = pDataSet->GetRasterBand(1);
507 // pBand->SetNoDataValue(m_dMissingValue);
508 // nRet = pBand->RasterIO(GF_Write, 0, 0, m_nXGridSize, m_nYGridSize, pdRaster, m_nXGridSize, m_nYGridSize, GDT_Float64, 0, 0, NULL);
509 //
510 // if (nRet == CE_Failure)
511 // return RTN_ERR_GRIDCREATE;
512 //
513 // GDALClose(pDataSet);
514 // delete[] pdRaster;
515 // // DEBUG CODE ===========================================================================================================
516
517 // Modify the wave breaking properties (wave height, wave dir, breaking depth, breaking distance) for coastline points within the shadow zone
518 for (int nCoast = 0; nCoast < static_cast<int>(m_VCoast.size()); nCoast++)
519 {
520 int const nNumProfiles = m_VCoast[nCoast].nGetNumProfiles();
521
522 for (int nProfile = 0; nProfile < nNumProfiles; nProfile++)
524 }
525
526 // // DEBUG CODE ===========================================================================================================
527 // strOutFile = m_strOutPath;
528 // strOutFile += "sea_wave_height_CHECKPOINT_3_";
529 // strOutFile += to_string(m_ulIter);
530 // strOutFile += ".tif";
531 //
532 // pDriver = GetGDALDriverManager()->GetDriverByName("gtiff");
533 // pDataSet = pDriver->Create(strOutFile.c_str(), m_nXGridSize, m_nYGridSize, 1, GDT_Float64, m_papszGDALRasterOptions);
534 // pDataSet->SetProjection(m_strGDALBasementDEMProjection.c_str());
535 // pDataSet->SetGeoTransform(m_dGeoTransform);
536 //
537 // nn = 0;
538 // pdRaster = new double[m_nXGridSize * m_nYGridSize];
539 // for (int nY = 0; nY < m_nYGridSize; nY++)
540 // {
541 // for (int nX = 0; nX < m_nXGridSize; nX++)
542 // {
543 // pdRaster[nn++] = m_pRasterGrid->m_Cell[nX][nY].dGetWaveHeight();
544 // }
545 // }
546 //
547 // pBand = pDataSet->GetRasterBand(1);
548 // pBand->SetNoDataValue(m_dMissingValue);
549 // nRet = pBand->RasterIO(GF_Write, 0, 0, m_nXGridSize, m_nYGridSize, pdRaster, m_nXGridSize, m_nYGridSize, GDT_Float64, 0, 0, NULL);
550 //
551 // if (nRet == CE_Failure)
552 // return RTN_ERR_GRIDCREATE;
553 //
554 // GDALClose(pDataSet);
555 // delete[] pdRaster;
556 // // DEBUG CODE ===========================================================================================================
557
558 // // DEBUG CODE ===========================================================================================================
559 // strOutFile = m_strOutPath;
560 // strOutFile += "sea_wave_angle_CHECKPOINT_3_";
561 // strOutFile += to_string(m_ulIter);
562 // strOutFile += ".tif";
563 //
564 // pDriver = GetGDALDriverManager()->GetDriverByName("gtiff");
565 // pDataSet = pDriver->Create(strOutFile.c_str(), m_nXGridSize, m_nYGridSize, 1, GDT_Float64, m_papszGDALRasterOptions);
566 // pDataSet->SetProjection(m_strGDALBasementDEMProjection.c_str());
567 // pDataSet->SetGeoTransform(m_dGeoTransform);
568 //
569 // nn = 0;
570 // pdRaster = new double[m_nXGridSize * m_nYGridSize];
571 // for (int nY = 0; nY < m_nYGridSize; nY++)
572 // {
573 // for (int nX = 0; nX < m_nXGridSize; nX++)
574 // {
575 // pdRaster[nn++] = m_pRasterGrid->m_Cell[nX][nY].dGetWaveAngle();
576 // }
577 // }
578 //
579 // pBand = pDataSet->GetRasterBand(1);
580 // pBand->SetNoDataValue(m_dMissingValue);
581 // nRet = pBand->RasterIO(GF_Write, 0, 0, m_nXGridSize, m_nYGridSize, pdRaster, m_nXGridSize, m_nYGridSize, GDT_Float64, 0, 0, NULL);
582 //
583 // if (nRet == CE_Failure)
584 // return RTN_ERR_GRIDCREATE;
585 //
586 // GDALClose(pDataSet);
587 // delete[] pdRaster;
588 // // DEBUG CODE ===========================================================================================================
589
590 // Interpolate these wave properties for all remaining coastline points. Do this in along-coastline sequence
591 for (int nCoast = 0; nCoast < static_cast<int>(m_VCoast.size()); nCoast++)
592 {
593 int const nCoastSize = m_VCoast[nCoast].nGetCoastlineSize();
594 int const nNumProfiles = m_VCoast[nCoast].nGetNumProfiles();
595
596 // Interpolate these wave properties for all remaining coastline points. Do this in along-coastline sequence, but do not do this for the end-of-coastline profile (which is the final one)
597 for (int n = 0; n < nNumProfiles - 1; n++)
599
600 // And do the same for the coastline cells
602
603 // Calculate wave energy at breaking for every point on the coastline
604 for (int nCoastPoint = 0; nCoastPoint < nCoastSize; nCoastPoint++)
605 {
606 // Equation 4 from Walkden & Hall, 2005
607 double const dBreakingWaveHeight = m_VCoast[nCoast].dGetBreakingWaveHeight(nCoastPoint);
608 double const dCoastPointWavePeriod = m_VCoast[nCoast].dGetCoastDeepWaterWavePeriod(nCoastPoint);
609
610 // TODO 080 Why do we get -ve dBreakingWaveHeight here?
611 if (bFPIsEqual(dBreakingWaveHeight, DBL_NODATA, TOLERANCE))
612 {
613 m_VCoast[nCoast].SetBreakingWaveHeight(nCoastPoint, 0);
614 }
615
616 else
617 {
618 double const dErosiveWaveForce = pow(dBreakingWaveHeight, WALKDEN_HALL_PARAM_1) * pow(dCoastPointWavePeriod, WALKDEN_HALL_PARAM_2);
619
620 // Calculate total wave energy at this coast point during this timestep
621 double const dWaveEnergy = dErosiveWaveForce * m_dTimeStep * 3600;
622 m_VCoast[nCoast].SetWaveEnergyAtBreaking(nCoastPoint, dWaveEnergy);
623 }
624 }
625 }
626
627 return RTN_OK;
628}
629
630//===============================================================================================================================
632//===============================================================================================================================
633double CSimulation::dCalcWaveAngleToCoastNormal(double const dCoastAngle, double const dWaveAngle, int const nSeaHand)
634{
635 double dWaveToNormalAngle = 0;
636
637 if (nSeaHand == LEFT_HANDED)
638 // Left-handed coast
639 dWaveToNormalAngle = fmod((dWaveAngle - dCoastAngle + 360), 360) - 90;
640
641 else
642 // Right-handed coast
643 dWaveToNormalAngle = fmod((dWaveAngle - dCoastAngle + 360), 360) - 270;
644
645 if ((dWaveToNormalAngle >= 90) || (dWaveToNormalAngle <= -90))
646 dWaveToNormalAngle = DBL_NODATA;
647
648 return dWaveToNormalAngle;
649}
650
651//===============================================================================================================================
653//===============================================================================================================================
654int CSimulation::nCalcWavePropertiesOnProfile(int const nCoast, int const nCoastSize, CGeomProfile *pProfile, vector<double> *pVdX, vector<double> *pVdY, vector<double> *pVdHeightX, vector<double> *pVdHeightY, vector<bool> *pVbBreaking)
655{
656 // Only do this for profiles without problems. Still do start- and end-of-coast profiles however
657 if (!pProfile->bOKIncStartAndEndOfCoast())
658 {
659 // if (m_nLogFileDetail >= LOG_FILE_ALL)
660 // LogStream << m_ulIter << ": Coast " << nCoast << ", profile " << nProfile << " has been marked invalid, will not calc wave properties on this profile" << endl;
661
662 return RTN_OK;
663 }
664
665 // Calculate some wave properties based on the wave period following Airy wave theory
666 double const dDeepWaterWavePeriod = pProfile->dGetProfileDeepWaterWavePeriod();
667
668 m_dC_0 = (m_dG * dDeepWaterWavePeriod) / (2 * PI); // Deep water (offshore) wave celerity (m/s)
669 m_dL_0 = m_dC_0 * dDeepWaterWavePeriod; // Deep water (offshore) wave length (m)
670
671 int const nSeaHand = m_VCoast[nCoast].nGetSeaHandedness();
672 int const nCoastPoint = pProfile->nGetCoastPoint();
673
674 // Get the flux orientation (the orientation of a line which is tangential to the coast) at adjacent coastline points. Note special treatment for the coastline end points
675 double const dFluxOrientationThis = m_VCoast[nCoast].dGetFluxOrientation(nCoastPoint);
676 double dFluxOrientationPrev = 0;
677 double dFluxOrientationNext = 0;
678
679 if (nCoastPoint == 0)
680 {
681 dFluxOrientationPrev = dFluxOrientationThis;
682 dFluxOrientationNext = m_VCoast[nCoast].dGetFluxOrientation(1);
683 }
684
685 else if (nCoastPoint == nCoastSize - 1)
686 {
687 dFluxOrientationPrev = m_VCoast[nCoast].dGetFluxOrientation(nCoastPoint - 2);
688 dFluxOrientationNext = dFluxOrientationThis;
689 }
690
691 else
692 {
693 dFluxOrientationPrev = m_VCoast[nCoast].dGetFluxOrientation(nCoastPoint - 1);
694 dFluxOrientationNext = m_VCoast[nCoast].dGetFluxOrientation(nCoastPoint + 1);
695 }
696
697 // Get the deep water wave orientation for this profile
698 double const dDeepWaterWaveAngle = pProfile->dGetProfileDeepWaterWaveAngle();
699
700 // Calculate the angle between the deep water wave direction and a normal to the coast tangent
701 double dWaveToNormalAngle = dCalcWaveAngleToCoastNormal(dFluxOrientationThis, dDeepWaterWaveAngle, nSeaHand);
702
703 // Are the waves off-shore?
704 if (bFPIsEqual(dWaveToNormalAngle, DBL_NODATA, TOLERANCE))
705 {
706 // They are so, do nothing (each cell under the profile has already been initialised with deep water wave height and wave direction)
707 // LogStream << m_ulIter << ": profile " << nProfile << " has sea to " << (m_VCoast[nCoast].nGetSeaHandedness() == RIGHT_HANDED ? "right" : "left") << " dWaveToNormalAngle = " << dWaveToNormalAngle << " which is off-shore" << endl;
708
709 return RTN_OK;
710 }
711
712 // LogStream << m_ulIter << ": profile = " << nProfile << " has sea to " << (m_VCoast[nCoast].nGetSeaHandedness() == RIGHT_HANDED ? "right" : "left") << " dWaveToNormalAngle = " << dWaveToNormalAngle << " which is " << (dWaveToNormalAngle < 0 ? "DOWN" : "UP") << "-coast" << endl;
713
714 // Calculate the angle between the deep water wave direction and a normal to the coast tangent for the previous coast point
715 double dWaveToNormalAnglePrev;
716
717 if (nCoastPoint > 0)
718 {
719 // Get the deep water wave orientation for the up-coast point
720 double const dPrevDeepWaterWaveAngle = m_VCoast[nCoast].dGetCoastDeepWaterWaveAngle(nCoastPoint - 1);
721
722 dWaveToNormalAnglePrev = dCalcWaveAngleToCoastNormal(dFluxOrientationPrev, dPrevDeepWaterWaveAngle, nSeaHand);
723 }
724
725 else
726 {
727 dWaveToNormalAnglePrev = dWaveToNormalAngle;
728 }
729
730 // if (dWaveToNormalAnglePrev == DBL_NODATA)
731 // LogStream << "\tPrevious profile, dWaveToNormalAnglePrev = " << dWaveToNormalAnglePrev << " which is off-shore" << endl;
732 // else
733 // LogStream << "\tPrevious profile, dWaveToNormalAnglePrev = " << dWaveToNormalAnglePrev << " which is " << (dWaveToNormalAnglePrev < 0 ? "DOWN" : "UP") << "-coast" << endl;
734
735 // Calculate the angle between the deep water wave direction and a normal to the coast tangent for the next coast point
736 double dWaveToNormalAngleNext;
737
738 if (nCoastPoint < nCoastSize - 1)
739 {
740 // Get the deep water wave orientation for the down-coast point
741 double const dNextDeepWaterWaveAngle = m_VCoast[nCoast].dGetCoastDeepWaterWaveAngle(nCoastPoint + 1);
742
743 dWaveToNormalAngleNext = dCalcWaveAngleToCoastNormal(dFluxOrientationNext, dNextDeepWaterWaveAngle, nSeaHand);
744 }
745
746 else
747 {
748 dWaveToNormalAngleNext = dWaveToNormalAngle;
749 }
750
751 // if (dWaveToNormalAngleNext == DBL_NODATA)
752 // LogStream << "\tNext profile, dWaveToNormalAngleNext = " << dWaveToNormalAngleNext << " which is off-shore" << endl;
753 // else
754 // LogStream << "\tNext profile, dWaveToNormalAngleNext = " << dWaveToNormalAngleNext << " which is " << (dWaveToNormalAngleNext < 0 ? "DOWN" : "UP") << "-coast" << endl;
755
756 // Following Ashton and Murray (2006), if we have high-angle waves then use the flux orientation of the previous (up-coast) profile, if transitioning from diffusive to antidiffusive use flux maximizing angle (45 degrees)
757 if ((dWaveToNormalAngle > 0) && (!bFPIsEqual(dWaveToNormalAnglePrev, DBL_NODATA, TOLERANCE)) && (dWaveToNormalAnglePrev > 0))
758 {
759 if (dWaveToNormalAngle > 45)
760 {
761 if (dWaveToNormalAnglePrev < 45)
762 {
763 dWaveToNormalAngle = 45;
764 // LogStream << "\tA1" << endl;
765 }
766
767 else
768 {
769 dWaveToNormalAngle = dWaveToNormalAnglePrev;
770 // LogStream << "\tA2" << endl;
771 }
772 }
773 }
774
775 else if ((dWaveToNormalAngle < 0) && (!bFPIsEqual(dWaveToNormalAngleNext, DBL_NODATA, TOLERANCE)) && (dWaveToNormalAngleNext < 0))
776 {
777 if (dWaveToNormalAngle < -45)
778 {
779 if (dWaveToNormalAngleNext > -45)
780 {
781 dWaveToNormalAngle = -45;
782 // LogStream << "\tB1" << endl;
783 }
784
785 else
786 {
787 dWaveToNormalAngle = dWaveToNormalAngleNext;
788 // LogStream << "\tB2" << endl;
789 }
790 }
791 }
792
793 else if ((dWaveToNormalAngle > 45) && (!bFPIsEqual(dWaveToNormalAnglePrev, DBL_NODATA, TOLERANCE)) && (dWaveToNormalAnglePrev > 0))
794 {
795 // The wave direction here has an up-coast (decreasing indices) component: so for high-angle waves use the orientation from the up-coast (previous) profile
796 // LogStream << "\tCCC" << endl;
797
798 dWaveToNormalAngle = dFluxOrientationPrev;
799 }
800
801 else if ((dWaveToNormalAngle < -45) && (!bFPIsEqual(dWaveToNormalAngleNext, DBL_NODATA, TOLERANCE)) && (dWaveToNormalAngleNext < 0))
802 {
803 // The wave direction here has a down-coast (increasing indices) component: so for high-angle waves use the orientation from the down-coast (next) profile
804 // LogStream << "\tDDD" << endl;
805
806 dWaveToNormalAngle = dFluxOrientationNext;
807 }
808
809 // Initialize the wave properties at breaking for this profile
810 bool bBreaking = false;
811
812 int const nProfileSize = pProfile->nGetNumCellsInProfile();
813 int nProfileBreakingDist = 0;
814
815 double dProfileBreakingWaveHeight = DBL_NODATA;
816 double dProfileBreakingWaveAngle = 0;
817 double dProfileBreakingDepth = 0;
818 double dProfileWaveHeight = DBL_NODATA;
819 double const dProfileDeepWaterWaveHeight = pProfile->dGetProfileDeepWaterWaveHeight();
820 double dProfileWaveAngle = DBL_NODATA;
821 double const dProfileDeepWaterWaveAngle = pProfile->dGetProfileDeepWaterWaveAngle();
822
823 vector<bool> VbWaveIsBreaking(nProfileSize, 0);
824
825 vector<double> VdWaveHeight(nProfileSize, 0);
826 vector<double> VdWaveSetupSurge(nProfileSize, 0);
827 // vector<double> VdStormSurge(nProfileSize, 0);
828 vector<double> const VdWaveSetupRunUp(nProfileSize, 0);
829 vector<double> VdWaveDirection(nProfileSize, 0);
830
832 {
833 // We are using CShore to propagate the waves
834 double const dCShoreTimeStep = 3600; // In seconds, not important because we are not using CShore to erode the profile, just to get the hydrodynamics
835 double const dSurgeLevel = CSHORE_SURGE_LEVEL;
836
837 // Set up vectors for the coastline-normal profile elevations. The length of this vector line is given by the number of cells 'under' the profile. Thus each point on the vector relates to a single cell in the grid. This assumes that all points on the profile vector are equally spaced (not quite true, depends on the orientation of the line segments which comprise the profile)
838 vector<double> VdProfileZ; // Initial (pre-erosion) elevation of both consolidated and unconsolidated sediment for cells 'under' the profile, in CShore units
839 vector<double> VdProfileDistXY; // Along-profile distance measured from the seaward limit, in CShore units
840 vector<double> VdProfileFrictionFactor; // Along-profile friction factor from seaward limit
841
842 // The elevation of each of these profile points is the elevation of the centroid of the cell that is 'under' the point. However we cannot always be confident that this is the 'true' elevation of the point on the vector since (unless the profile runs planview N-S or W-E) the vector does not always run exactly through the centroid of the cell
843 int nRet = nGetThisProfileElevationsForCShore(nCoast, pProfile, nProfileSize, &VdProfileDistXY, &VdProfileZ, &VdProfileFrictionFactor);
844
845 if (nRet != RTN_OK)
846 {
847 // Could not create the profile elevation vectors
848 LogStream << m_ulIter << ": could not create CShore profile elevation vectors for profile " << pProfile->nGetCoastID() << endl;
849
850 return nRet;
851 }
852
853 // assert(static_cast<int>(VdProfileDistXY.size()) == nProfileSize);
854 // assert(static_cast<int>(VdProfileZ.size()) == nProfileSize);
855 // assert(static_cast<int>(VdProfileFrictionFactor.size()) == nProfileSize);
856
857 if (VdProfileDistXY.empty())
858 {
859 // The profile elevation vector was created, but was not populated
860 LogStream << m_ulIter << ": could not populate CShore profile elevation vector for profile " << pProfile->nGetCoastID() << endl;
861
863 }
864
865 // Constrain the wave to normal angle to be between -80 and 80 degrees, this is a requirement of CShore
866 dWaveToNormalAngle = tMax(dWaveToNormalAngle, -80.0);
867 dWaveToNormalAngle = tMin(dWaveToNormalAngle, 80.0);
868
869 int const nProfileDistXYSize = static_cast<int>(VdProfileDistXY.size());
870 vector<double> VdFreeSurfaceStd(nProfileDistXYSize, 0); // This is converted to Hrms by Hrms = sqr(8)*FreeSurfaceStd
871 vector<double> VdSinWaveAngleRadians(nProfileDistXYSize, 0); // This is converted to deg by asin(VdSinWaveAngleRadians)*(180/pi)
872 vector<double> VdFractionBreakingWaves(nProfileDistXYSize, 0); // Is 0 if no wave breaking, and 1 if all waves breaking
873
874 // Now define the other values that CShore requires
875 int const nILine = 1; // This is the number of cross-shore lines i.e. the number of CoastalME profiles. Only one at a time, at present
876 int const nIProfl = 0; // 0 for fixed bottom profile, 1 for profile evolution computation
877 int const nIPerm = 0; // 0 for impermeable bottom, 1 for permeable bottom of stone structure
878 int const nIOver = 0; // 0 for no wave overtopping and overflow on crest, 1 for wave overtopping and overflow
879 int const nIWCInt = 0; // 0 for no wave/current interaction, 1 for wave/current interaction in frequency dispersion, momentum and wave action equations
880 int const nIRoll = 0; // 0 for no roller effects, 1 for roller effects in governing equations
881 int const nIWind = 0; // 0 for no wind effects, 1 for wind shear stresses on momentum equations
882 int const nITide = 0; // 0 for no tidal effect on currents, 1 for longshore and cross-shore tidal currents
883 int const nILab = 0; // 0 for field data set, 1 for laboratory data set
884 int const nNWave = 1; // Number of waves at x = 0 starting from time = 0
885 int const nNSurge = 1; // Number of water levels at x = 0 from time = 0
886 double const dDX = VdProfileDistXY.back() / (static_cast<double>(VdProfileDistXY.size() - 1));
887 double const dWaveInitTime = 0; // CShore wave start time
888 double const dSurgeInitTime = 0; // CShore surge start time
889
890#if defined CSHORE_FILE_INOUT || CSHORE_BOTH
891 // Move to the CShore folder
892 nRet = chdir(CSHORE_DIR.c_str());
893
894 if (nRet != RTN_OK)
895 return nRet;
896
897#endif
898
899#if defined CSHORE_FILE_INOUT
900 // We are communicating with CShore using ASCII files, so create an input file for this profile which will be read by CShore
901 nRet = nCreateCShoreInfile(nCoast, nProfile, nILine, nIProfl, nIPerm, nIOver, nIWCInt, nIRoll, nIWind, nITide, nILab, nNWave, nNSurge, dDX, dCShoreTimeStep, dWaveInitTime, dDeepWaterWavePeriod, dProfileDeepWaterWaveHeight, dWaveToNormalAngle, dSurgeInitTime, dSurgeLevel, &VdProfileDistXY, &VdProfileZ, &VdProfileFrictionFactor);
902
903 if (nRet != RTN_OK)
904 return nRet;
905
906 // Set the error flag: this will be changed to 0 within CShore if CShore returns correctly
907 nRet = -1;
908
909 // Run CShore for this profile
910 CShore(&nRet);
911 // CShore(&nRet, &m_ulIter, &nProfile, &nProfile);
912
913 // Check return code for error
914 if (nRet != 0)
915 {
916 string strErr;
917
918 switch (nRet)
919 {
920 case -1:
921 strErr = to_string(m_ulIter) + ": CShore WARNING 1: negative depth at the first node ";
922 break;
923
924 case 2:
925 strErr = to_string(m_ulIter) + ": CShore WARNING 2: negative value at end of landward marching computation ";
926 break;
927
928 case 3:
929 strErr = to_string(m_ulIter) + ": CShore WARNING 3: large energy gradients at the first node: small waves with short period at sea boundary ";
930 break;
931
932 case 4:
933 strErr = to_string(m_ulIter) + ": CShore WARNING 4: zero energy at the first node ";
934 break;
935
936 case 5:
937 strErr = to_string(m_ulIter) + ": CShore WARNING 5: at end of landward marching computation, insufficient water depth ";
938 break;
939
940 case 7:
941 strErr = to_string(m_ulIter) + ": CShore WARNING 7: did not reach convergence ";
942 break;
943 }
944
945 strErr += "(coast " + to_string(nCoast) + " profile " + to_string(pProfile->nGetCoastID()) + " profile length " + to_string(nOutSize) + ")\n";
946
947 // OK, give up for this profile
948 // LogStream << strErr;
949 //
950 // return RTN_ERR_CSHORE_ERROR;
951 }
952
953 // Fetch the CShore results by reading files written by CShore
954 string strOSETUP = "OSETUP";
955 string strOYVELO = "OYVELO";
956 string strOPARAM = "OPARAM";
957
958 nRet = nReadCShoreOutput(nProfile, &strOSETUP, 4, 4, &VdProfileDistXY, &VdFreeSurfaceStd);
959
960 if (nRet != RTN_OK)
961 return nRet;
962
963 nRet = nReadCShoreOutput(nProfile, &strOYVELO, 4, 2, &VdProfileDistXY, &VdSinWaveAngleRadians);
964
965 if (nRet != RTN_OK)
966 return nRet;
967
968 nRet = nReadCShoreOutput(nProfile, &strOPARAM, 4, 3, &VdProfileDistXY, &VdFractionBreakingWaves);
969
970 if (nRet != RTN_OK)
971 return nRet;
972
973 // Read surge outputs
974 // VdTSurg = {dSurgeInitTime, dCShoreTimeStep}, // Ditto
975 // VdSWLin = {dSurgeLevel, dSurgeLevel}, // Ditto
976
977 // Clean up the CShore outputs
978#ifdef _WIN32
979 nRet = system("./clean.bat");
980#else
981
983 {
984 string strCommand = "./save_CShore_output.sh ";
985 strCommand += to_string(m_ulIter);
986 strCommand += " ";
987 strCommand += to_string(nCoast);
988 strCommand += " ";
989 strCommand += to_string(nProfile);
990
991 nRet = system(strCommand.c_str());
992
993 if (nRet != RTN_OK)
994 return nRet;
995 }
996
997 nRet = system("./clean.sh");
998#endif
999
1000 if (nRet != RTN_OK)
1001 return nRet;
1002
1003 // And return to the CoastalME folder
1004 nRet = chdir(m_strCMEDir.c_str());
1005
1006 if (nRet != RTN_OK)
1007 return nRet;
1008
1009#endif
1010
1011#if defined CSHORE_ARG_INOUT || CSHORE_BOTH
1012 // We are communicating with CShore by passing arguments
1013 int nOutSize = 0; // CShore will return the size of the output vectors
1014
1015 // Set the error flag: this will be changed within CShore if there is a problem
1016 nRet = 0;
1017
1018 vector<double> VdInitTime = {dWaveInitTime, dCShoreTimeStep}; // Size is nNwave+1, value 1 is for the start of the CShore run, value 2 for end of CShore run
1019 vector<double> VdTPIn = {dDeepWaterWavePeriod, dDeepWaterWavePeriod}; // Ditto
1020 vector<double> VdHrmsIn = {dProfileDeepWaterWaveHeight, dProfileDeepWaterWaveHeight}; // Ditto
1021 vector<double> VdWangIn = {dWaveToNormalAngle, dWaveToNormalAngle}; // Ditto
1022 vector<double> VdTSurg = {dSurgeInitTime, dCShoreTimeStep}; // Ditto
1023 vector<double> VdSWLin = {dSurgeLevel, dSurgeLevel}; // Ditto
1024 vector<double> VdFPInp = VdProfileFrictionFactor; // Set the value for wave friction at every point of the normal profile
1025 vector<double> VdXYDistFromCShoreOut(CSHOREARRAYOUTSIZE, 0); // Output from CShore
1026 vector<double> VdFreeSurfaceStdOut(CSHOREARRAYOUTSIZE, 0); // Ditto
1027 vector<double> VdWaveSetupSurgeOut(CSHOREARRAYOUTSIZE, 0); // Ditto
1028 // vector<double> VdStormSurgeOut(CSHOREARRAYOUTSIZE, 0); // Ditto
1029 vector<double> const VdWaveSetupRunUpOut(CSHOREARRAYOUTSIZE, 0); // Ditto
1030 vector<double> VdSinWaveAngleRadiansOut(CSHOREARRAYOUTSIZE, 0); // Ditto
1031 vector<double> VdFractionBreakingWavesOut(CSHOREARRAYOUTSIZE, 0); // Ditto
1032
1033 // Call CShore using the argument-passing wrapper
1034 // long lIter = static_cast<long>(m_ulIter); // Bodge to get round compiler 'invalid conversion' error
1035
1036 CShoreWrapper(&nILine, /* In_ILINE */
1037 &nIProfl, /* In_IPROFL */
1038 &nIPerm, /* In_IPERM */
1039 &nIOver, /* In_IOVER */
1040 &nIWCInt, /* In_IWCINT */
1041 &nIRoll, /* In_IROLL */
1042 &nIWind, /* In_IWIND */
1043 &nITide, /* In_ITIDE */
1044 &nILab, /* In_ILAB */
1045 &nNWave, /* In_NWAVE */
1046 &nNSurge, /* In_NSURG */
1047 &dDX, /* In_DX */
1048 &m_dBreakingWaveHeightDepthRatio, /* In_GAMMA */
1049 &VdInitTime[0], /* In_TWAVE */
1050 &VdTPIn[0], /* In_TPIN */
1051 &VdHrmsIn[0], /* In_HRMSIN */
1052 &VdWangIn[0], /* In_WANGIN */
1053 &VdTSurg[0], /* In_TSURG */
1054 &VdSWLin[0], /* In_SWLIN */
1055 &nProfileDistXYSize, /* In_NBINP */
1056 &VdProfileDistXY[0], /* In_XBINP */
1057 &VdProfileZ[0], /* In_ZBINP */
1058 &VdFPInp[0], /* In_FBINP */
1059 &nRet, /* Out_IError */
1060 &nOutSize, /* Out_nOutSize */
1061 &VdXYDistFromCShoreOut[0], /* Out_XYDist */
1062 &VdFreeSurfaceStdOut[0], /* Out_FreeSurfaceStd */
1063 &VdWaveSetupSurgeOut[0], /* Out_WaveSetupSurge */
1064 &VdSinWaveAngleRadiansOut[0], /* Out_SinWaveAngleRadians */
1065 &VdFractionBreakingWavesOut[0]); /* Out_FractionBreakingWaves */
1066
1067 // OK, now check for warnings and errors
1068 if (nOutSize < 2)
1069 {
1070 // CShore sometimes returns only one row of results, which contains data only for the seaward point of the profile. This happens when all other (more coastward) points give an invalid result during CShore's calculations. This is a problem. We don't want to abandon the simulation just because of this, so instead we just put some dummy data into the second row, and carry on with these two rows. The profile will get ignored later, since it is too small to be useful
1072 LogStream << m_ulIter << ": " << WARN << "for coast " << nCoast << " profile " << pProfile->nGetCoastID() << ", only " << nOutSize << " CShore output rows, abandoning this profile" << endl;
1073
1074 // // Set dummy data in the second row
1075 // VdXYDistFromCShoreOut[1] = 1e-5; // Dummy data, must not be the same as VdXYDistFromCShoreOut[0] tho', or get crash in linear interpolation routine
1076 // VdFreeSurfaceStdOut[1] = VdFreeSurfaceStdOut[0];
1077 // VdSinWaveAngleRadiansOut[1] = VdSinWaveAngleRadiansOut[0];
1078 // VdFractionBreakingWavesOut[1] = VdFractionBreakingWavesOut[0];
1079 // VdWaveSetupSurgeOut[1] = VdWaveSetupSurgeOut[0];
1080 // // VdStormSurgeOut[1] = VdStormSurgeOut[0];
1081 // // VdWaveSetupRunUpOut[1] = VdWaveSetupRunUpOut[0];
1082 //
1083 // // And increase the expected number of rows
1084 // nOutSize = 2;
1085
1086 return RTN_ERR_CSHORE_ERROR;
1087 }
1088
1089 if (nRet != RTN_OK)
1090 {
1091 string strErr;
1092
1093 switch (nRet)
1094 {
1095 case -1:
1096 strErr = to_string(m_ulIter) + ": CShore WARNING 1: negative depth at the first node ";
1097 break;
1098
1099 case 2:
1100 strErr = to_string(m_ulIter) + ": CShore WARNING 2: negative value at end of landward marching computation ";
1101 break;
1102
1103 case 3:
1104 strErr = to_string(m_ulIter) + ": CShore WARNING 3: large energy gradients at the first node: small waves with short period at sea boundary ";
1105 break;
1106
1107 case 4:
1108 strErr = to_string(m_ulIter) + ": CShore WARNING 4: zero energy at the first node ";
1109 break;
1110
1111 case 5:
1112 strErr = to_string(m_ulIter) + ": CShore WARNING 5: at end of landward marching computation, insufficient water depth ";
1113 break;
1114
1115 case 7:
1116 strErr = to_string(m_ulIter) + ": CShore WARNING 7: did not reach convergence ";
1117 break;
1118 }
1119
1120 strErr += "(coast " + to_string(nCoast) + " profile " + to_string(pProfile->nGetCoastID()) + " profile length " + to_string(nOutSize) + ")\n";
1121 LogStream << strErr;
1122
1123 // OK, give up for this profile
1124 // return RTN_ERR_CSHORE_ERROR;
1125 }
1126
1127 // LogStream << m_ulIter << ": interpolating profile " << nProfile << endl;
1128
1129 // All OK, so interpolate the CShore output, and convert from the CShore convention (cross-shore distance has its origin at the seaward end) to the CoastalME convention (origin at the shoreline)
1130 InterpolateCShoreOutput(&VdProfileDistXY, nOutSize, nProfileSize, &VdXYDistFromCShoreOut, &VdFreeSurfaceStdOut, &VdWaveSetupSurgeOut, &VdSinWaveAngleRadiansOut, &VdFractionBreakingWavesOut, &VdFreeSurfaceStd, &VdWaveSetupSurge, &VdSinWaveAngleRadians, &VdFractionBreakingWaves);
1131#endif
1132
1133#if defined CSHORE_BOTH
1134#if !defined _WIN32
1135
1137 {
1138 string strCommand = "./save_CShore_output.sh ";
1139 strCommand += to_string(m_ulIter);
1140 strCommand += " ";
1141 strCommand += to_string(nCoast);
1142 strCommand += " ";
1143 strCommand += to_string(nProfile);
1144
1145 nRet = system(strCommand.c_str());
1146
1147 if (nRet != RTN_OK)
1148 return nRet;
1149 }
1150
1151#endif
1152
1153 // Return to the CoastalME folder
1154 nRet = chdir(m_strCMEDir.c_str());
1155
1156 if (nRet != RTN_OK)
1157 return nRet;
1158
1159#endif
1160
1161 // Do more safety checks, then convert the CShore output to wave height and wave direction, and update wave profile attributes
1162 for (int nProfilePoint = (nProfileSize - 1); nProfilePoint >= 0; nProfilePoint--)
1163 {
1164 int const nX = pProfile->pPtiGetCellInProfile(nProfilePoint)->nGetX();
1165 int const nY = pProfile->pPtiGetCellInProfile(nProfilePoint)->nGetY();
1166
1167 // Safety check
1168 if (nProfilePoint > static_cast<int>(VdFreeSurfaceStd.size()) - 1)
1169 continue;
1170
1171 // Safety check: deal with NaN values
1172 if (isnan(VdFreeSurfaceStd[nProfilePoint]))
1173 VdFreeSurfaceStd[nProfilePoint] = 0;
1174
1175 VdWaveHeight[nProfilePoint] = sqrt(8) * VdFreeSurfaceStd[nProfilePoint];
1176
1177 // Another safety check: deal with NaN values
1178 if (isnan(VdSinWaveAngleRadians[nProfilePoint]))
1179 {
1180 VdSinWaveAngleRadians[nProfilePoint] = 0;
1181 VdWaveHeight[nProfilePoint] = 0;
1182 }
1183
1184 // More safety checks: constrain to the interval -1 to +1 to keep asin() happy
1185 if (VdSinWaveAngleRadians[nProfilePoint] < -1)
1186 VdSinWaveAngleRadians[nProfilePoint] = -1;
1187
1188 if (VdSinWaveAngleRadians[nProfilePoint] > 1)
1189 VdSinWaveAngleRadians[nProfilePoint] = 1;
1190
1191 double const dAlpha = asin(VdSinWaveAngleRadians[nProfilePoint]) * (180 / PI);
1192
1193 if (nSeaHand == LEFT_HANDED)
1194 VdWaveDirection[nProfilePoint] = dKeepWithin360(dAlpha + 90 + dFluxOrientationThis);
1195
1196 else
1197 VdWaveDirection[nProfilePoint] = dKeepWithin360(dAlpha + 270 + dFluxOrientationThis);
1198
1199 // Yet another safety check: deal with NaN values
1200 if (isnan(VdFractionBreakingWaves[nProfilePoint]))
1201 {
1202 VdFractionBreakingWaves[nProfilePoint] = 0;
1203 VdWaveHeight[nProfilePoint] = 0;
1204 }
1205
1206 // if ((VdFractionBreakingWaves[nProfilePoint] >= 0.10) && (! bBreaking)) // Sometimes is possible that waves break again
1207 if ((VdFractionBreakingWaves[nProfilePoint] >= 0.10) && (m_dDepthOfClosure >= m_pRasterGrid->m_Cell[nX][nY].dGetSeaDepth()) && (!bBreaking))
1208 {
1209 bBreaking = true;
1210 // assert(VdWaveHeight[nProfilePoint] >= 0);
1211 dProfileBreakingWaveHeight = VdWaveHeight[nProfilePoint];
1212 dProfileBreakingWaveAngle = VdWaveDirection[nProfilePoint];
1213 dProfileBreakingDepth = m_pRasterGrid->m_Cell[nX][nY].dGetSeaDepth(); // Water depth for the cell 'under' this point in the profile
1214 nProfileBreakingDist = nProfilePoint + 1; // At the nearest point nProfilePoint = 0, so, plus one
1215
1216 // LogStream << m_ulIter << ": CShore breaking at [" << nX << "][" << nY << "] = {" << dGridCentroidXToExtCRSX(nX) << ", " << dGridCentroidYToExtCRSY(nY) << "} nProfile = " << nProfile << ", nProfilePoint = " << nProfilePoint << ", dBreakingWaveHeight = " << dBreakingWaveHeight << ", dBreakingWaveAngle = " << dBreakingWaveAngle << ", dProfileBreakingDepth = " << dProfileBreakingDepth << ", nProfileBreakingDist = " << nProfileBreakingDist << endl;
1217 }
1218
1219 VbWaveIsBreaking[nProfilePoint] = bBreaking;
1220 }
1221
1222 if (dProfileBreakingWaveHeight >= dProfileDeepWaterWaveHeight)
1223 {
1224 dProfileBreakingWaveHeight = DBL_NODATA; // checking poorly conditions profiles problems for cshore
1225 }
1226 }
1227
1229 {
1230 // We are using COVE's linear wave theory to propagate the waves
1231 double const dDepthLookupMax = m_dWaveDepthRatioForWaveCalcs * dProfileDeepWaterWaveHeight;
1232
1233 // Go landwards along the profile, calculating wave height and wave angle for every inundated point on the profile (don't do point zero, this is on the coastline) until the waves start to break after breaking wave height is assumed to decrease linearly to zero at the shoreline and wave angle is equalt to wave angle at breaking
1234 for (int nProfilePoint = (nProfileSize - 1); nProfilePoint >= 0; nProfilePoint--)
1235 {
1236 int const nX = pProfile->pPtiGetCellInProfile(nProfilePoint)->nGetX();
1237 int const nY = pProfile->pPtiGetCellInProfile(nProfilePoint)->nGetY();
1238
1239 // Safety check
1240 if (!m_pRasterGrid->m_Cell[nX][nY].bIsInContiguousSea())
1241 continue;
1242
1243 double const dSeaDepth = m_pRasterGrid->m_Cell[nX][nY].dGetSeaDepth(); // Water depth for the cell 'under' this point in the profile
1244
1245 if (dSeaDepth > dDepthLookupMax)
1246 {
1247 // Sea depth is too large relative to wave height to feel the bottom, so do nothing since each cell under the profile has already been initialised with deep water wave height and wave direction
1248 dProfileWaveHeight = dProfileDeepWaterWaveHeight;
1249 dProfileWaveAngle = dProfileDeepWaterWaveAngle;
1250 }
1251
1252 else
1253 {
1254 if (!bBreaking)
1255 {
1256 // Start calculating wave properties using linear wave theory
1257 double const dL = m_dL_0 * sqrt(tanh((2 * PI * dSeaDepth) / m_dL_0)); // Wavelength (m) in intermediate-shallow waters
1258 double const dC = m_dC_0 * tanh((2 * PI * dSeaDepth) / dL); // Wave speed (m/s) set by dSeaDepth, dL and m_dC_0
1259 double const dk = 2 * PI / dL; // Wave number (1/m)
1260 double const dn = ((2 * dSeaDepth * dk) / (sinh(2 * dSeaDepth * dk)) + 1) / 2; // Shoaling factor
1261 double const dKs = sqrt(m_dC_0 / (dn * dC * 2)); // Shoaling coefficient
1262 double const dAlpha = (180 / PI) * asin((dC / m_dC_0) * sin((PI / 180) * dWaveToNormalAngle)); // Calculate angle between wave direction and the normal to the coast tangent
1263 double const dKr = sqrt(cos((PI / 180) * dWaveToNormalAngle) / cos((PI / 180) * dAlpha)); // Refraction coefficient
1264 dProfileWaveHeight = dProfileDeepWaterWaveHeight * dKs * dKr; // Calculate wave height, based on the previous (more seaward) wave height
1265
1266 if (nSeaHand == LEFT_HANDED)
1267 dProfileWaveAngle = dKeepWithin360(dAlpha + 90 + dFluxOrientationThis);
1268
1269 else
1270 dProfileWaveAngle = dKeepWithin360(dAlpha + 270 + dFluxOrientationThis);
1271
1272 // Test to see if the wave breaks at this depth
1273 if (dProfileWaveHeight < (dSeaDepth * m_dBreakingWaveHeightDepthRatio))
1274 {
1275 dProfileBreakingWaveHeight = dProfileWaveHeight;
1276 dProfileBreakingWaveAngle = dProfileWaveAngle;
1277 }
1278 }
1279
1280 else
1281 {
1282 // It does
1283 bBreaking = true;
1284
1285 // Wave has already broken
1286 dProfileWaveAngle = dProfileBreakingWaveAngle; // Wave orientation remains equal to wave orientation at breaking
1287
1288 dProfileWaveHeight = dProfileBreakingWaveHeight * (nProfilePoint / nProfileBreakingDist); // Wave height decreases linearly to zero at shoreline
1289 // dProfileWaveHeight = dSeaDepth * m_dBreakingWaveHeightDepthRatio; // Wave height is limited by depth
1290 dProfileBreakingDepth = dSeaDepth;
1291 nProfileBreakingDist = nProfilePoint;
1292 }
1293 }
1294
1295 // Save current wave attributes
1296 VdWaveDirection[nProfilePoint] = dProfileWaveAngle;
1297 VdWaveHeight[nProfilePoint] = dProfileWaveHeight;
1298 VbWaveIsBreaking[nProfilePoint] = bBreaking;
1299 }
1300 }
1301
1302 // Go landwards along the profile, fetching the calculated wave height and wave angle for every inundated point on this profile
1303 for (int nProfilePoint = (nProfileSize - 1); nProfilePoint >= 0; nProfilePoint--)
1304 {
1305 int nX = pProfile->pPtiGetCellInProfile(nProfilePoint)->nGetX();
1306 int nY = pProfile->pPtiGetCellInProfile(nProfilePoint)->nGetY();
1307
1308 // Safety check
1309 if (!m_pRasterGrid->m_Cell[nX][nY].bIsInContiguousSea())
1310 continue;
1311
1312 // Get the wave attributes calculated for this profile: wave height, wave angle, and whether is in the active zone
1313 double const dWaveHeight = VdWaveHeight[nProfilePoint];
1314 double const dWaveAngle = VdWaveDirection[nProfilePoint];
1315
1316 bBreaking = VbWaveIsBreaking[nProfilePoint];
1317
1318 // And store the wave properties for this point in the all-profiles vectors
1319 pVdX->push_back(nX);
1320 pVdY->push_back(nY);
1321 pVdHeightX->push_back(dWaveHeight * sin(dWaveAngle * PI / 180));
1322 pVdHeightY->push_back(dWaveHeight * cos(dWaveAngle * PI / 180));
1323 pVbBreaking->push_back(bBreaking);
1324 }
1325
1326 // Obtain the profile nodes near the coast
1327 int const nX = pProfile->pPtiGetCellInProfile(nProfileSize - 2)->nGetX();
1328 int const nY = pProfile->pPtiGetCellInProfile(nProfileSize - 2)->nGetY();
1329 int const nX1 = pProfile->pPtiGetCellInProfile(nProfileSize - 1)->nGetX();
1330 int const nY1 = pProfile->pPtiGetCellInProfile(nProfileSize - 1)->nGetY();
1331
1332 // Calculate the horizontal distance between the profile points
1333 double dXDist = tAbs(dGridCentroidXToExtCRSX(nX1) - dGridCentroidXToExtCRSX(nX));
1334 double dYDist = tAbs(dGridCentroidYToExtCRSY(nY1) - dGridCentroidYToExtCRSY(nY));
1335
1336 // Safety check
1337 if (bFPIsEqual(dXDist, 0.0, TOLERANCE))
1338 dXDist = 1e-3;
1339
1340 // Safety check
1341 if (bFPIsEqual(dYDist, 0.0, TOLERANCE))
1342 dYDist = 1e-3;
1343
1344 double const dDiffProfileDistXY = hypot(dXDist, dYDist);
1345
1346 // Compute the beach slope
1347 double const dtanBeta = tan(tAbs(m_pRasterGrid->m_Cell[nX1][nY1].dGetSeaDepth() - m_pRasterGrid->m_Cell[nX][nY].dGetSeaDepth()) / dDiffProfileDistXY);
1348
1349 // Compute the wave run-up using NIELSEN & HANSLOW (1991) & DHI (2004)
1350 int nValidPointsWaveHeight = 0;
1351 int nValidPointsWaveSetup = 0;
1352
1353 for (int nPoint = 0; nPoint < static_cast<int>(VdWaveHeight.size()); nPoint++)
1354 {
1355 if (VdWaveHeight[nPoint] > 1e-4)
1356 {
1357 nValidPointsWaveHeight += 1;
1358 }
1359
1360 else
1361 {
1362 break;
1363 }
1364 }
1365
1366 nValidPointsWaveHeight -= 1;
1367
1368 for (int nPoint = 0; nPoint < static_cast<int>(VdWaveSetupSurge.size()); nPoint++)
1369 {
1370 if (tAbs(VdWaveSetupSurge[nPoint]) < 1) // limiting the absolute value of setup + surge if cshore run fails
1371 {
1372 nValidPointsWaveSetup += 1;
1373 }
1374
1375 else
1376 {
1377 break;
1378 }
1379 }
1380
1381 nValidPointsWaveSetup -= 1;
1382
1383 double dWaveHeight = 0;
1384
1385 // Safety checks
1386 if ((nValidPointsWaveHeight >= 0) && (!bFPIsEqual(VdWaveHeight[nValidPointsWaveHeight], DBL_NODATA, TOLERANCE)))
1387 {
1388 dWaveHeight = VdWaveHeight[nValidPointsWaveHeight];
1389 }
1390
1391 // TODO 060 Remove these 'magic numbers'
1392 double dRunUp = 0;
1393
1394 if (m_nRunUpEquation == 0)
1395 {
1396 // Compute the run-up using Nielsen & Hanslow (1991) & DHI (2004)
1397 dRunUp = 0.36 * pow(9.81, 0.5) * dtanBeta * pow(dWaveHeight, 0.5) * dDeepWaterWavePeriod;
1398 }
1399
1400 else if (m_nRunUpEquation == 1)
1401 {
1402 // Compute the run-up using MASE 1989
1403 double const dS0 = 2 * PI * dWaveHeight / (9.81 * dDeepWaterWavePeriod * dDeepWaterWavePeriod);
1404 dRunUp = 1.86 * dWaveHeight * pow(pow(dtanBeta / dS0, 0.5), 0.71);
1405 }
1406
1407 else if (m_nRunUpEquation == 2)
1408 {
1409 // Compute the run-up using STOCKDON (2006)
1410 double const dS0 = 2 * PI * dWaveHeight / (9.81 * dDeepWaterWavePeriod * dDeepWaterWavePeriod);
1411 // dRunUp = 1.1 * ((0.35 * dWaveHeight * (pow((1 / dS0) * dWaveHeight * dWaveHeight, 0.5))) + (((((1 / dS0) * dWaveHeight * dWaveHeight) * (0.563 * dWaveHeight * dWaveHeight + 0.0004)), 0.5)) / 2);
1412
1413 double const dH0OverL0 = (1 / dS0) * dWaveHeight;
1414 double const dTmp1 = 0.35 * dWaveHeight * pow(dH0OverL0, 0.5);
1415 double const dTmp2 = pow(dH0OverL0 * ((0.563 * dWaveHeight * dWaveHeight) + 0.0004), 0.5);
1416 dRunUp = 1.1 * (dTmp1 + (dTmp2 / 2));
1417 }
1418
1419 if ((tAbs(dRunUp) < 1e-4) || (isnan(dRunUp)))
1420 {
1421 dRunUp = 0;
1422 }
1423
1424 double dWaveSetupSurge = 0;
1425
1426 // Safety checks
1427 if ((nValidPointsWaveSetup >= 0) && (!bFPIsEqual(VdWaveSetupSurge[nValidPointsWaveSetup], DBL_NODATA, TOLERANCE)))
1428 {
1429 dWaveSetupSurge = VdWaveSetupSurge[nValidPointsWaveSetup];
1430 }
1431
1432 if ((tAbs(dWaveSetupSurge) < 1e-4) || (isnan(dWaveSetupSurge)))
1433 {
1434 dWaveSetupSurge = 0;
1435 }
1436
1437 // Update wave attributes along the coastline object. Wave height at the coast is always calculated (i.e. whether or not waves are breaking)
1438 // cout << "Wave Height at the coast is " << VdWaveHeight[nProfileSize - 1] << endl;
1439 m_VCoast[nCoast].SetCoastWaveHeight(nCoastPoint, dWaveHeight);
1440 m_VCoast[nCoast].SetWaveSetupSurge(nCoastPoint, dWaveSetupSurge);
1441 m_VCoast[nCoast].SetRunUp(nCoastPoint, dRunUp);
1442
1443 if (nProfileBreakingDist > 0)
1444 {
1445 // This coast point is in the active zone, so set breaking wave height, breaking wave angle, and depth of breaking for the coast point
1446 m_VCoast[nCoast].SetBreakingWaveHeight(nCoastPoint, dProfileBreakingWaveHeight);
1447 m_VCoast[nCoast].SetBreakingWaveAngle(nCoastPoint, dProfileBreakingWaveAngle);
1448 m_VCoast[nCoast].SetDepthOfBreaking(nCoastPoint, dProfileBreakingDepth);
1449 m_VCoast[nCoast].SetBreakingDistance(nCoastPoint, nProfileBreakingDist);
1450
1451 // LogStream << m_ulIter << ": nProfile = " << nProfile << ", nCoastPoint = " << nCoastPoint << " in active zone, dBreakingWaveHeight = " << dBreakingWaveHeight << endl;
1452 }
1453
1454 else
1455 {
1456 // This coast point is not in the active zone, so breaking wave height, breaking wave angle, and depth of breaking are all meaningless
1457 m_VCoast[nCoast].SetBreakingWaveHeight(nCoastPoint, DBL_NODATA);
1458 m_VCoast[nCoast].SetBreakingWaveAngle(nCoastPoint, DBL_NODATA);
1459 m_VCoast[nCoast].SetDepthOfBreaking(nCoastPoint, DBL_NODATA);
1460 m_VCoast[nCoast].SetBreakingDistance(nCoastPoint, INT_NODATA);
1461
1462 // LogStream << m_ulIter << ": nProfile = " << nProfile << ", nCoastPoint = " << nCoastPoint << " NOT in active zone" << endl;
1463 }
1464
1465 return RTN_OK;
1466}
1467
1468#if defined CSHORE_FILE_INOUT
1469//===============================================================================================================================
1471//===============================================================================================================================
1472int CSimulation::nCreateCShoreInfile(int const nCoast, int const nProfile, int const nILine, int const nIProfl, int const nIPerm, int const nIOver, int const nIWcint, int const nIRoll, int const nIWind, int const nITide, int const nILab, int const nWave, int const nSurge, double const dX, double const dTimestep, double const dWaveInitTime, double const dWavePeriod, double const dHrms, double const dWaveAngle, double const dSurgeInitTime, double const dSurgeLevel, vector<double> const *pVdXdist, vector<double> const *pVdBottomElevation, vector<double> const *pVdWaveFriction)
1473{
1474 // Create the CShore input file
1475 ofstream CShoreOutStream;
1476 CShoreOutStream.open(CSHORE_INFILE.c_str(), ios::out | ios::app);
1477
1478 if (CShoreOutStream.fail())
1479 {
1480 // Error, cannot open file for writing
1481 LogStream << m_ulIter << ": " << ERR << "cannot write to CShore input file '" << CSHORE_INFILE << "'" << endl;
1483 }
1484
1485 // And write to the file
1486 CShoreOutStream << 3 << endl; // Number of comment lines
1487 CShoreOutStream << "------------------------------------------------------------" << endl;
1488 CShoreOutStream << "CShore input file created by CoastalME for iteration " << m_ulIter << ", coast " << nCoast << ", profile " << nProfile << endl;
1489 CShoreOutStream << "------------------------------------------------------------" << endl;
1490 CShoreOutStream << nILine << " -> ILINE" << endl;
1491 CShoreOutStream << nIProfl << " -> IPROFL" << endl;
1492 CShoreOutStream << nIPerm << " -> IPERM" << endl;
1493 CShoreOutStream << nIOver << " -> IOVER" << endl;
1494 CShoreOutStream << nIWcint << " -> IWCINT" << endl;
1495 CShoreOutStream << nIRoll << " -> IROLL" << endl;
1496 CShoreOutStream << nIWind << " -> IWIND" << endl;
1497 CShoreOutStream << nITide << " -> ITIDE" << endl;
1498 CShoreOutStream << fixed;
1499 CShoreOutStream << setw(11) << setprecision(4) << dX << " -> DX" << endl;
1500 CShoreOutStream << setw(11) << m_dBreakingWaveHeightDepthRatio << " -> GAMMA" << endl;
1501 CShoreOutStream << setw(11) << nILab << " -> ILAB" << endl;
1502 CShoreOutStream << setw(11) << nWave << " -> NWAVE" << endl;
1503 CShoreOutStream << setw(11) << nSurge << " -> NSURGE" << endl;
1504
1505 // Line 18 of infile
1506 CShoreOutStream << setw(11) << setprecision(2) << dWaveInitTime; // TWAVE(1) in CShore
1507 CShoreOutStream << setw(11) << setprecision(4) << dWavePeriod; // TPIN(1) in CShore
1508 CShoreOutStream << setw(11) << dHrms; // HRMS(1) in CShore
1509 CShoreOutStream << setw(11) << dWaveAngle << endl; // WANGIN(1) in CShore
1510
1511 // Line 19 of infile
1512 CShoreOutStream << setw(11) << setprecision(2) << dTimestep; // TWAVE(2) in CShore
1513 CShoreOutStream << setw(11) << setprecision(4) << dWavePeriod; // TPIN(2) in CShore
1514 CShoreOutStream << setw(11) << dHrms; // HRMS(2) in CShore
1515 CShoreOutStream << setw(11) << dWaveAngle << endl; // WANGIN(2) in CShore
1516
1517 // Line 20 of infile
1518 CShoreOutStream << setw(11) << setprecision(2) << dSurgeInitTime; // TSURG(1) in CShore
1519 CShoreOutStream << setw(11) << setprecision(4) << dSurgeLevel << endl; // SWLIN(1) in CShore
1520
1521 // Line 21 of infile
1522 CShoreOutStream << setw(11) << setprecision(2) << dTimestep; // TSURG(2) in CShore
1523 CShoreOutStream << setw(11) << setprecision(4) << dSurgeLevel << endl; // SWLIN(2) in CShore
1524
1525 // Line 22 of infile
1526 CShoreOutStream << setw(8) << pVdXdist->size() << " -> NBINP" << endl;
1527
1528 CShoreOutStream << fixed << setprecision(4);
1529
1530 for (unsigned int i = 0; i < pVdXdist->size(); i++)
1531 // These are BINP(J,1), ZBINP(J,1), FBINP(J-1,1) in CShore
1532 CShoreOutStream << setw(11) << pVdXdist->at(i) << setw(11) << pVdBottomElevation->at(i) << setw(11) << pVdWaveFriction->at(i) << endl;
1533
1534 CShoreOutStream << endl;
1535
1536 // File written, so close it
1537 CShoreOutStream.close();
1538
1539 return RTN_OK;
1540}
1541#endif
1542
1543//===============================================================================================================================
1545//===============================================================================================================================
1546int CSimulation::nGetThisProfileElevationsForCShore(int const nCoast, CGeomProfile *pProfile, int const nProfSize, vector<double> *VdDistXY, vector<double> *VdVZ, vector<double> *VdFricF)
1547{
1548 bool bIsBehindIntervention = false;
1549
1550 int nX1 = 0;
1551 int nY1 = 0;
1552
1553 double dXDist;
1554 double dYDist;
1555 double dProfileDistXY = 0;
1556 double dProfileFricFact;
1557 double dPrevDist = -1;
1558
1559 for (int i = nProfSize - 1; i >= 0; i--)
1560 {
1561 int const nX = pProfile->pPtiVGetCellsInProfile()->at(i).nGetX();
1562 int const nY = pProfile->pPtiVGetCellsInProfile()->at(i).nGetY();
1563
1564 // Calculate the horizontal distance relative to the most seaward point
1565 if (i == nProfSize - 1)
1566 dProfileDistXY = 0;
1567
1568 else
1569 {
1572 dProfileDistXY = dProfileDistXY + hypot(dXDist, dYDist);
1573 }
1574
1575 // Before we store the X-Y distance, must check that it is not the same as the previously-stored distance (if it is, we get zero-divide errors in CShore). If they are the same, add on a small distance
1576 if (bFPIsEqual(dProfileDistXY, dPrevDist, TOLERANCE))
1577 dProfileDistXY += 0.1; // TODO 084 Improve this
1578
1579 // Update the cell indexes, the initial cell is now the previous one
1580 nX1 = nX;
1581 nY1 = nY;
1582
1583 // Get the number of the highest layer with non-zero thickness
1584 int const nTopLayer = m_pRasterGrid->m_Cell[nX][nY].nGetTopNonZeroLayerAboveBasement();
1585
1586 // Safety checks
1587 if (nTopLayer == INT_NODATA)
1588 return RTN_ERR_NO_TOP_LAYER;
1589
1590 if (nTopLayer == NO_NONZERO_THICKNESS_LAYERS)
1591 // TODO 009 We are down to basement, decide what to do
1592 return RTN_OK;
1593
1594 // Get the elevation for both consolidated and unconsolidated sediment on this cell
1595 double const dTopElev = m_pRasterGrid->m_Cell[nX][nY].dGetSedimentTopElev() + m_pRasterGrid->m_Cell[nX][nY].dGetInterventionHeight();
1596 double const VdProfileZ = dTopElev - m_dThisIterSWL;
1597
1598 // Check that landward elevation is greater than SWL
1599 if (i == 0)
1600 {
1601 if (VdProfileZ < 0)
1602 {
1603 VdVZ->push_back(0.1); // TODO 053 Set it to a small +ve elevation (compared with SWL). However there must be a better way of doing this
1604
1605 // Could not create the profile elevation vectors
1606 LogStream << m_ulIter << ": " << WARN << "for coast " << nCoast << ", profile " << pProfile->nGetCoastID() << ", elevation at the landward end is " << dTopElev << " m. This is lower than this-iteration SWL (" << m_dThisIterSWL << " m). For CShore, changing the landward elevation for profile " << pProfile->nGetCoastID() << " to " << m_dThisIterSWL + 0.1 << "m" << endl;
1607 }
1608
1609 else
1610 {
1611 VdVZ->push_back(VdProfileZ);
1612 }
1613 }
1614
1615 else
1616 {
1617 VdVZ->push_back(VdProfileZ);
1618 }
1619
1620 // Now store the X-Y plane distance from the start of the profile
1621 VdDistXY->push_back(dProfileDistXY);
1622
1623 // Get the landform type at each point along the profile
1624 double const dInterventionHeight = m_pRasterGrid->m_Cell[nX][nY].dGetInterventionHeight();
1625
1626 // Modify default friction factor if a structural intervention is found, otherwise use the default
1627 if (dInterventionHeight > 0 || bIsBehindIntervention)
1628 {
1629 // Use an arbitrarily high value if structure is present
1630 dProfileFricFact = 100 * CSHORE_FRICTION_FACTOR;
1631 bIsBehindIntervention = true;
1632 }
1633
1634 else
1635 dProfileFricFact = CSHORE_FRICTION_FACTOR;
1636
1637 // Store the friction factor
1638 VdFricF->push_back(dProfileFricFact);
1639
1640 // For next time round
1641 dPrevDist = dProfileDistXY;
1642 }
1643
1644 return RTN_OK;
1645}
1646
1647#if defined CSHORE_FILE_INOUT
1648//===============================================================================================================================
1650//===============================================================================================================================
1651int CSimulation::nReadCShoreOutput(int const nProfile, string const *strCShoreFilename, int const nExpectedColumns, int const nCShorecolumn, vector<double> const *pVdProfileDistXYCME, vector<double> *pVdInterpolatedValues)
1652{
1653 // Read in the first column (contains XY distance relative to seaward limit) and CShore column from the CShore output file
1654 ifstream InStream;
1655 InStream.open(strCShoreFilename->c_str(), ios::in);
1656
1657 // Did it open OK?
1658 if (!InStream.is_open())
1659 {
1660 // Error: cannot open CShore file for input
1661 LogStream << m_ulIter << ": " << ERR << "for profile " << nProfile << ", cannot open " << *strCShoreFilename << " for input" << endl;
1662
1664 }
1665
1666 // Opened OK, so set up the vectors to hold the CShore output data
1667 vector<double> VdXYDistCShore;
1668 vector<double> VdValuesCShore;
1669
1670 // And read in the data
1671 int n = -1;
1672 int nExpectedRows = 0;
1673 string strLineIn;
1674
1675 while (getline(InStream, strLineIn))
1676 {
1677 n++;
1678
1679 if (n == 0)
1680 {
1681 // Read in the header line
1682 vector<string> VstrItems = VstrSplit(&strLineIn, SPACE);
1683
1684 if (!bIsStringValidInt(VstrItems[1]))
1685 {
1686 string strErr = ERR + "invalid integer for number of expected rows '" + VstrItems[1] + "' in " + *strCShoreFilename + "\n";
1687 cerr << strErr;
1688 LogStream << strErr;
1689
1691 }
1692
1693 // Get the number of expected rows
1694 nExpectedRows = stoi(VstrItems[1]);
1695 }
1696
1697 else
1698 {
1699 // Read in a data line
1700 vector<string> VstrItems = VstrSplit(&strLineIn, SPACE);
1701
1702 int nCols = static_cast<int>(VstrItems.size());
1703
1704 if (nCols != nExpectedColumns)
1705 {
1706 // Error: did not read the expected number of CShore output columns
1707 LogStream << m_ulIter << ": " << ERR << "for profile " << nProfile << ", expected " << nExpectedColumns << " CShore output columns but read " << nCols << " columns from header section of file " << *strCShoreFilename << endl;
1708
1710 }
1711
1712 // Number of columns is OK
1713 VdXYDistCShore.push_back(strtod(VstrItems[0].c_str(), NULL));
1714 VdValuesCShore.push_back(strtod(VstrItems[nCShorecolumn - 1].c_str(), NULL));
1715 }
1716 }
1717
1718 // Check that we have read nExpectedRows from the file
1719 int nReadRows = static_cast<int>(VdXYDistCShore.size());
1720
1721 if (nReadRows != nExpectedRows)
1722 {
1723 // Error: did not get nExpectedRows CShore output rows
1724 LogStream << m_ulIter << ": " << ERR << "for profile " << nProfile << ", expected " << nExpectedRows << " CShore output rows, but read " << nReadRows << " rows from file " << *strCShoreFilename << endl;
1725
1727 }
1728
1729 if (nReadRows < 2)
1730 {
1731 // CShore sometimes returns only one row, which contains data for the seaward point of the profile. This happens when all other (more coastward) points give an invalid result during CShore's calculations. This is a problem. We don't want to abandon the simulation just because of this, so instead we just duplicate the row, so that the profile will later get marked as invalid
1733 LogStream << m_ulIter << ": " << WARN << "for profile " << nProfile << ", only " << nReadRows << " CShore output rows in file " << *strCShoreFilename << endl;
1734
1735 // Duplicate the data
1736 VdXYDistCShore.push_back(VdXYDistCShore[0]);
1737 VdValuesCShore.push_back(VdValuesCShore[0]);
1738
1739 // And increase the expected number of rows
1740 nReadRows++;
1741 }
1742
1743 // The output is OK, so change the origin of the across-shore distance from the CShore convention to the one used here (i.e. with the origin at the shoreline)
1744 vector<double> VdXYDistCShoreTmp(nReadRows, 0);
1745
1746 for (int i = 0; i < nReadRows; i++)
1747 VdXYDistCShoreTmp[i] = VdXYDistCShore[nReadRows - 1] - VdXYDistCShore[i];
1748
1749 // Reverse the XY-distance and value vectors (i.e. first point is at the shoreline and must be in strictly ascending order)
1750 reverse(VdXYDistCShoreTmp.begin(), VdXYDistCShoreTmp.end());
1751
1752 // Similarly, reverse the CShore output
1753 reverse(VdValuesCShore.begin(), VdValuesCShore.end());
1754
1755 // Using a simple linear interpolation approach
1756 vector<double> VdDistXYCopy(pVdProfileDistXYCME->begin(), pVdProfileDistXYCME->end());
1757
1758 // assertVdXYDistCShoreTmp.size() == VdValuesCShore.size());
1759 *pVdInterpolatedValues = VdInterpolateCShoreProfileOutput(&VdXYDistCShoreTmp, &VdValuesCShore, &VdDistXYCopy);
1760
1761 return RTN_OK;
1762}
1763#endif
1764
1765#if defined CSHORE_ARG_INOUT || CSHORE_BOTH
1766//===============================================================================================================================
1768//===============================================================================================================================
1769void CSimulation::InterpolateCShoreOutput(vector<double> const *pVdProfileDistXYCME, int const nOutSize, int const nProfileSize, vector<double> *pVdXYDistFromCShore, vector<double> *pVdFreeSurfaceStdCShore, vector<double> *pVdWaveSetupSurgeCShore, vector<double> *pVdSinWaveAngleRadiansCShore, vector<double> *pVdFractionBreakingWavesCShore, vector<double> *pVdFreeSurfaceStdCME, vector<double> *pVdWaveSetupSurgeCME, vector<double> *pVdSinWaveAngleRadiansCME, vector<double> *pVdFractionBreakingWavesCME)
1770{
1771 // // DEBUG CODE ========================================================================================================
1772 // LogStream << m_ulIter << ": nOutSize = " << nOutSize << " nProfileSize = " << nProfileSize << " pVdProfileDistXYCME->size() = " << pVdProfileDistXYCME->size() << " pVdXYDistFromCShore->size() = " << pVdXYDistFromCShore->size() << endl;
1773 // // DEBUG CODE ========================================================================================================
1774
1775 // Sometimes the profile length returned from CShore is shorter than the CoastalME profile length
1776 bool bTooShort = false;
1777 int nTooShort = 0;
1778
1779 if (nOutSize < nProfileSize)
1780 {
1781 bTooShort = true;
1782 nTooShort = nProfileSize - nOutSize;
1783
1784 // LogStream << m_ulIter << ": CShore PROFILE IS TOO SHORT" << endl;
1785 }
1786
1787 // // DEBUG CODE ========================================================================================================
1788 // for (int n = 0; n < static_cast<int>(pVdProfileDistXYCME->size()); n++)
1789 // LogStream << "pVdProfileDistXYCME[" << n << "] = " << pVdProfileDistXYCME->at(n) << endl;
1790 //
1791 // LogStream << endl;
1792 //
1793 // LogStream << "ORIGINAL" << endl;
1794 // for (int n = 0; n < nOutSize; n++)
1795 // LogStream << "pVdXYDistFromCShore[" << n << "] = " << pVdXYDistFromCShore->at(n) << " pVdFreeSurfaceStdCShore[" << n << "] = " << pVdFreeSurfaceStdCShore->at(n) << " pVdWaveSetupSurgeCShore[" << n << "] = " << pVdWaveSetupSurgeCShore->at(n) << " pVdSinWaveAngleRadiansCShore[" << n << "] = " << pVdSinWaveAngleRadiansCShore->at(n) << " pVdFractionBreakingWavesCShore[" << n << "] = " << pVdFractionBreakingWavesCShore->at(n) << endl;
1796 //
1797 // LogStream << endl;
1798 // // DEBUG CODE ========================================================================================================
1799
1800 if (bTooShort)
1801 {
1802 // Add extrapolated value(s) to the end of the valid part of each profile vector so that we have nProfileSize valid values
1803 double dLastDiff = pVdXYDistFromCShore->at(nOutSize - 1) - pVdXYDistFromCShore->at(nOutSize - 2);
1804
1805 for (int n = 0; n < nTooShort; n++)
1806 pVdXYDistFromCShore->at(nOutSize + n) = pVdXYDistFromCShore->at(nOutSize - 1 + n) + dLastDiff;
1807
1808 for (int n = 0; n < nTooShort; n++)
1809 pVdFreeSurfaceStdCShore->at(nOutSize + n) = pVdFreeSurfaceStdCShore->at(nOutSize - 1);
1810
1811 for (int n = 0; n < nTooShort; n++)
1812 pVdWaveSetupSurgeCShore->at(nOutSize + n) = pVdWaveSetupSurgeCShore->at(nOutSize - 1);
1813
1814 // TODO 007 Do same for pVdStormSurgeCShore and pVdWaveSetupRunUpCShore ?
1815
1816 for (int n = 0; n < nTooShort; n++)
1817 pVdSinWaveAngleRadiansCShore->at(nOutSize + n) = pVdSinWaveAngleRadiansCShore->at(nOutSize - 1);
1818
1819 dLastDiff = pVdFractionBreakingWavesCShore->at(nOutSize - 1) - pVdFractionBreakingWavesCShore->at(nOutSize - 2);
1820
1821 for (int n = 0; n < nTooShort; n++)
1822 pVdFractionBreakingWavesCShore->at(nOutSize + n) = tMin(pVdFractionBreakingWavesCShore->at(nOutSize - 1 + n) + dLastDiff, 1.0);
1823
1824 // // DEBUG CODE ========================================================================================================
1825 // LogStream << "EXTENDED" << endl;
1826 // for (int n = 0; n < nProfileSize; n++)
1827 // LogStream << "pVdXYDistFromCShore[" << n << "] = " << pVdXYDistFromCShore->at(n) << " pVdFreeSurfaceStdCShore[" << n << "] = " << pVdFreeSurfaceStdCShore->at(n) << " pVdWaveSetupSurgeCShore[" << n << "] = " << pVdWaveSetupSurgeCShore->at(n) << " pVdSinWaveAngleRadiansCShore[" << n << "] = " << pVdSinWaveAngleRadiansCShore->at(n) << " pVdFractionBreakingWavesCShore[" << n << "] = " << pVdFractionBreakingWavesCShore->at(n) << " " << (n > (nOutSize-1) ? "EXTRAPOLATED" : "") << endl;
1828 //
1829 // LogStream << endl;
1830 // // DEBUG CODE ========================================================================================================
1831 }
1832
1833 vector<double> VdXYDistCShoreTmp(nProfileSize);
1834 copy(pVdXYDistFromCShore->begin(), pVdXYDistFromCShore->begin() + nProfileSize, begin(VdXYDistCShoreTmp));
1835
1836 // The CShore cross-shore distance has its origin at the seaward end, but the CoastalME convention has its origin at the shoreline. So we mst reverse the other vectors to conform with the CoastalMEME convention
1837 vector<double> VdFreeSurfaceStdCShoreTmp(nProfileSize);
1838 reverse_copy(pVdFreeSurfaceStdCShore->begin(), pVdFreeSurfaceStdCShore->begin() + nProfileSize, begin(VdFreeSurfaceStdCShoreTmp));
1839
1840 vector<double> VdWaveSetupSurgeCShoreTmp(nProfileSize);
1841 reverse_copy(pVdWaveSetupSurgeCShore->begin(), pVdWaveSetupSurgeCShore->begin() + nProfileSize, begin(VdWaveSetupSurgeCShoreTmp));
1842
1843 vector<double> VdSinWaveAngleRadiansCShoreTmp(nProfileSize);
1844 reverse_copy(pVdSinWaveAngleRadiansCShore->begin(), pVdSinWaveAngleRadiansCShore->begin() + nProfileSize, begin(VdSinWaveAngleRadiansCShoreTmp));
1845
1846 vector<double> VdFractionBreakingWavesCShoreTmp(nProfileSize);
1847 reverse_copy(pVdFractionBreakingWavesCShore->begin(), pVdFractionBreakingWavesCShore->begin() + nProfileSize, begin(VdFractionBreakingWavesCShoreTmp));
1848
1849 // // DEBUG CODE ========================================================================================================
1850 // LogStream << "REVERSED" << endl;
1851 // for (int n = 0; n < nProfileSize; n++)
1852 // LogStream << "VdXYDistCShoreTmp[" << n << "] = " << VdXYDistCShoreTmp.at(n) << " VdFreeSurfaceStdCShoreTmp[" << n << "] = " << VdFreeSurfaceStdCShoreTmp.at(n) << " VdWaveSetupSurgeCShoreTmp[" << n << "] = " << VdWaveSetupSurgeCShoreTmp.at(n) << " VdSinWaveAngleRadiansCShoreTmp[" << n << "] = " << VdSinWaveAngleRadiansCShoreTmp.at(n) << " VdFractionBreakingWavesCShoreTmp[" << n << "] = " << VdFractionBreakingWavesCShoreTmp.at(n) << endl;
1853 //
1854 // LogStream << endl;
1855 // // DEBUG CODE ========================================================================================================
1856
1857 // Finally we linearly interpolate the CShore output vectors to each point on the CoastalME profile. Input parameters x_old, y_old, x_new and the routine returns y_new
1858 *pVdFreeSurfaceStdCME = VdInterpolateCShoreProfileOutput(&VdXYDistCShoreTmp, pVdFreeSurfaceStdCShore, pVdProfileDistXYCME); // was &VdDistXYCopy);
1859 *pVdWaveSetupSurgeCME = VdInterpolateCShoreProfileOutput(&VdXYDistCShoreTmp, pVdWaveSetupSurgeCShore, pVdProfileDistXYCME);
1860 // *pVdStormSurgeCME = VdInterpolateCShoreProfileOutput(&VdXYDistCShoreTmp, pVdStormSurgeCShore, pVdProfileDistXYCME);
1861 // *pVdWaveSetupRunUpCME = VdInterpolateCShoreProfileOutput(&VdXYDistCShoreTmp, pVdWaveSetupRunUpCShore, pVdProfileDistXYCME);
1862 *pVdSinWaveAngleRadiansCME = VdInterpolateCShoreProfileOutput(&VdXYDistCShoreTmp, pVdSinWaveAngleRadiansCShore, pVdProfileDistXYCME);
1863 *pVdFractionBreakingWavesCME = VdInterpolateCShoreProfileOutput(&VdXYDistCShoreTmp, pVdFractionBreakingWavesCShore, pVdProfileDistXYCME);
1864
1865 // LogStream << "INTERPOLATED" << endl;
1866 // for (int n = 0; n < nProfileSize; n++)
1867 // LogStream << "pVdProfileDistXYCME[" << n << "] = " << pVdProfileDistXYCME->at(n) << " pVdFreeSurfaceStdCME[" << n << "] = " << pVdFreeSurfaceStdCME->at(n) << " pVdWaveSetupSurgeCME[" << n << "] = " << pVdWaveSetupSurgeCME->at(n) << " pVdSinWaveAngleRadiansCME[" << n << "] = " << pVdSinWaveAngleRadiansCME->at(n) << " pVdFractionBreakingWavesCME[" << n << "] = " << pVdFractionBreakingWavesCME->at(n) << endl;
1868 //
1869 // LogStream << "================================================ " << endl;
1870}
1871#endif
1872
1873//===============================================================================================================================
1875//===============================================================================================================================
1877{
1878 CGeomProfile *pProfile = m_VCoast[nCoast].pGetProfile(nProfile);
1879
1880 // Only do this for profiles without problems, including the start and end-of-coast profile
1881 if (!pProfile->bOKIncStartAndEndOfCoast())
1882 return;
1883
1884 bool bModfiedWaveHeightisBreaking = false;
1885 bool bProfileIsinShadowZone = false;
1886 int const nThisCoastPoint = pProfile->nGetCoastPoint();
1887 int const nProfileSize = pProfile->nGetNumCellsInProfile();
1888 int nThisBreakingDist = m_VCoast[nCoast].nGetBreakingDistance(nThisCoastPoint);
1889 double dThisBreakingWaveHeight = m_VCoast[nCoast].dGetBreakingWaveHeight(nThisCoastPoint); // This could be DBL_NODATA
1890 double dThisBreakingWaveAngle = m_VCoast[nCoast].dGetBreakingWaveAngle(nThisCoastPoint);
1891 double dThisBreakingDepth = m_VCoast[nCoast].dGetDepthOfBreaking(nThisCoastPoint);
1892
1893 // Traverse the profile landwards, checking if any profile cell is within the shadow zone
1894 for (int nProfilePoint = (nProfileSize - 1); nProfilePoint >= 0; nProfilePoint--)
1895 {
1896 int const nX = pProfile->pPtiGetCellInProfile(nProfilePoint)->nGetX();
1897 int const nY = pProfile->pPtiGetCellInProfile(nProfilePoint)->nGetY();
1898
1899 // If there is any cell profile within the shadow zone and waves are breaking then modify wave breaking properties otherwise continue
1900 if (m_pRasterGrid->m_Cell[nX][nY].bIsinAnyShadowZone())
1901 {
1902 bProfileIsinShadowZone = true;
1903
1904 // Check if the new wave height is breaking
1905 double const dWaveHeight = m_pRasterGrid->m_Cell[nX][nY].dGetWaveHeight();
1906
1907 // Check that wave height at the given point is lower than maximum real wave height. If breaking wave height is expected that no good wave height are obtained, so, do not take it
1908 if (dWaveHeight > (m_dDepthOfClosure * m_dBreakingWaveHeightDepthRatio) && (!bModfiedWaveHeightisBreaking) && (!bFPIsEqual(dThisBreakingWaveHeight, DBL_NODATA, TOLERANCE)))
1909 {
1910 // It is breaking
1911 bModfiedWaveHeightisBreaking = true;
1912
1913 dThisBreakingWaveHeight = m_dDepthOfClosure * m_dBreakingWaveHeightDepthRatio;
1914 dThisBreakingWaveAngle = m_pRasterGrid->m_Cell[nX][nY].dGetWaveAngle();
1915 dThisBreakingDepth = m_dDepthOfClosure;
1916 nThisBreakingDist = nProfilePoint;
1917 }
1918 }
1919 }
1920
1921 // Update breaking wave properties along coastal line object (Wave height, dir, distance). TODO 010 Update the active zone cells
1922 if (bProfileIsinShadowZone && bModfiedWaveHeightisBreaking) // Modified wave height is still breaking
1923 {
1924 // This coast point is in the active zone, so set breaking wave height, breaking wave angle, and depth of breaking for the coast point TODO 007 Where does the 0.78 come from? TODO 011 Should it be an input variable or a named constant?
1925 if (dThisBreakingWaveHeight > dThisBreakingDepth * 0.78)
1926 {
1927 dThisBreakingWaveHeight = dThisBreakingDepth * 0.78; // Likely CShore output wave height is not adequately reproduced due to input profile and wave properties. TODO 007 Info needed. Does something need to be changed then?
1928 }
1929
1930 // assert(dThisBreakingWaveHeight >= 0);
1931 m_VCoast[nCoast].SetBreakingWaveHeight(nThisCoastPoint, dThisBreakingWaveHeight);
1932 m_VCoast[nCoast].SetBreakingWaveAngle(nThisCoastPoint, dThisBreakingWaveAngle);
1933 m_VCoast[nCoast].SetDepthOfBreaking(nThisCoastPoint, dThisBreakingDepth);
1934 m_VCoast[nCoast].SetBreakingDistance(nThisCoastPoint, nThisBreakingDist);
1935
1936 // LogStream << m_ulIter << ": nProfile = " << nProfile << ", nCoastPoint = " << nCoastPoint << " in active zone, dBreakingWaveHeight = " << dBreakingWaveHeight << endl;
1937 }
1938
1939 else if (bProfileIsinShadowZone && (!bModfiedWaveHeightisBreaking))
1940 {
1941 // This coast point is no longer in the active zone
1942 m_VCoast[nCoast].SetBreakingWaveHeight(nThisCoastPoint, DBL_NODATA);
1943 m_VCoast[nCoast].SetBreakingWaveAngle(nThisCoastPoint, DBL_NODATA);
1944 m_VCoast[nCoast].SetDepthOfBreaking(nThisCoastPoint, DBL_NODATA);
1945 m_VCoast[nCoast].SetBreakingDistance(nThisCoastPoint, INT_NODATA);
1946
1947 // LogStream << m_ulIter << ": nProfile = " << nProfile << ", nCoastPoint = " << nCoastPoint << " NOT in active zone" << endl;
1948 }
1949
1950 return;
1951}
1952
1953//===============================================================================================================================
1955//===============================================================================================================================
1956void CSimulation::InterpolateWavePropertiesBetweenProfiles(int const nCoast, int const nCount)
1957{
1958 CGeomProfile *pProfile = m_VCoast[nCoast].pGetProfileWithDownCoastSeq(nCount);
1959
1960 // Only do this for profiles without problems, including the start-of-coast profile (but not the end-of-coast profile)
1961 // if (!pProfile->bOKIncStartOfCoast())
1962 // return;
1963
1964 int const nThisCoastPoint = pProfile->nGetCoastPoint();
1965
1966 // For the breaking wave stuff, to go into the in-between coastline points
1967 int const nThisBreakingDist = m_VCoast[nCoast].nGetBreakingDistance(nThisCoastPoint);
1968 double const dThisBreakingWaveHeight = m_VCoast[nCoast].dGetBreakingWaveHeight(nThisCoastPoint); // This could be DBL_NODATA
1969 double const dThisBreakingWaveAngle = m_VCoast[nCoast].dGetBreakingWaveAngle(nThisCoastPoint);
1970 double const dThisBreakingDepth = m_VCoast[nCoast].dGetDepthOfBreaking(nThisCoastPoint);
1971 double const dThisWaveSetupSurge = m_VCoast[nCoast].dGetWaveSetupSurge(nThisCoastPoint);
1972 // double dThisStormSurge = m_VCoast[nCoast].dGetStormSurge(nThisCoastPoint);
1973 double const dThisRunUp = m_VCoast[nCoast].dGetRunUp(nThisCoastPoint);
1974
1975 // Get the next profile along the coast, in the down-coast direction. If this next profile has a problem, go to the one after that, etc
1976 CGeomProfile const *pTmpProfile = pProfile;
1977 CGeomProfile *pNextProfile;
1978
1979 while (true)
1980 {
1981 pNextProfile = pTmpProfile->pGetDownCoastAdjacentProfile();
1982
1983 if (pNextProfile->bOKIncStartAndEndOfCoast())
1984 {
1985 // The next profile is OK
1986 break;
1987 }
1988
1989 // The next profile is not OK, so prepare to the one after that
1990 pTmpProfile = pNextProfile;
1991
1992 if (pTmpProfile->bEndOfCoast())
1993 {
1994 // Uh-oh, we've reached the down-coast end of the coast without finding an OK down-coast profile. So give up
1995 return;
1996 }
1997 }
1998
1999 // The next profile is OK
2000 int const nNextCoastPoint = pNextProfile->nGetCoastPoint();
2001 int const nDistBetween = nNextCoastPoint - nThisCoastPoint;
2002
2003 // Safety check
2004 if (nDistBetween <= 0)
2005 // Nothing to do
2006 return;
2007
2008 int const nNextBreakingDist = m_VCoast[nCoast].nGetBreakingDistance(nNextCoastPoint);
2009 double const dNextBreakingWaveHeight = m_VCoast[nCoast].dGetBreakingWaveHeight(nNextCoastPoint); // This could be DBL_NODATA
2010 double const dNextBreakingWaveAngle = m_VCoast[nCoast].dGetBreakingWaveAngle(nNextCoastPoint);
2011 double const dNextBreakingDepth = m_VCoast[nCoast].dGetDepthOfBreaking(nNextCoastPoint);
2012 double const dNextWaveSetupSurge = m_VCoast[nCoast].dGetWaveSetupSurge(nNextCoastPoint);
2013 double const dNextRunUp = m_VCoast[nCoast].dGetRunUp(nNextCoastPoint);
2014
2015 // OK, fill coast point between profiles for setupsurge and runup
2016 for (int n = nThisCoastPoint; n <= nNextCoastPoint; n++)
2017 {
2018 // Fill first wave setup and surge
2019 int const nDist = n - nThisCoastPoint;
2020 double const dThisWeight = (nDistBetween - nDist) / static_cast<double>(nDistBetween);
2021 double const dNextWeight = 1 - dThisWeight;
2022 double dWaveSetupSurge = 0;
2023 double dRunUp = 0;
2024
2025 dWaveSetupSurge = (dThisWeight * dThisWaveSetupSurge) + (dNextWeight * dNextWaveSetupSurge);
2026 m_VCoast[nCoast].SetWaveSetupSurge(n, dWaveSetupSurge);
2027
2028 dRunUp = (dThisWeight * dThisRunUp) + (dNextWeight * dNextRunUp);
2029 m_VCoast[nCoast].SetRunUp(n, dRunUp);
2030 }
2031
2032 // int const nNextProfile = pNextProfile->nGetCoastID();
2033
2034 // If both this profile and the next profile are not in the active zone, then do no more
2035 if ((bFPIsEqual(dThisBreakingWaveHeight, DBL_NODATA, TOLERANCE)) && (bFPIsEqual(dNextBreakingWaveHeight, DBL_NODATA, TOLERANCE)))
2036 {
2037 // if (m_nLogFileDetail >= LOG_FILE_HIGH_DETAIL)
2038 // LogStream << m_ulIter << ": both profile " << pProfile->nGetCoastID() << " at coast point " << nThisCoastPoint << ", and profile " << nNextProfile << " at coast point " << nNextCoastPoint << ", are not in the active zone" << endl;
2039
2040 // Set the breaking wave height, breaking wave angle, and depth of breaking to DBL_NODATA
2041 for (int n = nThisCoastPoint; n < nNextCoastPoint; n++)
2042 {
2043 m_VCoast[nCoast].SetBreakingWaveHeight(nThisCoastPoint, DBL_NODATA);
2044 m_VCoast[nCoast].SetBreakingWaveAngle(nThisCoastPoint, DBL_NODATA);
2045 m_VCoast[nCoast].SetDepthOfBreaking(nThisCoastPoint, DBL_NODATA);
2046 m_VCoast[nCoast].SetBreakingDistance(nThisCoastPoint, DBL_NODATA);
2047 }
2048
2049 return;
2050 }
2051
2052 // OK, at least one of the two profiles is in the active zone
2053 if (bFPIsEqual(dThisBreakingWaveHeight, DBL_NODATA, TOLERANCE))
2054 {
2055 // The next profile must be in the active zone, so use values from the next profile
2056 for (int n = nThisCoastPoint; n < nNextCoastPoint; n++)
2057 {
2058 // Set the breaking wave height, breaking wave angle, and depth of breaking for this coast point TODO 056 This assert sometimes fails: why?
2059 // assert(dNextBreakingWaveHeight >= 0);
2060 m_VCoast[nCoast].SetBreakingWaveHeight(n, dNextBreakingWaveHeight);
2061 m_VCoast[nCoast].SetBreakingWaveAngle(n, dNextBreakingWaveAngle);
2062 m_VCoast[nCoast].SetDepthOfBreaking(n, dNextBreakingDepth);
2063 m_VCoast[nCoast].SetBreakingDistance(n, nNextBreakingDist);
2064 }
2065
2066 return;
2067 }
2068
2069 if (bFPIsEqual(dNextBreakingWaveHeight, DBL_NODATA, TOLERANCE))
2070 {
2071 // This profile must be in the active zone, so use values from this profile
2072 for (int n = nThisCoastPoint + 1; n <= nNextCoastPoint; n++)
2073 {
2074 // Set the breaking wave height, breaking wave angle, and depth of breaking for this coast point TODO 056 This assert sometimes fails: why?
2075 // assert(dThisBreakingWaveHeight >= 0);
2076 m_VCoast[nCoast].SetBreakingWaveHeight(n, dThisBreakingWaveHeight);
2077 m_VCoast[nCoast].SetBreakingWaveAngle(n, dThisBreakingWaveAngle);
2078 m_VCoast[nCoast].SetDepthOfBreaking(n, dThisBreakingDepth);
2079 m_VCoast[nCoast].SetBreakingDistance(n, nThisBreakingDist);
2080 }
2081
2082 return;
2083 }
2084
2085 // The values for both this profile point and the next profile point are fine, so do a weighted interpolation between this profile and the next profile
2086 for (int n = nThisCoastPoint + 1; n < nNextCoastPoint; n++)
2087 {
2088 int const nDist = n - nThisCoastPoint;
2089
2090 double dBreakingWaveHeight = DBL_NODATA;
2091 double dBreakingWaveAngle = DBL_NODATA;
2092 double dBreakingDepth = DBL_NODATA;
2093 double dBreakingDist = DBL_NODATA;
2094
2095 if ((dNextBreakingDepth > 0) && (dThisBreakingDepth > 0))
2096 {
2097 double const dThisWeight = (nDistBetween - nDist) / static_cast<double>(nDistBetween);
2098 double const dNextWeight = 1 - dThisWeight;
2099
2100 dBreakingWaveHeight = (dThisWeight * dThisBreakingWaveHeight) + (dNextWeight * dNextBreakingWaveHeight),
2101 dBreakingWaveAngle = (dThisWeight * dThisBreakingWaveAngle) + (dNextWeight * dNextBreakingWaveAngle),
2102 dBreakingDepth = (dThisWeight * dThisBreakingDepth) + (dNextWeight * dNextBreakingDepth),
2103 dBreakingDist = (dThisWeight * nThisBreakingDist) + (dNextWeight * nNextBreakingDist);
2104 }
2105
2106 else if (dNextBreakingDepth > 0)
2107 {
2108 dBreakingWaveHeight = dNextBreakingWaveHeight,
2109 dBreakingWaveAngle = dNextBreakingWaveAngle,
2110 dBreakingDepth = dNextBreakingDepth,
2111 dBreakingDist = nNextBreakingDist;
2112 }
2113
2114 else if (dThisBreakingDepth > 0)
2115 {
2116 dBreakingWaveHeight = dThisBreakingWaveHeight,
2117 dBreakingWaveAngle = dThisBreakingWaveAngle,
2118 dBreakingDepth = dThisBreakingDepth,
2119 dBreakingDist = nThisBreakingDist;
2120 }
2121
2122 // Set the breaking wave height, breaking wave angle, and depth of breaking for this coast point TODO 056 This assert sometimes fails: why?
2123 // assert(dBreakingWaveHeight >= 0);
2124 m_VCoast[nCoast].SetBreakingWaveHeight(n, dBreakingWaveHeight);
2125 m_VCoast[nCoast].SetBreakingWaveAngle(n, dBreakingWaveAngle);
2126 m_VCoast[nCoast].SetDepthOfBreaking(n, dBreakingDepth);
2127 m_VCoast[nCoast].SetBreakingDistance(n, nRound(dBreakingDist));
2128 }
2129}
2130
2131//===============================================================================================================================
2133//===============================================================================================================================
2135{
2136 int const nCoastPoints = m_VCoast[nCoast].nGetCoastlineSize();
2137
2138 // Initialize all vectors pairs (x,y) for each variable
2139 vector<int> nVCoastWaveHeightX;
2140 vector<double> dVCoastWaveHeightY;
2141
2142 // Search all coast points for non NaN values and store them into temporary variables for later interporlation
2143 for (int n = 0; n < nCoastPoints; n++)
2144 {
2145 double const dCoastWaveHeight = m_VCoast[nCoast].dGetCoastWaveHeight(n);
2146
2147 if (!bFPIsEqual(dCoastWaveHeight, DBL_NODATA, TOLERANCE))
2148 {
2149 nVCoastWaveHeightX.push_back(n);
2150 dVCoastWaveHeightY.push_back(dCoastWaveHeight);
2151 }
2152 }
2153
2154 // Interpolate all coast points. Check first that x,y have more than 3 points and are both of equal size
2155 if ((nVCoastWaveHeightX.size() >= 3) && (nVCoastWaveHeightX.size() == dVCoastWaveHeightY.size()))
2156 {
2157 for (int n = 0; n < nCoastPoints; n++)
2158 {
2159 double const dInterpCoastWaveHeight = dGetInterpolatedValue(&nVCoastWaveHeightX, &dVCoastWaveHeightY, n, false);
2160 m_VCoast[nCoast].SetCoastWaveHeight(n, dInterpCoastWaveHeight);
2161 }
2162 }
2163
2164 else
2165 {
2166 for (int n = 0; n < nCoastPoints; n++)
2167 {
2168 m_VCoast[nCoast].SetCoastWaveHeight(n, 0);
2169 }
2170 }
2171}
2172
2173//===============================================================================================================================
2175//===============================================================================================================================
2177{
2178 int const nCoastSize = m_VCoast[nCoast].nGetCoastlineSize();
2179
2180 for (int nCoastPoint = 0; nCoastPoint < nCoastSize; nCoastPoint++)
2181 {
2182 double dXDiff;
2183 double dYDiff;
2184
2185 if (nCoastPoint == 0)
2186 {
2187 // For the point at the start of the coastline: use the straight line from 'this' point to the next point
2188 CGeom2DPoint const PtThis = *m_VCoast[nCoast].pPtGetCoastlinePointExtCRS(nCoastPoint); // In external CRS
2189 CGeom2DPoint const PtAfter = *m_VCoast[nCoast].pPtGetCoastlinePointExtCRS(nCoastPoint + 1); // In external CRS
2190
2191 dXDiff = PtAfter.dGetX() - PtThis.dGetX();
2192 dYDiff = PtAfter.dGetY() - PtThis.dGetY();
2193 }
2194
2195 else if (nCoastPoint == nCoastSize - 1)
2196 {
2197 // For the point at the end of the coastline: use the straight line from the point before to 'this' point
2198 CGeom2DPoint const PtBefore = *m_VCoast[nCoast].pPtGetCoastlinePointExtCRS(nCoastPoint - 1); // In external CRS
2199 CGeom2DPoint const PtThis = *m_VCoast[nCoast].pPtGetCoastlinePointExtCRS(nCoastPoint); // In external CRS
2200
2201 dXDiff = PtThis.dGetX() - PtBefore.dGetX();
2202 dYDiff = PtThis.dGetY() - PtBefore.dGetY();
2203 }
2204
2205 else
2206 {
2207 // For coastline points not at the start or end of the coast: start with a straight line which links the coastline points before and after 'this' coastline point
2208 CGeom2DPoint const PtBefore = *m_VCoast[nCoast].pPtGetCoastlinePointExtCRS(nCoastPoint - 1); // In external CRS
2209 CGeom2DPoint const PtAfter = *m_VCoast[nCoast].pPtGetCoastlinePointExtCRS(nCoastPoint + 1); // In external CRS
2210
2211 dXDiff = PtAfter.dGetX() - PtBefore.dGetX();
2212 dYDiff = PtAfter.dGetY() - PtBefore.dGetY();
2213 }
2214
2215 // Calculate angle between line and north point, measured clockwise (the azimuth)
2216 if (bFPIsEqual(dYDiff, 0.0, TOLERANCE))
2217 {
2218 // The linking line runs either W-E or E-W
2219 if (dXDiff > 0)
2220 // It runs W to E
2221 m_VCoast[nCoast].SetFluxOrientation(nCoastPoint, 90);
2222
2223 else
2224 // It runs E to W
2225 m_VCoast[nCoast].SetFluxOrientation(nCoastPoint, 270);
2226 }
2227
2228 else if (bFPIsEqual(dXDiff, 0.0, TOLERANCE))
2229 {
2230 // The linking line runs N-S or S-N
2231 if (dYDiff > 0)
2232 // It runs S to N
2233 m_VCoast[nCoast].SetFluxOrientation(nCoastPoint, 0);
2234
2235 else
2236 // It runs N to S
2237 m_VCoast[nCoast].SetFluxOrientation(nCoastPoint, 180);
2238 }
2239
2240 else
2241 {
2242 // The linking line runs neither W-E nor N-S so we have to work a bit harder to find the angle between it and the azimuth
2243 double dAzimuthAngle;
2244
2245 if (dXDiff > 0)
2246 dAzimuthAngle = (180 / PI) * (PI * 0.5 - atan(dYDiff / dXDiff));
2247
2248 else
2249 dAzimuthAngle = (180 / PI) * (PI * 1.5 - atan(dYDiff / dXDiff));
2250
2251 m_VCoast[nCoast].SetFluxOrientation(nCoastPoint, dAzimuthAngle);
2252 }
2253
2254 // LogStream << m_ulIter << ": nCoastPoint = " << nCoastPoint << " FluxOrientation = " << m_VCoast[nCoast].dGetFluxOrientation(nCoastPoint) << endl;
2255 }
2256}
2257
2258//===============================================================================================================================
2260//===============================================================================================================================
2262{
2263 vector<int> VnPolygonD50Count(m_nNumPolygonGlobal + 1, 0); // TODO 044
2264 vector<double> VdPolygonD50(m_nNumPolygonGlobal + 1, 0); // TODO 044
2265
2266 for (int nX = 0; nX < m_nXGridSize; nX++)
2267 {
2268 for (int nY = 0; nY < m_nYGridSize; nY++)
2269 {
2270 if (m_pRasterGrid->m_Cell[nX][nY].bIsInContiguousSea())
2271 {
2272 // This is a sea cell, first get polygon ID
2273 int const nPolyID = m_pRasterGrid->m_Cell[nX][nY].nGetPolygonID();
2274
2275 // Is it in the active zone?
2276 bool const bActive = m_pRasterGrid->m_Cell[nX][nY].bIsInActiveZone();
2277
2278 if (bActive)
2279 {
2280 // It is in the active zone. Does it have unconsolidated sediment on it? Test this using the UnconD50 value: if dGetUnconsD50() returns DBL_NODATA, there is no unconsolidated sediment
2281 double const dTmpd50 = m_pRasterGrid->m_Cell[nX][nY].dGetUnconsD50();
2282
2283 if (!bFPIsEqual(dTmpd50, DBL_NODATA, TOLERANCE))
2284 {
2285 // It does have unconsolidated sediment, so which polygon is this cell in?
2286 if (nPolyID != INT_NODATA)
2287 {
2288 VnPolygonD50Count[nPolyID]++;
2289 VdPolygonD50[nPolyID] += dTmpd50;
2290 }
2291 }
2292 }
2293
2294 // // Now fill in wave calc holes
2295 // if (m_pRasterGrid->m_Cell[nX][nY].dGetWaveHeight() == 0)
2296 // {
2297 // if (nID == INT_NODATA)
2298 // m_pRasterGrid->m_Cell[nX][nY].SetWaveHeight(m_dAllCellsDeepWaterWaveHeight);
2299 // }
2300 //
2301 // if (m_pRasterGrid->m_Cell[nX][nY].dGetWaveAngle() == 0)
2302 // {
2303 // if (nID == INT_NODATA)
2304 // m_pRasterGrid->m_Cell[nX][nY].SetWaveAngle(m_dAllCellsDeepWaterWaveAngle);
2305 // }
2306
2307 // Next look at the cell's N-S and W-E neighbours
2308 int nXTmp;
2309 int nYTmp;
2310 int nActive = 0;
2311 int nShadow = 0;
2312 int nShadowNum = 0;
2313 int nDownDrift = 0;
2314 int nDownDriftNum = 0;
2315 int nCoast = 0;
2316 int nRead = 0;
2317 double dWaveHeight = 0;
2318 double dWaveAngle = 0;
2319
2320 // North
2321 nXTmp = nX;
2322 nYTmp = nY - 1;
2323
2324 if (bIsWithinValidGrid(nXTmp, nYTmp))
2325 {
2326 if (m_pRasterGrid->m_Cell[nXTmp][nYTmp].bIsInContiguousSea())
2327 {
2328 nRead++;
2329 dWaveHeight += m_pRasterGrid->m_Cell[nXTmp][nYTmp].dGetWaveHeight();
2330 dWaveAngle += (m_pRasterGrid->m_Cell[nXTmp][nYTmp].dGetWaveAngle());
2331
2332 if (m_pRasterGrid->m_Cell[nXTmp][nYTmp].bIsInActiveZone())
2333 nActive++;
2334
2335 int nTmp = m_pRasterGrid->m_Cell[nXTmp][nYTmp].nGetShadowZoneNumber();
2336
2337 if (nTmp != 0)
2338 {
2339 nShadow++;
2340 nShadowNum = nTmp;
2341 }
2342
2343 nTmp = m_pRasterGrid->m_Cell[nXTmp][nYTmp].nGetDownDriftZoneNumber();
2344
2345 if (nTmp > 0)
2346 {
2347 nDownDrift++;
2348 nDownDriftNum = nTmp;
2349 }
2350 }
2351
2352 else
2353 nCoast++;
2354 }
2355
2356 // East
2357 nXTmp = nX + 1;
2358 nYTmp = nY;
2359
2360 if (bIsWithinValidGrid(nXTmp, nYTmp))
2361 {
2362 if (m_pRasterGrid->m_Cell[nXTmp][nYTmp].bIsInContiguousSea())
2363 {
2364 nRead++;
2365 dWaveHeight += m_pRasterGrid->m_Cell[nXTmp][nYTmp].dGetWaveHeight();
2366 dWaveAngle += (m_pRasterGrid->m_Cell[nXTmp][nYTmp].dGetWaveAngle());
2367
2368 if (m_pRasterGrid->m_Cell[nXTmp][nYTmp].bIsInActiveZone())
2369 nActive++;
2370
2371 int nTmp = m_pRasterGrid->m_Cell[nXTmp][nYTmp].nGetShadowZoneNumber();
2372
2373 if (nTmp != 0)
2374 {
2375 nShadow++;
2376 nShadowNum = nTmp;
2377 }
2378
2379 nTmp = m_pRasterGrid->m_Cell[nXTmp][nYTmp].nGetDownDriftZoneNumber();
2380
2381 if (nTmp > 0)
2382 {
2383 nDownDrift++;
2384 nDownDriftNum = nTmp;
2385 }
2386 }
2387
2388 else
2389 nCoast++;
2390 }
2391
2392 // South
2393 nXTmp = nX;
2394 nYTmp = nY + 1;
2395
2396 if (bIsWithinValidGrid(nXTmp, nYTmp))
2397 {
2398 if (m_pRasterGrid->m_Cell[nXTmp][nYTmp].bIsInContiguousSea())
2399 {
2400 nRead++;
2401 dWaveHeight += m_pRasterGrid->m_Cell[nXTmp][nYTmp].dGetWaveHeight();
2402 dWaveAngle += (m_pRasterGrid->m_Cell[nXTmp][nYTmp].dGetWaveAngle());
2403
2404 if (m_pRasterGrid->m_Cell[nXTmp][nYTmp].bIsInActiveZone())
2405 nActive++;
2406
2407 int nTmp = m_pRasterGrid->m_Cell[nXTmp][nYTmp].nGetShadowZoneNumber();
2408
2409 if (nTmp != 0)
2410 {
2411 nShadow++;
2412 nShadowNum = nTmp;
2413 }
2414
2415 nTmp = m_pRasterGrid->m_Cell[nXTmp][nYTmp].nGetDownDriftZoneNumber();
2416
2417 if (nTmp > 0)
2418 {
2419 nDownDrift++;
2420 nDownDriftNum = nTmp;
2421 }
2422 }
2423
2424 else
2425 nCoast++;
2426 }
2427
2428 // West
2429 nXTmp = nX - 1;
2430 nYTmp = nY;
2431
2432 if (bIsWithinValidGrid(nXTmp, nYTmp))
2433 {
2434 if (m_pRasterGrid->m_Cell[nXTmp][nYTmp].bIsInContiguousSea())
2435 {
2436 nRead++;
2437 dWaveHeight += m_pRasterGrid->m_Cell[nXTmp][nYTmp].dGetWaveHeight();
2438 dWaveAngle += (m_pRasterGrid->m_Cell[nXTmp][nYTmp].dGetWaveAngle());
2439
2440 if (m_pRasterGrid->m_Cell[nXTmp][nYTmp].bIsInActiveZone())
2441 nActive++;
2442
2443 int nTmp = m_pRasterGrid->m_Cell[nXTmp][nYTmp].nGetShadowZoneNumber();
2444
2445 if (nTmp != 0)
2446 {
2447 nShadow++;
2448 nShadowNum = nTmp;
2449 }
2450
2451 nTmp = m_pRasterGrid->m_Cell[nXTmp][nYTmp].nGetDownDriftZoneNumber();
2452
2453 if (nTmp > 0)
2454 {
2455 nDownDrift++;
2456 nDownDriftNum = nTmp;
2457 }
2458 }
2459
2460 else
2461 nCoast++;
2462 }
2463
2464 if (nRead > 0)
2465 {
2466 // Calculate the average of neighbours
2467 dWaveHeight /= nRead;
2468 dWaveAngle /= nRead;
2469 dWaveAngle = dKeepWithin360(dWaveAngle);
2470
2471 // If this sea cell has four active-zone neighbours, then it must also be in the active zone: give it wave height and orientation which is the average of its neighbours
2472 if (nActive == 4)
2473 {
2474 m_pRasterGrid->m_Cell[nX][nY].SetInActiveZone(true);
2475 m_pRasterGrid->m_Cell[nX][nY].SetWaveHeight(dWaveHeight);
2476 m_pRasterGrid->m_Cell[nX][nY].SetWaveAngle(dWaveAngle);
2477 }
2478
2479 // If this sea cell has a wave height which is the same as its deep-water wave height, but its neighbours have a different average wave height, then give it the average of its neighbours
2480 double const dDeepWaterWaveHeight = m_pRasterGrid->m_Cell[nX][nY].dGetCellDeepWaterWaveHeight();
2481
2482 if ((bFPIsEqual(m_pRasterGrid->m_Cell[nX][nY].dGetWaveHeight(), dDeepWaterWaveHeight, TOLERANCE)) && (!bFPIsEqual(dDeepWaterWaveHeight, dWaveHeight, TOLERANCE)))
2483 {
2484 m_pRasterGrid->m_Cell[nX][nY].SetWaveHeight(dWaveHeight);
2485 }
2486
2487 // If this sea cell has a wave orientation which is the same as its deep-water wave orientation, but its neighbours have a different average wave orientation, then give it the average of its neighbours
2488 double const dDeepWaterWaveAngle = m_pRasterGrid->m_Cell[nX][nY].dGetCellDeepWaterWaveAngle();
2489
2490 if ((bFPIsEqual(m_pRasterGrid->m_Cell[nX][nY].dGetWaveAngle(), dDeepWaterWaveAngle, TOLERANCE)) && (!bFPIsEqual(dDeepWaterWaveAngle, dWaveAngle, TOLERANCE)))
2491 {
2492 m_pRasterGrid->m_Cell[nX][nY].SetWaveAngle(dWaveAngle);
2493 }
2494
2495 // Is this sea cell is not already marked as in a shadow zone (note could be marked as in a shadow zone but not yet processed: a -ve number)?
2496 int const nShadowZoneCode = m_pRasterGrid->m_Cell[nX][nY].nGetShadowZoneNumber();
2497
2498 if (nShadowZoneCode <= 0)
2499 {
2500 // If the cell has four neighbours which are all in a shadow zone, or four neighbours some of which are shadow zone and the remainder downdrift zone, or four neighbours some of which are shadow zone and the remainder coast; then it should also be in the shadow zone: give it the average of its neighbours
2501 if ((nShadow == 4) || (nShadow + nDownDrift == 4) || (nShadow + nCoast == 4))
2502 {
2503 m_pRasterGrid->m_Cell[nX][nY].SetShadowZoneNumber(nShadowNum);
2504 m_pRasterGrid->m_Cell[nX][nY].SetWaveHeight(dWaveHeight);
2505 m_pRasterGrid->m_Cell[nX][nY].SetWaveAngle(dWaveAngle);
2506 }
2507 }
2508
2509 // If this sea cell is not marked as in a downdrift zone but has four neighbours which are in a downdrift zone, then it should also be in the downdrift zone: give it the average of its neighbours
2510 int const nDownDriftZoneCode = m_pRasterGrid->m_Cell[nX][nY].nGetDownDriftZoneNumber();
2511
2512 if ((nDownDriftZoneCode == 0) && (nDownDrift == 4))
2513 {
2514 m_pRasterGrid->m_Cell[nX][nY].SetDownDriftZoneNumber(nDownDriftNum);
2515 m_pRasterGrid->m_Cell[nX][nY].SetWaveHeight(dWaveHeight);
2516 m_pRasterGrid->m_Cell[nX][nY].SetWaveAngle(dWaveAngle);
2517 }
2518 }
2519 }
2520 }
2521 }
2522
2523 // Calculate the average d50 for every polygon TODO 044
2524 for (int nCoast = 0; nCoast < static_cast<int>(m_VCoast.size()); nCoast++)
2525 {
2526 for (int nPoly = 0; nPoly < m_VCoast[nCoast].nGetNumPolygons(); nPoly++)
2527 {
2528 CGeomCoastPolygon *pPolygon = m_VCoast[nCoast].pGetPolygon(nPoly);
2529 int const nPolyID = pPolygon->nGetCoastID();
2530
2531 if (VnPolygonD50Count[nPolyID] > 0)
2532 VdPolygonD50[nPolyID] /= VnPolygonD50Count[nPolyID];
2533
2534 pPolygon->SetAvgUnconsD50(VdPolygonD50[nPolyID]);
2535 }
2536 }
2537}
Contains CGeom2DPoint definitions.
int nGetY(void) const
Returns the CGeom2DIPoint object's integer Y coordinate.
Definition 2di_point.cpp:55
int nGetX(void) const
Returns the CGeom2DIPoint object's integer X coordinate.
Definition 2di_point.cpp:49
Geometry class used to represent 2D point objects with floating-point coordinates.
Definition 2d_point.h:27
double dGetY(void) const
Returns the CGeom2DPoint object's double Y coordinate.
Definition 2d_point.cpp:53
double dGetX(void) const
Returns the CGeom2DPoint object's double X coordinate.
Definition 2d_point.cpp:47
Geometry class used for coast polygon objects.
void SetAvgUnconsD50(double const)
Set the average d50 for unconsolidated sediment on this polygon.
int nGetCoastID(void) const
Get the coast ID, this is the same as the down-coast sequence of polygons.
Geometry class used to represent coast profile objects.
Definition profile.h:37
double dGetProfileDeepWaterWaveAngle(void) const
Returns the deep-water wave orientation for this profile.
Definition profile.cpp:602
double dGetProfileDeepWaterWaveHeight(void) const
Returns the deep-water wave height for this profile.
Definition profile.cpp:590
int nGetCoastPoint(void) const
Returns the coast point at which the profile starts.
Definition profile.cpp:87
CGeomProfile * pGetDownCoastAdjacentProfile(void) const
Definition profile.cpp:478
CGeom2DIPoint * pPtiGetCellInProfile(int const)
Returns a single cell in the profile.
Definition profile.cpp:511
int nGetNumCellsInProfile(void) const
Returns the number of cells in the profile.
Definition profile.cpp:525
void SetCShoreProblem(bool const)
Sets a switch to indicate whether this profile has a CShore problem.
Definition profile.cpp:135
bool bOKIncStartAndEndOfCoast(void) const
Returns true if this is a problem-free profile (however it could be a start-of-coast or an end-of-coa...
Definition profile.cpp:254
int nGetCoastID(void) const
Returns the profile's coast ID.
Definition profile.cpp:75
bool bStartOfCoast(void) const
Returns the switch to indicate whether this is a start-of-coast profile.
Definition profile.cpp:117
vector< CGeom2DIPoint > * pPtiVGetCellsInProfile(void)
Returns all cells in the profile.
Definition profile.cpp:504
bool bEndOfCoast(void) const
Returns the switch to indicate whether this is an end-of-coast profile.
Definition profile.cpp:129
double dGetProfileDeepWaterWavePeriod(void) const
Returns the deep-water wave Period for this profile.
Definition profile.cpp:614
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:574
double m_dAllCellsDeepWaterWaveHeight
Deep water wave height (m) for all sea cells.
Definition simulation.h:754
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 > *)
double m_dG
Gravitational acceleration (m**2/sec)
Definition simulation.h:811
double m_dThisIterSWL
The still water level for this timestep (this includes tidal changes and any long-term SWL change)
Definition simulation.h:712
CGeomRasterGrid * m_pRasterGrid
Pointer to the raster grid object.
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 ...
int m_nXGridSize
The size of the grid in the x direction.
Definition simulation.h:457
double m_dL_0
Deep water wave length (m)
Definition simulation.h:745
ofstream LogStream
vector< CRWCoast > m_VCoast
The coastline objects.
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:385
int m_nYGridSize
The size of the grid in the y direction.
Definition simulation.h:460
int nDoAllShadowZones(void)
Finds wave shadow zones and modifies waves in and near them. Note that where up-coast and down-coast ...
static double dKeepWithin360(double const)
Constrains the supplied angle to be within 0 and 360 degrees.
void CalcCoastTangents(int const)
Calculates tangents to a coastline: the tangent is assumed to be the orientation of energy/sediment f...
double m_dWaveDepthRatioForWaveCalcs
Start depth for wave calculations.
Definition simulation.h:748
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 *)
void InterpolateWavePropertiesBetweenProfiles(int const, int const)
Interpolates wave properties from profiles to the coastline points between two profiles....
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(...
int m_nRunUpEquation
The run-up equation used TODO 007.
Definition simulation.h:577
double m_dBreakingWaveHeightDepthRatio
Breaking wave height-to-depth ratio.
Definition simulation.h:751
int nReadCShoreOutput(int const, string const *, int const, int const, vector< double > const *, vector< double > *)
double dGridCentroidYToExtCRSY(int const) const
Definition gis_utils.cpp:98
void CalcD50AndFillWaveCalcHoles(void)
Calculates an average d50 for each polygon. Also fills in 'holes' in active zone and wave calcs i....
string m_strCMEDir
The CME folder.
vector< double > VdInterpolateCShoreProfileOutput(vector< double > const *, vector< double > const *, vector< double > const *)
Returns a linearly interpolated vector of doubles, to make CShore profile output compatible with CME....
int m_nWavePropagationModel
The wave propagation model used. Possible values are WAVE_MODEL_CSHORE and WAVE_MODEL_COVE.
Definition simulation.h:550
double m_dAllCellsDeepWaterWaveAngle
Deep water wave angle for all sea cells.
Definition simulation.h:757
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_dTimeStep
The length of an iteration (a time step) in hours.
Definition simulation.h:676
int nGetThisProfileElevationsForCShore(int const, CGeomProfile *, int const, vector< double > *, vector< double > *, vector< double > *)
Get profile horizontal distance and bottom elevation vectors in CShore units.
bool bIsWithinValidGrid(int const, int const) const
void ModifyBreakingWavePropertiesWithinShadowZoneToCoastline(int const, int const)
Modifies the wave breaking properties at coastline points of profiles within the shadow zone.
double m_dC_0
Deep water wave speed (m/s)
Definition simulation.h:742
unsigned long m_ulIter
The number of the current iteration (time step)
Definition simulation.h:598
double dGridCentroidXToExtCRSX(int const) const
Definition gis_utils.cpp:87
int m_nNumPolygonGlobal
Number of global (all coasts) polygons.
Definition simulation.h:523
double dGetInterpolatedValue(vector< double > const *, vector< double > const *, double, bool)
int nDoAllPropagateWaves(void)
Simulates wave propagation along all coastline-normal profiles, on all coasts.
double m_dDepthOfClosure
Depth of closure (in m) TODO 007 can be calculated using Hallermeier, R.J. (1978) or Birkemeier (1985...
Definition simulation.h:817
void InterpolateWaveHeightToCoastPoints(int const)
Linearly interpolates wave properties from profiles to the coastline cells between two profiles.
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:2443
static double dCalcWaveAngleToCoastNormal(double const, double const, int const)
Calculates the angle between the wave direction and a normal to the coastline tangent....
This file contains global definitions for CoastalME.
int const NO_NONZERO_THICKNESS_LAYERS
Definition cme.h:769
int const INT_NODATA
Definition cme.h:476
T tMin(T a, T b)
Definition cme.h:1255
int const WAVE_MODEL_COVE
Definition cme.h:782
double const TOLERANCE
Definition cme.h:813
int const RTN_ERR_NO_TOP_LAYER
Definition cme.h:740
string const CSHORE_DIR
Definition cme.h:889
int const RTN_ERR_CSHORE_FILE_INPUT
Definition cme.h:748
string const ERR
Definition cme.h:892
double const CSHORE_FRICTION_FACTOR
Definition cme.h:810
int const WAVE_MODEL_CSHORE
Definition cme.h:783
int const LOG_FILE_MIDDLE_DETAIL
Definition cme.h:491
double const PI
Definition cme.h:795
int const RTN_ERR_CSHORE_ERROR
Definition cme.h:753
bool bFPIsEqual(const T d1, const T d2, const T dEpsilon)
Definition cme.h:1293
double const WALKDEN_HALL_PARAM_2
Definition cme.h:803
T tMax(T a, T b)
Definition cme.h:1242
int const RTN_ERR_READING_CSHORE_FILE_OUTPUT
Definition cme.h:749
double const CSHORE_SURGE_LEVEL
Definition cme.h:811
string const WARN
Definition cme.h:893
int const RTN_OK
Definition cme.h:694
bool const SAVE_CSHORE_OUTPUT
Definition cme.h:459
double const DBL_NODATA
Definition cme.h:824
int const CSHOREARRAYOUTSIZE
The size of the arrays output by CShore. If you change this, then you must also set the same value on...
Definition cme.h:472
double const WALKDEN_HALL_PARAM_1
Definition cme.h:802
T tAbs(T a)
Definition cme.h:1267
int const RTN_ERR_CSHORE_EMPTY_PROFILE
Definition cme.h:747
string const CSHORE_INFILE
Definition cme.h:890
int const LEFT_HANDED
Definition cme.h:512
char const SPACE
Definition cme.h:453
Contains CRWCoast definitions.
void CShore(int *)
void CShoreWrapper(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[], double[], double[], double[], double[], double[], int const *, double[], double[], double[], int *, int *, double[], double[], double[], double[], double[])
Contains CSimulation definitions.
bool bIsStringValidInt(string &str)
Checks to see if a string can be read as a valid integer, from https://stackoverflow....
int nRound(double const d)
Version of the above that returns an int.