-
Notifications
You must be signed in to change notification settings - Fork 2
/
calv46sn.cpp
1693 lines (1553 loc) · 86.8 KB
/
calv46sn.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/* Copyright 2017 Lambert Rubash
This file is part of TopNetCpp, a translation and enhancement of
Fortran TopNet.
TopNetCpp 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.
TopNetCpp 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.
You should have received a copy of the GNU General Public License
along with TopNetCpp. If not, see <http://www.gnu.org/licenses/>.
*/
#include "topnet.hh"
#include <sstream>
#include <iomanip>
#include "snow.hh"
using namespace std;
ofstream oFile[51], extraFiles[13];
ofstream snowcontrol3_File; // unit 10
// This version,V3, has lakes and snow modelling added to it
// ******************************************************
// * subroutine calcts
// ******************************************************
int calcts( double **Si, const vector<vector<double> > &Sp, double **Rp, const valarray<int> &ll,
const int Nsub, const valarray<int> &Nka, const valarray<double> &tl, double **atb,
double **pka, const valarray<int> &nd, double **cl, double **pd,
const double units, const int ipsub, const int ipatb, const bool reinit,
const bool modwrt, const int stim, vector<vector<double> > &bRain, const long int interval,
const int m, const int mi, const int mps, const int mpe,
bool &ok, const valarray<double> &xlat, const valarray<double> &xlong, const double stdlon,
const valarray<double> &elevtg, double **bdtBar, const int sDate, int &sHour,
const valarray<double> &temper, const valarray<double> &dewp, const valarray<double> &tRange, const int Neq,
const int Nout, const int nBout, const valarray<int> &iBout, const valarray<double> &wind2m,
double **bTmin, double **bTmax, double **bTdew, const valarray<double> &bXlat,
const valarray<double> &bXlon, int *ntdh, const int ndump, const int maxInt,
const int maxSlp, const int maxA, const int maxC, const int maxChn,
const int idebugoutput, const int idebugbasin, const int idebugcase)
{
istringstream line;
char dateStr[11];
double daysInYear = -1.0;
const int offset = 18; // Offset needed for file name array.
const int withdrawalOffset = 39; // Offset for withdrawal by user type file output code.
const int NuserTypes = 6;
ostringstream TcodeStr, AcodeStr;
int ncode;
int userTypeCode[] = {constant_definitions::SoilMoistureIrrigationUseCode, // 1
constant_definitions::FixedDemandIrrigationUseCode, // 2
constant_definitions::DownstreamReservoirReleaseUseCode, // 3
constant_definitions::PWSUseCode, // 4
constant_definitions::NonPWSMandIUseCode, // 5
constant_definitions::DairyUseCode, // 6
constant_definitions::RanchUseCode, // 7
constant_definitions::PoultryUseCode, // 8
constant_definitions::ParkGolfCemeteryUseCode, // 9
constant_definitions::InstreamFlowUseCode, // 10
constant_definitions::DiversionUseCode}; // 11
// The Model
int ievap_method, n_irrig, n_drainage, kk, i_irrig, i_drainage, kp;
double wt0, wt1, wt2, wt12, depth_irrig_dem, rate_irrig, art_drainage;
double dep, qinst_out_0, dr_out_d, dr_out_0;
double zbm_d, acsem_d, s1_d, sr_d, cv_d, bal_d, qinst_out_d;
double zbar_d;
double sumr_d, sumq_d, sumae_d, s0_d, sumpe_d, dth1;
double art_drainage_out; // artificial drainage
double scalefactor;
// DGT 5/27/12 allocating arrays that seemed to be allocated dynamically before
double *groundwater_to_take; // groundwater_to_take(maxSlp)
double evap_mm, qlat_mm;
const int nreg = 6; // how many regions can we model within a sub-basin?
// use these regions for irrigation and drainage
//flood comment out next line
static vector<int> irr_l(Nip1);
static array<double,Nip1> rirr;
static int istep, j, n;
// locals
int js;
// ET variables
static double albedo, elevsb, rlapse;
static double PETsngl, elev;
static double temper_temp, dewp_temp, trange_temp;
// snow
static double ddf;
static int irad;
// REACH ROUTING ( ROUTE )
static int jr, jr1;
// locals:
static double prec, pet=0.0;
static double *snowst;
int jsub, i;
static int jj;
static double smin;
// This is special for CALCTS
static int iyear, month = 0, iday, ihr, imm, isec, ihh, iss;
static double hour1;
static double **sr, **cv;
static double *q0, **s0;
static double **tdh;
static double **aciem;
static double **acsem;
static double **zbm;
static double **sumr;
static double **sumae;
static double **sumpe;
static double **sumq;
static double **sumie;
static double **sumse;
static double **sumqb;
static double **sumce;
static double sumad;
static double **sumsle;
static double **sumr1;
static double **sumqv;
static double *qb;
static double *zr;
static double *ak0fzrdt;
static double *logoqm;
static double *qvmin;
static double *dth;
static int ngut, natball;
// SNOWUEB
static int ndepletionpoints, ipflag, mQgOption;
static int nintstep, isurftmpoption, nstepday, isnow_method;
static double bca, bcc, a, b, timestep;
double surfacewaterinput;
static double snowevaporation, areafractionsnow;
static double ta, p, ws, qsiobs, qnetob;
static string inFile, outFileName, pFile, svFile, bcFile;
static string aFile, dfcFile, outFileDet;
static string inLine;
static array<double,Nsv> snowsitev;
//double **snowstatev;
static array<double,Niv> snowforcing;
static vector<double> snowparam(Npar);
static array<int,7> snowcontrol;
static array<double,12> dtbar;
static double *cump;
static double *cume;
static double *cummr;
static double *errmbal;
static double *w1;
// arrays to keep records of surface temperature and snowpack average
// temperature. This is for the fourth model (Modified force restore approach)
static double **dfc;
#if TRACE
static int ncalls = 0;
double tm0 = static_cast<double>(clock())/static_cast<double>(CLOCKS_PER_SEC);
string save_caller = caller;
if (ncalls < MAX_TRACE) {
traceFile << setw(30) << caller << " -> calcts(" << ncalls << ")" << std::endl;
}
caller = "calcts";
#endif
string resultsFileNames[] = {"results/annual_upwelling_mm.txt", // 18
"results/annual_recharge_mm.txt", // 19
"results/annual_precipitation_mm.txt", // 20
"results/annual_potential_evapotranspiration_mm.txt", // 21
"results/annual_return_flow_cmy.txt", // 22, cubic meters per year
"results/monthly_upwelling_mm.txt", // 23
"results/monthly_recharge_mm.txt", // 24
"results/monthly_precipitation_mm.txt", // 25
"results/monthly_potentialevap_mm.txt", // 26
"results/monthly_return_flow_cmm.txt", // 27, cubic meters per month
"results/annual_mean_upwelling_mm.txt", // 28
"results/annual_mean_recharge_mm.txt", // 29
"results/annual_mean_precipitation_mm.txt", // 30
"results/annual_mean_potentialevap_mm.txt", // 31
"results/annual_mean_return_flow_cmy.txt", // 32, cubic meters per year
"results/annual_artificial_drainage_cmy.txt", // 33
"results/monthly_artificial_drainage_cmm.txt", // 34
"results/annual_mean_artificial_drainage_cmy.txt", // 35
"results/annual_evaporation_mm.txt", // 36
"results/monthly_evaporation_mm.txt", // 37
"results/annual_mean_evaporation_mm.txt"}; // 38
string extraFileNames[] = { "results/upwelling_irrigated_natural_drainage_mm.txt", // 0
"results/upwelling_irrigated_tile_drainage_mm.txt", // 1
"results/upwelling_irrigated_ditch_drainage_mm.txt", // 2
"results/upwelling_unirrigated_natural_drainage_mm.txt", // 3
"results/upwelling_unirrigated_tile_drainage_mm.txt", // 4
"results/upwelling_unirrigated_ditch_drainage_mm.txt", // 5
"results/recharge_irrigated_natural_drainage_mm.txt", // 6
"results/recharge_irrigated_tile_drainage_mm.txt", // 7
"results/recharge_irrigated_ditch_drainage_mm.txt", // 8
"results/recharge_unirrigated_natural_drainage_mm.txt", // 9
"results/recharge_unirrigated_tile_drainage_mm.txt", // 10
"results/recharge_unirrigated_ditch_drainage_mm.txt", // 11
"results/precip_minus_et_mm.txt"}; // 12
// cms calculation for e.g. user withdrawal is dividided by 864000 which implies that user withdrawal
// stored in the User structure is in units of cubic meters per day. See the
// 'write_output_tables.cpp' file for examples.
string drainageWithdrawalTotalsFileNames[NuserTypes], drainageWithdrawalAnnAveFileNames[NuserTypes];
build_topnet_to_client_index(); // Early initialization of TopNet to Client drainage number conversion array.
// annual, monthly and average files must be opened and initialized here because their
// write routines cannot look for timestep = 0
strftime(dateStr, 11,"%Y %m %d", timeinfo);
// --------------------------------------------------------------------------------
for (i = 0; i < 13; ++i) {
extraFiles[i].open(extraFileNames[i]);
if (!extraFiles[i].is_open()) {
cerr << "Failed to open " << extraFileNames[i] << '\n';
exit(EXIT_FAILURE);
} else {
cerr << "Extra File # " << dec << setw(2) << i << " " << extraFileNames[i] << " opened" << endl;
}
extraFiles[i] << left << dec << setw(11) << "Date";
for (jsub = 0; jsub < Nsub; jsub++) {
extraFiles[i] << setw(12) << "Drainage_" << left << dec << setw(5)
<< input_structures::Drainage[index_to_real_DID[jsub]].RealDrainageID;
}
extraFiles[i] << endl;
}
// --------------------------------------------------------------------------------
for (i = offset; i < 28; ++i) {
oFile[i].open(resultsFileNames[i-offset]);
if (!oFile[i].is_open()) {
cerr << "Failed to open " << resultsFileNames[i-offset] << '\n';
exit(EXIT_FAILURE);
} else {
cerr << "File # " << dec << setw(2) << i << " " << resultsFileNames[i-offset] << " opened" << endl;
}
oFile[i] << dec << setw(11) << "Date";
for (jsub = 0; jsub < Nsub; jsub++) {
oFile[i] << setw(10) << "Drainage_" << left << dec << setw(5)
<< input_structures::Drainage[index_to_real_DID[jsub]].RealDrainageID;;
}
oFile[i] << '\n';
}
for (i = 28; i < withdrawalOffset; ++i) {
oFile[i].open(resultsFileNames[i-offset]);
if (!oFile[i].is_open()) {
cerr << "Failed to open " << resultsFileNames[i-offset] << '\n';
exit(EXIT_FAILURE);
} else {
cerr << "File # " << dec << setw(2) << i << " " << resultsFileNames[i-offset] << " opened" << endl;
}
if (i == 34 || i == 37) {
oFile[i] << setw(11) << "Months";
} else {
oFile[i] << setw(11) << "Years";
}
for (jsub = 0; jsub < Nsub; jsub++) {
oFile[i] << setw(10) << "Drainage_" << left << dec << setw(5)
<< input_structures::Drainage[index_to_real_DID[jsub]].RealDrainageID;;
}
oFile[i] << '\n';
}
i = withdrawalOffset;
for (ncode = 0; ncode < NuserTypes; ++ncode) {
TcodeStr << "results/DrainageWithdrawalByUserType_" << userTypeCode[ncode] << "_cmy.txt"; // cubic meters per year
AcodeStr << "results/AnnAveDrainageWithdrawalByUserType_" << userTypeCode[ncode] << "_cmy.txt"; // cubic meters per year
oFile[i+ncode].open(TcodeStr.str());
oFile[i+ncode+NuserTypes].open(AcodeStr.str());
if (!oFile[i+ncode].is_open()) {
cerr << "Failed to open " << TcodeStr.str() << '\n';
exit(EXIT_FAILURE);
} else {
cerr << "File # " << dec << setw(2) << i+ncode << " " << TcodeStr.str() << " opened" << endl;
}
if (!oFile[i+ncode+NuserTypes].is_open()) {
cerr << "Failed to open " << AcodeStr.str() << '\n';
exit(EXIT_FAILURE);
} else {
cerr << "File # " << dec << setw(2) << i+ncode+NuserTypes << " " << AcodeStr.str() << " opened" << endl;
}
drainageWithdrawalTotalsFileNames[ncode] = TcodeStr.str();
drainageWithdrawalAnnAveFileNames[ncode] = AcodeStr.str();
TcodeStr.clear();
TcodeStr.str("");
AcodeStr.clear();
AcodeStr.str("");
oFile[i+ncode] << setw(11) << "Year end";
if (i+ncode+NuserTypes >= 23 &&i+ncode+NuserTypes < 28) {
oFile[i+ncode+NuserTypes] << setw(11) << "Months";
} else {
oFile[i+ncode+NuserTypes] << setw(11) << "Years";
}
for (jsub = 0; jsub < Nsub; jsub++) {
oFile[i+ncode] << setw(10) << "Drainage_" << left << dec << setw(5)
<< input_structures::Drainage[index_to_real_DID[jsub]].RealDrainageID;
oFile[i+ncode+NuserTypes] << setw(10) << "Drainage_" << left << dec << setw(5)
<< input_structures::Drainage[index_to_real_DID[jsub]].RealDrainageID;
}
oFile[i+ncode] << '\n';
oFile[i+ncode+NuserTypes] << '\n';
}
//----------------------------------------------------------------------------------------------------
oFile[16].open("results/recharge_mm.txt");
if (!oFile[16].is_open()) {
cerr << "Failed to open results/recharge_mm.txt\n";
exit(EXIT_FAILURE);
} else {
cerr << "results/recharge_mm.txt opened" << endl;
}
oFile[16] << left << dec << setw(11) << "Date";
for (jsub = 0; jsub < Nsub; jsub++) {
oFile[16] << setw(12) << "Drainage_" << left << dec << setw(5)
<< input_structures::Drainage[index_to_real_DID[jsub]].RealDrainageID;
}
oFile[16] << endl;
//----------------------------------------------------------------------------------------------------
// Allocate and initialize all the now dynamic arrays
vector<double> baseflow(maxSlp,0.0);
vector<double> totalrunoff(maxSlp,0.0);
vector<double> tempave(maxSlp,0.0);
vector<double> surfro(maxSlp,0.0);
vector<double> canstore(maxSlp,0.0);
vector<double> soilstore(maxSlp,0.0);
vector<double> tiled(maxSlp,0.0);
vector<double> ditchd(maxSlp,0.0);
valarray<double> volume_irrig_sup(0.0,maxSlp);
valarray<double> ArtDrainage(0.0,maxSlp);
valarray<double> annual_ArtDrainage(0.0,maxSlp);
valarray<double> monthly_ArtDrainage(0.0,maxSlp);
valarray<double> annAve_ArtDrainage(0.0,maxSlp);
valarray<double> irrigated_natural_upwelling(0.0,maxSlp);
valarray<double> irrigated_tile_upwelling(0.0,maxSlp);
valarray<double> irrigated_ditch_upwelling(0.0,maxSlp);
valarray<double> unirrigated_natural_upwelling(0.0,maxSlp);
valarray<double> unirrigated_tile_upwelling(0.0,maxSlp);
valarray<double> unirrigated_ditch_upwelling(0.0,maxSlp);
double upwelling_temporary;
valarray<double> upwelling(0.0,maxSlp);
valarray<double> annual_upwelling(0.0,maxSlp);
valarray<double> monthly_upwelling(0.0,maxSlp);
valarray<double> annAve_upwelling(0.0,maxSlp);
valarray<double> irrigated_natural_recharge(0.0,maxSlp);
valarray<double> irrigated_tile_recharge(0.0,maxSlp);
valarray<double> irrigated_ditch_recharge(0.0,maxSlp);
valarray<double> unirrigated_natural_recharge(0.0,maxSlp);
valarray<double> unirrigated_tile_recharge(0.0,maxSlp);
valarray<double> unirrigated_ditch_recharge(0.0,maxSlp);
double recharge_temporary;
double precip_minus_et_temporary;
valarray<double> precip_minus_et(0.0,maxSlp);
valarray<double> recharge(0.0,maxSlp);
valarray<double> annual_recharge(0.0,maxSlp);
valarray<double> monthly_recharge(0.0,maxSlp);
valarray<double> annAve_recharge(0.0,maxSlp);
valarray<double> returnflow(0.0,maxSlp);
valarray<double> GroundWaterReturnFlow(0.0,maxSlp);
valarray<double> SurfaceWaterReturnFlow(0.0,maxSlp);
valarray<double> annual_returnflow(0.0,maxSlp);
valarray<double> monthly_returnflow(0.0,maxSlp);
valarray<double> annAve_returnflow(0.0,maxSlp);
valarray<double> precip_for_watermgmt(0.0,maxSlp);
valarray<double> annual_precipitation(0.0,maxSlp);
valarray<double> monthly_precipitation(0.0,maxSlp);
valarray<double> annAve_precipitation(0.0,maxSlp);
valarray<double> potentialevap(0.0,maxSlp);
valarray<double> annual_potentialevap(0.0,maxSlp);
valarray<double> monthly_potentialevap(0.0,maxSlp);
valarray<double> annAve_potentialevap(0.0,maxSlp);
valarray<double> evap_for_watermgmt(0.0,maxSlp);
valarray<double> annual_evaporation(0.0,maxSlp);
valarray<double> monthly_evaporation(0.0,maxSlp);
valarray<double> annAve_evaporation(0.0,maxSlp);
vector<vector<double> > BasinWithdrawalByUser(NuserTypes,vector<double>(Nsub)); // first 6 user types
vector<vector<double> > annAve_basin_withdrawal(NuserTypes,vector<double>(Nsub)); // first 6 user types
int countFullYears = 0;
vector<double> vol_irrig_demand(maxSlp,0.0);
groundwater_to_take = new double[maxSlp];
for (i = 0; i < maxSlp; i++) {
groundwater_to_take[i] = 0.0;
}
vector<array<double,12> > bdtBarR(maxSlp); // maxSlp rows, 12 columns
vector<vector<double> > zbar(maxSlp,vector<double>(nreg)); // maxSlp rows, nreg columns
snowst = new double[maxSlp];
q0 = new double[maxSlp];
qb = new double[maxSlp];
qvmin = new double[maxSlp];
zr = new double[maxSlp];
ak0fzrdt = new double[maxSlp];
logoqm = new double[maxSlp];
dth = new double[maxSlp];
cump = new double[maxSlp];
cume = new double[maxSlp];
cummr = new double[maxSlp];
errmbal = new double[maxSlp];
w1 = new double[maxSlp];
sr = new double*[maxSlp];
cv = new double*[maxSlp];
s0 = new double*[maxSlp];
aciem = new double*[maxSlp];
acsem = new double*[maxSlp];
zbm = new double*[maxSlp];
sumr = new double*[maxSlp];
sumae = new double*[maxSlp];
sumpe = new double*[maxSlp];
sumq = new double*[maxSlp];
sumie = new double*[maxSlp];
sumse = new double*[maxSlp];
sumqb = new double*[maxSlp];
sumce = new double*[maxSlp];
sumsle = new double*[maxSlp];
sumr1 = new double*[maxSlp];
sumqv = new double*[maxSlp];
for (j = 0; j < maxSlp; j++) {
sr[j] = new double[nreg];
cv[j] = new double[nreg];
s0[j] = new double[nreg];
aciem[j] = new double[nreg];
acsem[j] = new double[nreg];
zbm[j] = new double[nreg];
sumr[j] = new double[nreg];
sumae[j] = new double[nreg];
sumpe[j] = new double[nreg];
sumq[j] = new double[nreg];
sumie[j] = new double[nreg];
sumse[j] = new double[nreg];
sumqb[j] = new double[nreg];
sumce[j] = new double[nreg];
sumsle[j] = new double[nreg];
sumr1[j] = new double[nreg];
sumqv[j] = new double[nreg];
}
double ***dr = new double**[maxSlp];
double ***qinst = new double**[maxSlp];
for (js = 0; js < maxSlp; ++js) {
dr[js] = new double*[nreg];
qinst[js] = new double*[nreg];
for (kk = 0; kk < nreg; ++kk) {
dr[js][kk] = new double[MAX_NTDH];
qinst[js][kk] = new double[MAX_NTDH];
}
}
vector<double> dr1(MAX_NTDH);
vector<double> qinst1(MAX_NTDH);
double***zbar_in = new double**[m+1];
for (int i = 0; i <= m; ++i) {
zbar_in[i] = new double*[Nsub];
for (int j = 0; j < Nsub; ++j) {
zbar_in[i][j] = new double[nreg];
}
}
tdh = new double*[MAX_NTDH];
for (j = 0; j <= MAX_NTDH; j++) {
tdh[j] = new double[maxSlp];
}
vector<vector<double> > snowstatev(maxSlp,vector<double>(Nxv)); // maxSlp rows, Nxv columns
// *******************************************************************************
// THIS SECTION IS RELEVANT TO THE WHOLE MODEL
// ASSUME EVERYTHING WILL TURN OUT OK
ok = true;
// Calculate sum of catchment areas above each reach
ngut = Nsub;
ievap_method = 2; //0=P-T, 1 and 2 are P-M
// initialise snowueb model
aFile = "snow.in";
ifstream snowFile(aFile.c_str()); // unit 1
if (!snowFile.is_open()) {
cerr << "Failed to open " << aFile << '\n';
exit(EXIT_FAILURE);
}
snowFile >> inFile >> outFileName >> outFileDet >> pFile >> svFile >> bcFile >> dfcFile;
outFileName = "results/" + outFileName;
snowFile >> irad;
getline(snowFile, inLine, '\n'); // Read the remainder of the line
// Flag to control radiation
// 0 is no measurements - radiation estimated from diurnal temperature range
// 1 is incoming shortwave radiation read from file (measured), incoming longwave estimated
// 2 is incoming shortwave and longwave radiation read from file (measured)
// 3 is net radiation read from file (measured)
snowFile >> ipflag; // Flag to control printing (1=print)
getline(snowFile, inLine, '\n'); // Read the remainder of the line
snowFile >> nintstep; // Number of internal time steps to use
getline(snowFile, inLine, '\n'); // Read the remainder of the line
snowFile >> isurftmpoption; // Surface temperature algorithm option
getline(snowFile, inLine, '\n'); // Read the remainder of the line
snowFile >> mQgOption; // ground heat input
getline(snowFile, inLine, '\n'); // Read the remainder of the line
snowFile.close();
ifstream snowparamFile(pFile.c_str()); // unit 88
if (!snowparamFile.is_open()) {
cerr << "Failed to open " << pFile << '\n';
exit(EXIT_FAILURE);
}
for (i = 0; i < Npar; i++) {
snowparamFile >> snowparam[i];
getline(snowparamFile, inLine, '\n'); // Read the remainder of the line
}
snowparamFile.close();
bcparm(dtbar, bca, bcc, bcFile); // warning we wont use the dtbar values we read here but use this read to get bca and bcc
// dtbar values used are read in hyData and averaged for nearby gages
ifstream dcFile(dfcFile.c_str()); // unit 88
if (!dcFile.is_open()) {
cerr << "Failed to open " << dfcFile << '\n';
exit(EXIT_FAILURE);
}
i = 0;
getline(dcFile, inLine, '\n'); // Read the first line.
while (inLine.length() > 1) {
line.str(inLine);
line >> a >> b;
i++;
getline(dcFile, inLine, '\n');
}
ndepletionpoints = i;
dfc = new double*[ndepletionpoints];
for (j = 0; j < ndepletionpoints; j++) {
dfc[j] = new double[2];
}
// rewind
dcFile.clear();
dcFile.seekg(0);
for (i = 0; i < ndepletionpoints; i++) {
dcFile >> dfc[i][0] >> dfc[i][1];
}
dcFile.close();
for (j = 0; j < maxSlp; j++) {
for (i = 0; i < Nxv; i++) {
snowstatev[j][i] = 0.0;
}
}
// month,day,year,hour,
timestep = 24;
snowcontrol[0] = irad;
snowcontrol[1] = ipflag;
snowcontrol[2] = 10;
snowcontrol[3] = nintstep;
snowcontrol[6] = isurftmpoption;
nstepday = 24/timestep*nintstep;
vector<vector<double> > snowsurfacetemp(maxSlp,vector<double>(nstepday)); // maxSlp columns
vector<vector<double> > snowaveragetemp(maxSlp,vector<double>(nstepday)); // maxSlp columns
for (j = 0; j < maxSlp; ++j) {
for (i = 0; i < nstepday; ++i) {
snowsurfacetemp[j][i] = -9999.0;
snowaveragetemp[j][i] = -9999.0;
}
}
for (j = 0; j < maxSlp; j++) {
snowsurfacetemp[j][nstepday-1] = 0.0;
snowaveragetemp[j][nstepday-1] = 0.0;
}
if (snowcontrol[1] >= 1) {
snowcontrol3_File.open(outFileName.c_str());
if (!snowcontrol3_File.is_open()) {
cerr << "Failed to open " << outFileName << '\n';
exit(EXIT_FAILURE);
}
}
// initialize variables for mass balance
for (i = 0; i < maxSlp; i++) {
w1[i] = snowstatev[i][1];
cump[i] = 0.0;
cume[i] = 0.0;
cummr[i] = 0.0;
errmbal[i] = 0.0;
}
// ***********************************************************************
// SUBCATCHMENT SECTION
// Landcare
natball = 0;
// implementing multiple raingauges for each sub-basin
// (see also top.f, mddatav4.f)
for (js = 0; js < Nsub; js++) { // for each subbasin
// Landcare
natball += Nka[js];
} // Now pass 'bRain' instead of 'rain' and use row js of bRain, rather than row links(2,js) of rain
// do 3551 it=1,m
// 3551 write(21,*)(bRain(js,it),js=1,Nsub)
// write subbasin output headers
if ( modwrt ) {
if (idebugoutput >= 1) {
lunmodFile << nBout << " " << m << '\n';
lunmodFile << " Basin";
lunmodFile << " TimeStep IrrDrainCat Afrac SWInput_mm Qlat_mm Qtot_mm Qb_mm";
lunmodFile << " Recharge_mm SatEx_mm InfEx_mm SurfRo_mm SatAfrac InfAfrac";
lunmodFile << " IntStore_mm WTDepth_mm SoilStore_mm Pet_mm Aet_mm";
lunmodFile << " Irrig_mm GWTake_mm IrrDem_mm Prec_mm SWE_mm Sublim_mm";
lunmodFile << " Tave_C Tdew_C Trange_C ErrClosure_mm\n";
// Landcare
if (mpe == -1) {
lundatFile << Nsub << '\n';
}
luntopFile << "Measured and modeled flows\n";
luntopFile << "Timestep (column 1), Measured reaches (columns 2 to " << dec << setw(3) << Neq+1;
luntopFile << ") Modeled reaches (columns " << dec << setw(4) << Neq+2 << " to " << dec << setw(4) << Neq+Nout+1 << ")\n";
luntopFile << "units ARE um/interval normalized by each basins own area";
luntopFile << " unless they have a -ve site No. in which case they are lake";
luntopFile << " levels in metres\n";
luntopFile << "The next two rows give the number of flows to be used for fitting";
luntopFile << " and printing, followed by the reach numbers to which they relate\n";
}
}
isnow_method = 2;
if (isnow_method != 2) {
// get any snow input parameters
ddf = -1.0;
ifstream snowinpFile("snowinp.txt");
if (!snowinpFile.is_open()) {
cerr << "Failed to open 'snowinp.txt'\n";
} else {
snowinpFile >> ddf;
snowinpFile.close();
}
} else {
cout << "Snow input parameters are not used in this version of TopNet (isnow_method == 2)\n";
}
ofstream snowoutFile("snowout.txt"); // unit 78
// ***********************************************************************
// INITIALISE FOR REACH ROUTING
// The routing cannot handle 0 slopes so fix those.
// The relationship used is
// S = C A^theta = .133351 (A[m^2])^(-.425) = 47.31513 (A[mm^2])^(-.425)
// This was fit as a "lower bound" to the scatter in a slope vs area
// plot for the Grey river, New Zealand.
for (jr1 = 1; jr1 <= model1::Nrch; jr1++) {
jr = ll[jr1+ngut-1];
// smin=47.31513 * area(jr)**(-0.425) Area not valid due to reachlogic above inoperable
Rp[0][jr-1] = max(smin, Rp[0][jr-1]);
}
for (i = 0; i < maxSlp; i++) {
snowst[i] = 0.0;
}
#ifdef ZBAR_IN
double t0 = (double)clock()/(double)CLOCKS_PER_SEC;
char cr;
ifstream zbarInFile("zbar_in.dat");
if (!zbarInFile) {
cerr << "Failed to open zbar_in.dat\n";
exit(1);
}
int in_istep, in_basin;
for (istep = 0; istep <= m; istep++) {
for (jsub = 0; jsub < Nsub; jsub++) {
zbarInFile >> in_istep >> in_basin; // discard the first two integers.
for (n = 0; n < nreg; n++) {
zbarInFile >> zbar_in[istep][jsub][n];
}
zbarInFile.get(cr); // read line end
}
}
zbarInFile.close();
double t1 = static_cast<double>(clock())/static_cast<double>(CLOCKS_PER_SEC);
cout << t1 - t0 << " seconds to read depth to water file \n";
#endif
if (nooksack) {
n_irrig = 1; // DGT warning. If this is changed the logic associated with
// averaging depth_irrig_dem will need to be adjusted
n_drainage = 2;
}
else {
n_irrig = 0;
n_drainage = 0;
}
// start of time loop +++++++++++++++++++++++++++++++++++++++++
cout << "Starting time loop\n";
// cout << "Period starting date/time for calculation of mean is: " << asctime(startAnnual);
for (istep = 0; istep <= m; istep++) {
#if TRACE
static int nsteps = 0;
if (nsteps < MAX_TRACE) {
traceFile << setw(10) << "timestep " << istep << '\n';
}
nsteps++;
#endif
if (timeinfo->tm_mday == 1) {
strftime(dateStr, 11,"%Y %m %d", timeinfo);
cout << "--------------------------- " << dec << setw(5) << istep << " --- " << dateStr << "-----------------------------\n";
}
// LOOP OVER SUBBASINS
// Work through basins in "stream order" order, i.e.,
// do basins that feed a first order channels before second order and so on.
// The order is already defined by array LL, and linkR(2,:) gives the basin number.
// The first ngut rows in linkR refer to basins (basin # are transformed by
// having ngut added to them). The next nrch rows of linkR tell which reach
// is fed by which basin and/or upstream reach. We want the basin numbers
// from column 2, rows ngut+1 to ngut+nrch in linkR and we want them in the order they
// should be processed in. The order is in LL(ngut+1) to LL(ngut+nrch). DGT
// introduced an additive factor of MAXCHN to avoid confusion between basin
// and reach numbers and we have to remove it.
#ifdef SS
if (timeinfo->tm_mon == 0 && timeinfo->tm_mday == 1) { // January 1st
cout << " Updating depth to water table. timestep " << setw(5) << istep << " " << asctime(timeinfo);
char cr;
ifstream zbarInFile("zbar_ss.dat");
if (!zbarInFile) {
cerr << "Failed to open zbar_ss.dat\n";
exit(1);
}
double MF_depth; // MODFLOW output
int sub;
for (jsub = 0; jsub < Nsub; jsub++) {
zbarInFile >> sub >> MF_depth;
//cerr << dec << setw(3) << sub;
for (n = 0; n < nreg; n++) {
zbar[jsub][n] = MF_depth;
//cerr << fixed << setw(12) << setprecision(6) << zbar[jsub][n];
}
zbarInFile.get(cr); // read line end
//cerr << '\n';
}
zbarInFile.close();
}
#endif // SS
#ifdef IRR
if (timeinfo->tm_mon == 8 && timeinfo->tm_mday == 30) { // September 30th
cout << " Updating depth to water table. timestep " << setw(5) << istep << " " << asctime(timeinfo);
char cr;
ifstream zbarInFile("zbar_irr.dat");
if (!zbarInFile) {
cerr << "Failed to open zbar_irr.dat\n";
exit(1);
}
double MF_depth; // MODFLOW output
int sub;
for (jsub = 0; jsub < Nsub; jsub++) {
zbarInFile >> sub >> MF_depth;
//cerr << dec << setw(3) << sub;
for (n = 0; n < nreg; n++) {
zbar[jsub][n] = MF_depth;
//cerr << fixed << setw(12) << setprecision(6) << zbar[jsub][n];
}
zbarInFile.get(cr); // read line end
//cerr << '\n';
}
zbarInFile.close();
}
#endif // IRR
#ifdef NIRR
if (timeinfo->tm_mon == 2 && timeinfo->tm_mday == 31) { // March 31st
cout << " Updating depth to water table. timestep " << setw(5) << istep << " " << asctime(timeinfo);
char cr;
ifstream zbarInFile("zbar_nirr.dat");
if (!zbarInFile) {
cerr << "Failed to open zbar_nirr.dat\n";
exit(1);
}
double MF_depth; // MODFLOW output
int sub;
for (jsub = 0; jsub < Nsub; jsub++) {
zbarInFile >> sub >> MF_depth;
//cerr << dec << setw(3) << sub;
for (n = 0; n < nreg; n++) {
zbar[jsub][n] = MF_depth;
//cerr << fixed << setw(12) << setprecision(6) << zbar[jsub][n];
}
zbarInFile.get(cr); // read line end
//cerr << '\n';
}
zbarInFile.close();
}
#endif // NIRR
for (jsub = 1; jsub <= Nsub; jsub++) {
js = jsub-1; // Fortran to C indexing
prec = bRain[js][max(istep, 1)-1];
jr = ll[js];
// snow
if (istep > 0) {
// Subroutine to compute ET DGT 17 May 1998
// parameters for ET - because at least albedo may be fitted, it is not possible
// to put this loop outside calcts - therefore leave as is for the time being RPI 18/3/2004
// debugging. To debug a particular watershed uncomment the lines below and enter its topnetID
// if (jsub == 87 .and. istep == 245)then
// js=jsub
// snowcontrol(1)=1
// else
// snowcontrol(1)=0
// endif
albedo = Sp[11][js];
rlapse = Sp[12][js];
elevsb = Sp[13][js];
temper_temp = temper[istep-1];
dewp_temp = dewp[istep-1];
trange_temp = tRange[istep-1];
if (ievap_method != 0)
dewp_temp = bTdew[js][istep-1]; // for use with penman-Monteith
if (ievap_method != 0)
trange_temp = min(30.0, max(0.0, bTmax[js][istep-1] - bTmin[js][istep-1])); // check for bad data
if (ievap_method != 0)
temper_temp = (bTmax[js][istep-1] + bTmin[js][istep-1])/2;
for (i = 0; i < 12; i++) {
for (j = 0; j < maxSlp; j++) {
bdtBarR[j][i] = bdtBar[i][j];
}
}
sHour = 0;
// Warning here - this code would have to be changed for time steps other than daily
// In general it would be better to inherit sHour from the calling program but coming in at 240000
// the integration of radiation across a day fails. A more general solution to this issue would involve
// reprogramming hyri to handle time steps that cross the day break.
// Setting sHour=0 also achieves compatibility with the convention that inputs are associated with
// measurements recorded at the end of the time step (240000 for daily) but ET and snow computations
// integrate from the start time over the time step.
elev = elevsb; // we are driving this using a basin average temperature, so it must be at basin average elev
etall(bXlat[js], bXlon[js], stdlon, elev, bdtBarR[js], PETsngl, temper_temp, dewp_temp,
trange_temp, elevsb, albedo, rlapse, sDate, sHour, interval, m, istep, iyear, month, iday, ihr,
imm, isec, hour1, bTmin[js][istep-1], bTmax[js][istep-1], wind2m[istep-1], ievap_method);
pet = PETsngl;
ta= (bTmin[js][istep-1] + bTmax[js][istep-1])/2.0;
p = bRain[js][istep-1]*3600.0/interval/1000.0; //mm/int *3600s/h / (s/int) / (1000mm/m) = m/hr for snowueb
ws = wind2m[istep-1];
// RH=1 !unknown for Nooksack, hope we don't need it!
qsiobs = 0; // unknown for Nooksack, hope we don't need it!
qnetob = 0; // unknown for Nooksack, hope we don't need it!
if (isnow_method == 2) {
snowcontrol[4] = iyear*10000 + month*100 + iday;
ihh = hour1;
imm = (hour1-ihh)*60;
iss = hour1*3600-ihh*3600-imm*60;
snowcontrol[5] = ihh*10000 + imm*100*iss;
snowforcing[0] = ta;
snowforcing[1] = p;
snowforcing[2] = ws;
snowforcing[3] = dewp_temp;
// snowforcing(4) = 237.3/(1/(log(RH)/17.27+Ta/(Ta+237.3))-1) ! from http://williams.best.vwh.net/ftp/avsig/avform.txt
snowforcing[4] = bTmax[js][istep-1] - bTmin[js][istep-1];
snowforcing[5] = qsiobs;
snowforcing[6] = qnetob;
snowsitev[0] = 0; // Sp(39,js) !0. !forest cover DGT 4/1/05
snowsitev[1] = Sp[13][js]; //elevsb
snowsitev[2] = bXlat[js]; //lat
snowsitev[3] = bXlon[js]; //lon
snowsitev[4] = stdlon; //stdlong
snowsitev[5] = Sp[13][js]; //elevtg is assumed = elevsb
snowsitev[6] = Sp[12][js]; //rlapse
snowsitev[7] = 0.0; //slope
snowsitev[8] = 0.0; //azimuth
snowueb(istep,jsub,snowsitev, snowstatev[js], snowparam, ndepletionpoints, dfc, snowcontrol, bdtBarR[js],
snowforcing, snowsurfacetemp[js], snowaveragetemp[js], timestep, nstepday,
surfacewaterinput, snowevaporation, // outputs (both in m/h)
areafractionsnow, js+1); // added js so that within snow one knows which element one is in for debugging
cump[js] = cump[js] + p*timestep;
cume[js] = cume[js] + snowevaporation*timestep;
cummr[js] = cummr[js] + surfacewaterinput*timestep;
errmbal[js] = w1[js] + cump[js] - cummr[js] - cume[js] - snowstatev[js][1];
if (snowcontrol[1] >= 2) {
snowcontrol3_File << dec << setw(2) << js;
snowcontrol3_File << dec << setw(5) << istep;
snowcontrol3_File << dec << setw(5) << iyear;
snowcontrol3_File << dec << setw(3) << month;
snowcontrol3_File << dec << setw(3) << iday;
snowcontrol3_File << fixed << setw(5) << setprecision(1) << hour1;
snowcontrol3_File << fixed << setw(18) << setprecision(9) << surfacewaterinput;
snowcontrol3_File << fixed << setw(18) << setprecision(9) << snowevaporation;
snowcontrol3_File << fixed << setw(18) << setprecision(9) << areafractionsnow;
snowcontrol3_File << fixed << setw(18) << setprecision(9) << snowstatev[js][1];
snowcontrol3_File << fixed << setw(18) << setprecision(9) << cump[js];
snowcontrol3_File << fixed << setw(18) << setprecision(9) << cume[js];
snowcontrol3_File << fixed << setw(18) << setprecision(9) << cummr[js] << '\n'; //swe
}
} else {
if (ddf >= 0) {
snow(snowoutFile, temper, elevtg[0], elevsb, rlapse, prec, ddf, snowst[js], interval, Nsub, m, js+1, istep, maxSlp, maxInt);
}
}
if ( modwrt && istep == m ) {
if (idebugoutput >= 1) { // DGT 8/17/05 debugout
lunpFile << "\n\n Output for sub-catchment " << dec << setw(3) << js+1 << '\n';
lunpFile << " ---------------------------\n";
}
}
}
prec = surfacewaterinput*1000.0*interval/3600.0; //mm/timestep=m/h*1000mm/m*h/timestep
wt0 = 1.0;
kk = 0;
depth_irrig_dem = 0.0;
// xIRR(:)=0
baseflow[js] = 0.0;
totalrunoff[js] = 0.0;
evap_for_watermgmt[js] = 0.0;
potentialevap[js] = 0.0; // DGT 10/21/12 initialize for averaging over categories
surfro[js] = 0.0;
canstore[js] = 0.0;
soilstore[js] = 0.0;
tiled[js] = 0.0;
ditchd[js] = 0.0;
zbm_d = 0;
acsem_d = 0;
s1_d = 0;
zbar_d = 0.0;
sr_d = 0;
cv_d = 0;
bal_d = 0;
sumr_d = 0;
sumq_d = 0;
sumae_d = 0;
s0_d = 0;
sumpe_d = 0;
qinst_out_d = 0;
dr_out_d = 0; // Initializing
art_drainage_out = 0.0; // Artificial drainage
#ifdef ZBAR_OUT
if (istep > 0) {
zbarFile << dec << setw(4) << istep;
zbarFile << dec << setw(4) << jsub;
// all nreg regions have the same value at this point in the program
dateZbarFile << dec << setw(4) << istep;
dateZbarFile << dec << setw(4) << jsub;
dateZbarFile << dec << setw(11) << dateStr;
}
#endif
// the previous timestep's call to watermgmt calculated groundwater_to_take
if (istep > 1) { //can't do any pumping on first timestep because we haven't called watermgmt yet: assume zero pumping when istep=1
dth1 = Sp[3][js]; // dgt 11/4/07 added this line to get dth1 from corresponding basin parameter
for (n = 0; n < nreg; n++) {
#ifdef ZBAR_IN
zbar[js][n] = zbar_in[istep][js][n];
#endif // ZBAR_IN
zbar[js][n] += groundwater_to_take[js]/(Sp[0][js]/1.0e6)/dth1; //RAW 18-Jul-2005 bug fix: Sp(1,JS)/1e6 was Sp(1,JS)*1e6
}
// m/timestep = m3/timestep/(m^2)/porosity
// DGT 11/4/07 changed the above from DTH to DTH1 because topmodel saturated zone works with that.
}
for (i_irrig = 1; i_irrig <= (n_irrig+1); i_irrig++) {
if (i_irrig == 1) {
wt1 = wt0*(1.0 - Sp[19][js]); //unirrigated fraction
rate_irrig = 0;
} else {
wt1 = wt0*Sp[19][js]; // irrigated fraction