-
Notifications
You must be signed in to change notification settings - Fork 16
/
EcoSLIM.f90
1910 lines (1647 loc) · 72.9 KB
/
EcoSLIM.f90
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
!--------------------------------------------------------------------
! **EcoSLIM** is a Lagrangian, particle-tracking that simulates advective
! and diffusive movement of water parcels. This code can be used to
! simulate age, diagnosing travel times, source water composition and
! flowpaths. It integrates seamlessly with **ParFlow-CLM**.
!
! Developed by: Reed Maxwell-August 2016 ([email protected])
!
! Contributors: Laura Condon ([email protected])
! Mohammad Danesh-Yazdi ([email protected])
! Lindsay Bearup ([email protected])
! Zach Perzan ([email protected])
!
! released under GNU LPGL, see LICENSE file for details
!
!--------------------------------------------------------------------
! MAIN FORTRAN CODE
!--------------------------------------------------------------------
! EcoSLIM_main.f90: The main fortran code performing particle tracking
! Use Makefile provided to build.
!
!
!
!--------------------------------------------------------------------
! INPUTS
!--------------------------------------------------------------------
! slimin.txt: Includes the domain's geometric information,
! ParFlow timesteps and their total number, and particles
! initial locations.
!
!--------------------------------------------------------------------
! SUBROUTINES
!--------------------------------------------------------------------
! pfb_read(arg1,...,arg5).f90: Reads a ParFlow .pfb output file and
! stores it in a matrix. Arguments
! in order are:
!
! - arg1: Name of the matrix in which
! ParFlow .pfb is stored
! - arg2: Corresponding .pfb file name,
! e.g., test.out.press.00100.pfb
! - arg3: Number of cells in x-direction
! - arg4: Number of cells in y-direction
! - arg5: Number of cells in z-direction
!
!--------------------------------------------------------------------
! OUTPUTS
!--------------------------------------------------------------------
! XXXX_log.txt: Reports the domain's geometric information,
! ParFlow's timesteps and their total number,
! and particles initial condition. XXXX is the name of
! the SLIM2 run already set in slimin.txt
!
! XXXX_particle.3D: Contains particles' trajectory information in
! space (i.e., X, Y, Z) and time (i.e., residence
! time). XXXX is the name of the SLIM2 run
! already set in slimin.txt
!
! XXXX_endparticle.txt: Contains the final X, Y, Z location of all
! particles as well as their travel time.
! XXXX is the name of the SLIM2 run already set
! in slimin.txt
!
!--------------------------------------------------------------------
! CODE STRUCTURE
!--------------------------------------------------------------------
! (1) Define variables
!
! (2) Read inputs, set up domain, write the log file, and
! initialize particles,
!
! (3) For each timestep, loop over all particles to find and
! update their new locations
!--------------------------------------------------------------------
program EcoSLIM
use ran_mod
implicit none
!--------------------------------------------------------------------
! (1) Define variables
!--------------------------------------------------------------------
real*8,allocatable::P(:,:)
! P = Particle array [np,attributes]
! np = Number of particles
! P(np,1) = X coordinate [L]
! P(np,2) = Y coordinate [L]
! P(np,3) = Z coordinate [L]
! P(np,4) = Particle residence time [T]
! P(np,5) = Saturated particle residence time [T]
! P(np,6) = Particle mass; assigned via preciptiation or snowmelt rate (Evap_Trans*density*volume*dT)
! P(np,7) = Particle source (1=IC, 2=rain, 3=snowmelt...)
! P(np,8) = Particle Status (1=active, 0=inactive)
! P(np,9) = concentration
! P(np,10) = Exit status (1=surface outflow, 2=ET, 3=subsurface outflow)
! P(np,11) = Particle Number (This is a unique integer identifier for the particle)
! P(np,12) = Partical Initial X coordinate [L]
! P(np,13) = Partical Initial Y coordinate [L]
! P(np,14) = Partical Initial Z coordinate [L]
! P(np,15) = Time that particle was added [T]
! P(np,16) = Length of flow path [L]
! P(np,17) = Length of saturated flow path [L]
! P(np,18:(17+nind)) = Length of flow path in indicator i [L]
! P(np,(18+nind):(17+nind*2)) = particle age in indicator i [T]
!@ RMM, why is this needed?
real*8,allocatable::PInLoc(:,:)
! PInLoc(np,1) = Particle initial X location
! PInLoc(np,2) = Particle initial Y location
! PInLoc(np,3) = Particle initial Z location
real*8,allocatable::Vx(:,:,:)
real*8,allocatable::Vy(:,:,:)
real*8,allocatable::Vz(:,:,:)
! Vx = Velocity x-direction [nx+1,ny,nz] -- ParFlow output
! Vy = Velocity y-direction [nx,ny+1,nz] -- ParFlow output
! Vz = Velocity z-direction [nx,ny,nz+1] -- ParFlow output
real*8,allocatable::C(:,:,:,:)
! Concentration array, in i,j,k with l (first index) as consituent or
! property. These are set by user at runtime using input
!real*8,allocatable::ET_grid(:,:,:,:)
! ET array, in i,j,k with l (first index) as consituent or property
CHARACTER*20,allocatable:: conc_header(:)
! name for variables written in the C array above. Dimensioned as l above.
real*8,allocatable::Time_Next(:)
! Vector of real times at which ParFlow dumps outputs
real*8,allocatable::dz(:), Zt(:)
! Delta Z values in the vertical direction
! Elevations in z-direction in local coordinates
real*8,allocatable::Sx(:,:) ! Sx: Slopes in x-direction (not used)
real*8,allocatable::Sy(:,:) ! Sy: Slopes in y-direction (not used)
real*8,allocatable::Saturation(:,:,:) ! Saturation (read from ParFlow)
real*8,allocatable::Porosity(:,:,:) ! Porosity (read from ParFlow)
real*8,allocatable::EvapTrans(:,:,:) ! CLM EvapTrans (read from ParFlow, [1/T] units)
real*8,allocatable::CLMvars(:,:,:) ! CLM Output (read from ParFlow, following single file
! CLM output as specified in the manual)
real*8,allocatable::Pnts(:,:), DEM(:,:) ! DEM and grid points for concentration output
real*8,allocatable::Ind(:,:,:) ! Flowpath indicator file
real*8,allocatable::BoundFlux(:,:,:) ! An array of flux across each domain boundary
real*8,allocatable::PSourceBool(:,:,:) ! Boolean particle source file
integer Ploc(3)
! Particle's location whithin a cell
integer nind, itemp
! Number of units in the indicator file, temporary integer to store current indicator
integer nx, nnx, ny, nny, nz, nnz, nztemp
! number of cells in the domain and cells+1 in x, y, and z directions
integer np_ic, np, np_active, np_active2, icwrite, jj, npnts, ncell, npout
! number of particles for intial pulse IC, total, and running active
integer nt, n_constituents
! number of timesteps ParFlow; numer of C vectors written for VTK output
integer pid
! Counter for particle ID number
real*8 pfdt, advdt(3), Ltemp
! ParFlow timestep value, advection timestep for each direction
! for each individual particle step; used to chose optimal particle timestep
! Temporary storage for step length calculations
integer pft1, pft2, tout1, pfnt, n_cycle
! parflow start and stop file numbers number of ParFlow timesteps
! flag specifying the number for the first output write (0= start with pft1)
! number of timestep cycles
real*8 Time_first, part_tstart
! initial timestep for Parflow ((pft1-1)*pfdt)
! initial time for particles (used for recording the particle insert time)
integer kk, ktemp
! Loop counter for the time steps (pfnt)
! temporary loop counter
integer pfkk, outkk
! Counter for the file numbers starts at pft1
! Counter for the output writing
integer ii
! Loop counter for the number of particles (np)
integer iflux_p_res, i_added_particles
! Number of particles per cell for flux input and
! number of particle added total for precip input per
! timestep
integer i, j, k, l, ik, ji, m, ij, nzclm, nCLMsoil
integer itime_loc
! Local indices / counters
integer*4 ir
character*200 runname, filenum, filenumout, pname, fname, vtk_file, DEMname, Indname, Boolname
! runname = SLIM runname
! filenum = ParFlow file number
! filenumout = File number for Ecoslim writing
! pname = ParFlow output runname
! fname = Full name of a ParFlow's output
! vtk_file = concentration file
! DEMname = DEM file name
! Indname = Indicator File name
! Boolname = Boolean particle source file name
integer nlines
! number of lines in slimin.txt
real*8 Clocx, Clocy, Clocz, Z, maxz
! The fractional location of each particle within it's grid cell
! Particle Z location
real*8 V_mult
! Multiplier for forward/backward particle tracking
! If V_mult = 1, forward tracking
! If V_mult = -1, backward tracking
logical clmtrans, clmfile, velfile, boolfile
! logical for mode of operation with CLM, will add particles with P-ET > 0
! will remove particles if ET > 0
! clmfile governs reading of the full CLM output, not just evaptrans
! velfile governs addition of particles to the domain via velocity pfb output
! boolfile governs reading of a PSourceBool file, which spatially limits particle addition
real*8 dtfrac
! fraction of dx/Vx (as well as dy/Vy and dz/Vz) assuring
! numerical stability while advecting a particle to a new
! location.
real*8 Xmin, Xmax, Ymin, Ymax, Zmin, Zmax
! Domain boundaries in local / grid coordinates. min values set to zero,
! DEM is read in later to output to Terrain Following Grid used by ParFlow.
real*8 dx, dy
! Domain's number of cells in x and y directions
real*8 Vpx, Vpy, Vpz
! Particle velocity in x, y, and z directions
real*8 particledt, delta_time
! The time it takes for a particle to displace from
! one location to another and the local particle from-to time
! for each PF timestep
real*8 mean_age, mean_comp, mean_mass, total_mass
! mean age and composition and mass of all particles in domain
real*8 local_flux, et_flux, water_vol, Zr, z1, z2, z3
! The local cell flux convergence
! The volumetric ET flux
! The availble water volume in a cell
! random variable
real*8 Xlow, Xhi, Ylow, Yhi, Zlow, Zhi
! Particles initial locations i.e., where they are injected
! into the domain.
integer, dimension(2):: xbd, ybd, zbd, vxbd, vybd, vzbd
! used to calculate flux at corners of domain
! xbd = (/ 1, nx /)
! vxbd = (/ 1, nx+1 /)
real*8, dimension(2):: fluxsign
! sign (+/-) of velocity at the each end of the domain that
! corresponds to flux into the domain
! fluxsign = (/ 1, -1 \)
real*8 velx, vely, velz
! float of an individual value of Vx, Vy, Vz
integer xi, yi, zi
! counters within xbd, etc and xfluxsign, etc
! density of water (M/L3), molecular diffusion (L2/T), fractionation
real*8 denh2o, moldiff, Efract !, ran1
! time history of ET, time (1,:) and mass for rain (2,:), snow (3,:),
! PET balance is water balance flux from PF accumulated over the domain at each
! timestep
real*8, allocatable::ET_age(:,:), ET_comp(:,:), water_balance(:,:), ET_mass(:)
real*8, allocatable::Surf_age(:,:), Surf_mass(:,:), Surf_comp(:,:)
real*8, allocatable::Subsurf_age(:,:), Subsurf_mass(:,:), Subsurf_comp(:,:)
integer, allocatable:: ET_np(:), Surf_np(:), Subsurf_np(:)
real*8 ET_dt, DR_Temp
! time interval for ET
! integer counters and operators.
! the first set are used for total run timing the latter for component timing
integer Total_time1, Total_time2, t1, t2, IO_time_read, IO_time_write, parallel_time
integer sort_time
!! integers for writing C or point based output
integer ipwrite, ibinpntswrite
! integers for writing gridded ET outputs
integer etwrite
!! IO control
!! ipwrite controls an ASCII, .3D particle file not recommended due to poor performance
!! this is left as a compiler option, currently disabled
!!!!!!!
!! icwrite controls VTK, binary grid based output where particle masses, concentrations,
!! ages are mapped to a grid and written every N timesteps. This is the most effiecient output
!! but loses some accuracy or flexbility becuase individual particle locations are aggregated
!! to the grid
!!!!!!!!
!! ibinpntswrite controls VTK, binary output of particle locations and attributes. This is much faster
!! than the .3D ASCII output but is still much slower than grid based output. It provides the most
!! information as particle locations are preserved
interface
SUBROUTINE vtk_write(time,x,conc_header,ixlim,iylim,izlim,icycle,n_constituents,Pnts,vtk_file)
real*8 :: time
REAL*8 :: x(:,:,:,:)
CHARACTER (LEN=20) :: conc_header(:)
INTEGER*4 :: ixlim
INTEGER*4 :: iylim
INTEGER*4 :: izlim
REAL*8 :: dx
REAL*8 :: dy
REAL*8 :: dz(izlim)
REAL*8 :: Pnts(:,:)
INTEGER :: icycle
INTEGER :: n_constituents
CHARACTER (LEN=200) :: vtk_file
end subroutine vtk_write
SUBROUTINE vtk_write_points(P,np_active, np,icycle,vtk_file,dx,dy,nx,ny,maxZ,DEM)
REAL*8 :: P(:,:)
INTEGER :: icycle
INTEGER*4 :: np_active
INTEGER*4 :: np
INTEGER :: n_constituents
REAL*8 :: dx
REAL*8 :: dy
INTEGER*4 :: nx
INTEGER*4 :: ny
REAL*8 :: maxZ
REAL*8 :: DEM(:,:)
CHARACTER (LEN=200) :: vtk_file
end subroutine vtk_write_points
end interface
!Set up timing
Total_time1 = 0
Total_time2 = 0
t1 = 0
t2 = 0
IO_time_read = 0
IO_time_write = 0
parallel_time = 0
sort_time = 0
call system_clock(Total_time1)
!--------------------------------------------------------------------
! (2) Read inputs, set up domain, write the log file, and
! initialize particles
!--------------------------------------------------------------------
! Note: The following file numbers refer to
!
! - #10: slimin.txt
! - #11: runname_log.txt
! - #12: runname_particle.3D (visualizes particles in VisIT)
! - #13: runname_endparticle.txt
call system_clock(T1)
! Count number of lines in slimin.txt to determine file version
! if nlines >= 31, then velfile option has been provided (even if false) and
! if nlines >= 32, then boolfile name has been provided (even if '')
nlines = 0
open (10,file='slimin.txt')
do
read (10,*,end=1010)
nlines = nlines + 1
end do
1010 close(10)
! open SLIM input .txt file
open (10,file='slimin.txt')
! read SLIM run name
read(10,*) runname
! read ParFlow run name
read(10,*) pname
! read DEM file name
read(10,*) DEMname
! open/create/write the output log.txt file. If doesn't exist, it's created.
open(11,file=trim(runname)//'_log.txt')
write(11,*) '### EcoSLIM Log File'
write(11,*)
write(11,*) 'run name:',trim(runname)
write(11,*)
write(11,*) 'ParFlow run name:',trim(pname)
write(11,*)
if (DEMname /= '') then
write(11,*) 'ParFlow DEM name:',trim(DEMname)
else
write(11,*) 'Not reading ParFlow DEM'
end if
write(11,*)
! read domain number of cells and number of particles to be injected
read(10,*) nx
read(10,*) ny
read(10,*) nz
! read in number of particles for IC (if np_ic = -1 then restart from a file)
read(10,*) np_ic
! read in the number of particles total
read(10,*) np
! check to make sure we don't assign more particles for IC than we have allocated
! in total
if (np_ic > np) then
write(11,*) ' warning NP_IC greater than NP total'
np = np_ic
end if
! write nx, ny, nz, and np in the log file
write(11,*) ' Grid information'
write(11,*) 'nx:',nx
write(11,*) 'ny:',ny
write(11,*) 'nz:',nz
write(11,*)
write(11,*) ' Particle IC Information'
write(11,*) 'np IC:',np_ic
if (np_ic == -1) &
write(11,*) 'Reading particle restart file:',trim(runname)//'_particle_restart.bin'
write(11,*) 'np:',np
! grid +1 variables
nnx=nx+1
nny=ny+1
nnz=nz+1
itemp=0
nCLMsoil = 10 ! number of CLM soil layers over the root zone
nzclm = 13+nCLMsoil ! CLM output is 13+nCLMsoil layers for different variables not domain NZ,
! e.g. 23 for 10 soil layers (default) and 17 for 4 soil layers (Noah soil
! layer setup)
! number of things written to C array, hard wired at 2 now for Testing
!n_constituents = 5
n_constituents = 9
!allocate arrays
allocate(PInLoc(np,3))
allocate(Sx(nx,ny),Sy(nx,ny), DEM(nx,ny))
allocate(dz(nz), Zt(0:nz))
allocate(Vx(nnx,ny,nz), Vy(nx,nny,nz), Vz(nx,ny,nnz))
allocate(Saturation(nx,ny,nz), Porosity(nx,ny,nz),EvapTrans(nx,ny,nz), Ind(nx,ny,nz), BoundFlux(nx,ny,nz))
allocate(PSourceBool(nx,ny,nz))
allocate(CLMvars(nx,ny,nzclm))
allocate(C(n_constituents,nx,ny,nz))
!allocate(ET_grid(1,nx,ny,nz))
allocate(conc_header(n_constituents))
!Intialize everything to Zero
Vx = 0.0d0
Vy = 0.0d0
Vz = 0.0d0
Saturation = 0.0D0
Porosity = 0.0D0
EvapTrans = 0.0d0
BoundFlux = 0.0d0
C = 0.0d0
!ET_grid=0.0d0
! read dx, dy as scalars
read(10,*) dx
read(10,*) dy
! read dz as an array
read(10,*) dz(1:nz)
! read in (constant for now) ParFlow dt
read(10,*) pfdt
! read in parflow start and stop times
read(10,*) pft1
read(10,*) pft2
read(10,*) tout1
read(10,*) part_tstart
read(10,*) n_cycle
pfnt=n_cycle*(pft2-pft1+1)
if (tout1 == 0) then
outkk = pft1
tout1 = pft1
else
outkk = tout1
end if
! set ET DT to ParFlow one and allocate ET arrays accordingly
ET_dt = pfdt
allocate(ET_age(pfnt,5), ET_comp(pfnt,3), ET_mass(pfnt), ET_np(pfnt))
allocate(water_balance(pfnt,2))
ET_age = 0.0d0
ET_mass = 0.0d0
ET_comp = 0.0d0
ET_np = 0
water_balance = 0.0D0
allocate(Surf_age(pfnt,5), Surf_comp(pfnt,3), Surf_mass(pfnt,5), Surf_np(pfnt))
Surf_age = 0.0d0
Surf_mass = 0.0d0
Surf_comp = 0.0d0
Surf_np = 0
allocate(Subsurf_age(pfnt,5), Subsurf_comp(pfnt,3), Subsurf_mass(pfnt,5), Subsurf_np(pfnt))
Subsurf_age = 0.0d0
Subsurf_mass = 0.0d0
Subsurf_comp = 0.0d0
Subsurf_np = 0
! clear out output particles
npout = 0
!! IO control, each value is a timestep interval, e.g. 1= every timestep, 2=every other, 0 = no writing
read(10,*) ipwrite ! controls an ASCII, .3D particle file not recommended due to poor performance
read(10,*) ibinpntswrite ! controls VTK, binary output of particle locations and attributes
read(10,*) etwrite ! controls ASCII ET output
read(10,*) icwrite ! controls VTK, binary grid based output where particle masses, concentrations,
! ages are mapped to a grid and written every N timesteps
! allocate and assign timesteps
allocate(Time_Next(pfnt))
do kk = 1, pfnt
Time_Next(kk) = float(kk+outkk-1)*pfdt
end do
Time_first = float(outkk-1)*pfdt
! read in velocity multiplier
read(10,*) V_mult
! do we read in clm evap trans?
read(10,*) clmtrans
! do we read in clm output file?
read(10,*) clmfile
! read in IC number of particles for flux
read(10,*) iflux_p_res
! read in density h2o
read(10,*) denh2o
! read in diffusivity
!moldiff = (1.15e-9)*3600.d0
read(10,*) moldiff
!saving Efract for a later time
!read(10,*) Efract
! fraction of dx/Vx
read(10,*) dtfrac
!wite out log file
write(11,*)
write(11,*) ' Grid Dimensions'
write(11,'("dx:",e12.5)') dx
write(11,'("dy:",e12.5)') dy
write(11,'("dz:",*(e12.5,", "))') dz(1:nz)
write(11,*)
write(11,*) 'Timestepping Information'
write(11,'("ParFlow delta-T, pfdt:",e12.5)') pfdt
write(11,'("ParFlow timesteps, pfnt:",i12)') pfnt
write(11,'("ParFlow start step, pft1:",i12)') pft1
write(11,'("ParFlow end step, pft2:",i12)') pft2
write(11,'("Output step start:",i12)') outkk
write(11,'("Time loops, cycles, n_cycle:",i12)') n_cycle
write(11,'("Total time steps:",i12)') pfnt
write(11,*)
write(11,*) 'V mult: ',V_mult,' for forward/backward particle tracking'
write(11,*) 'CLM Trans: ',clmtrans,' adds / removes particles based on LSM fluxes'
write(11,*)
write(11,*) 'Physical Constants'
write(11,*) 'denh2o: ',denh2o,' M/L^3'
write(11,*) 'Molecular Diffusivity: ',moldiff,' '
!write(11,*) 'Fractionation: ',Efract,' '
write(11,*)
write(11,*) 'Numerical Stability Information'
write(11,'("dtfrac: ",e12.5," fraction of dx/Vx")') dtfrac
read(10,*) nind !read in the number of indicators
read(10,*) Indname
write(11,*)
write(11,*) 'Indicator File'
write(11,*) nind, 'Indicators'
! allocate P, Sx, dz, Vx, Vy, Vz, Saturation, and Porosity arrays
if (nind<=0) then
allocate(P(np,17))
P(1:np,1:6) = 0.0 ! clear out all particle attributes
P(1:np,7:9) = 1.0 ! make all particles active to start with and original from 1 = GW/IC
P(1:np,10) = 0.0 ! no exit
P(1:np,11:14) = 0.0 ! Clear initial ID and location
P(1:np,15) = 0.0 ! Set initial start time to zero
P(1:np,16:17) = 0.0 ! Set initial lengths to zero
write(11,*) 'Zero indicators setting up regular P array'
else
allocate(P(np,(17+nind*2)))
P(1:np,1:6) = 0.0 ! clear out all particle attributes
P(1:np,7:9) = 1.0 ! make all particles active to start with and original from 1 = GW/IC
P(1:np,10) = 0.0 ! no exit
P(1:np,11:14) = 0.0 ! Clear initial ID and location
P(1:np,15) = 0.0 ! Set initial start time to zero
P(1:np,16:17) = 0.0 ! Set initial lengths to zero
P(1:np, 18:(17+nind*2))=0.0 !initizling indicator file times and lengths to zero
write(11,*) 'Adding columns for',nind, 'indicators'
!Read in the indictor file
end if
! Read in velfile and boolfile, if that line has been provided
if (nlines >= 31) then
read(10,*) velfile !bool of whether or not to add particles using velocity fluxes
write(11,*)
write(11,*) 'velfile: ',velfile,' whether or not to add particles from velocity fluxes'
if (nlines >= 32) then
read(10,*) Boolname
write(11,*)
if (Boolname /= '') then
boolfile = .TRUE.
write(11,*) 'Using boolean source files: ', Boolname
else
boolfile = .FALSE.
write(11,*) 'Not using boolean source file'
end if
end if
! If neither of these lines are provided, turn velfile and boolfile off
else
velfile = .FALSE.
boolfile = .FALSE.
write(11,*)
write(11,*) 'Did not find velfile line in slimin.txt. Setting velfile = False'
write(11,*)
write(11,*) 'Did not find PSourceBool line in slimin.txt. Setting boolfile = False'
end if
! end of SLIM input
close(10)
call system_clock(T2)
IO_time_read = IO_time_read + (T2-T1)
!! set up domain boundaries
Xmin = 0.0d0
Ymin = 0.0d0
Zmin = 0.0d0
Xmax = float(nx)*dx
Ymax = float(ny)*dy
Zmax = 0.0d0
do k = 1, nz
Zmax = Zmax + dz(k)
end do
DEM = 0.0d0
! read in DEM
if (DEMname /= '') then
fname = trim(adjustl(DEMname))
call pfb_read(DEM,fname,nx,ny,nztemp)
end if ! DEM
Ind = 1.0d0
! read in Indicator file
if (nind>0) then
if (Indname /= '') then
fname = trim(adjustl(Indname))
call pfb_read(Ind,fname,nx,ny,nz)
write(11,*) 'Read Indicator File:', fname
!write(11,*) 'IndREAD:', nx, ny, nztemp
else
write(11,*) 'WARNING: indicator flag >0 but no indicator file provided'
end if ! Ind
end if
!! hard wire DEM @RMM, to do, need to make this input
!do i = 1, nx
! do j = 1, ny
!DEM(i,j) = 0.0D0 + float(i)*dx*0.05
!end do
!end do
write(11,*)
write(11,*) '## Domain Info'
write(11,'("Xmin:",e12.5," Xmax:",e12.5)') Xmin, Xmax
write(11,'("Ymin:",e12.5," Ymax:",e12.5)') Ymin, Ymax
write(11,'("Zmin:",e12.5," Zmax:",e12.5)') Zmin, Zmax
!! DEM set to zero but will be read in as input
!! Set up grid locations for file output
npnts=nnx*nny*nnz
ncell=nx*ny*nz
allocate(Pnts(npnts,3))
Pnts=0
m=1
! Need the maximum height of the model and elevation locations
Z = 0.0d0
Zt(0) = 0.0D0
do ik = 1, nz
Z = Z + dz(ik)
Zt(ik) = Z
!print*, Z, dz(ik), Zt(ik), ik
end do
maxZ=Z
!! candidate loops for OpenMP
do k=1,nnz
do j=1,nny
do i=1,nnx
Pnts(m,1)=DBLE(i-1)*dx
Pnts(m,2)=DBLE(j-1)*dy
! This is a simple way of handling the maximum edges
if (i <= nx) then
ii=i
else
ii=nx
endif
if (j <= ny) then
jj=j
else
jj=ny
endif
! This step translates the DEM
! The specified initial heights in the pfb (z1) are ignored and the
! offset is computed based on the model thickness
Pnts(m,3)=(DEM(ii,jj)-maxZ)+Zt(k-1)
m=m+1
end do
end do
end do
! Read porosity values from ParFlow .pfb file
fname=trim(adjustl(pname))//'.out.porosity.pfb'
call pfb_read(Porosity,fname,nx,ny,nz)
! Read the in initial Saturation from ParFlow
kk = 0
pfkk=pft1-1
!print*, pfkk
write(filenum,'(i5.5)') pfkk
fname=trim(adjustl(pname))//'.out.satur.'//trim(adjustl(filenum))//'.pfb'
call pfb_read(Saturation,fname,nx,ny,nz)
! By default, add particles everywhere
PSourceBool = 1.0d0
! If we are limiting particles according to boolfile, read in Boolname to PSourceBool
if (boolfile) then
fname=trim(adjustl(Boolname))//'.'//trim(adjustl(filenum))//'.pfb'
call pfb_read(PSourceBool,fname,nx,ny,nz)
end if
! Intialize random seed
ir = -3333
!Initialize np_active2
np_active = 0
!Initialize Ltemp
Ltemp = 0.0
!! Define initial particles' locations and mass
!!
if (np_ic == 0) then
np_active = 0
pid = np_active
end if
if (np_ic > 0) then
np_active = 0
pid = 0
PInLoc=0.0d0
!call srand(333)
do i = 1, nx
do j = 1, ny
do k = 1, nz
if (np_active < np) then ! check if we have particles left
if (Saturation(i,j,k) > 0.0) then ! check if we are in the active domain
if (PSourceBool(i,j,k) == 1.0d0) then ! check if we should add particles here
do ij = 1, np_ic
np_active = np_active + 1
pid=pid + 1
ii = np_active
P(ii,11)=float(pid) !Saving a particle ID number
! assign X, Y, Z locations randomly to each cell
! assign X, Y, Z locations randomly to each cell
P(ii,1) = float(i-1)*dx +ran1(ir)*dx
PInLoc(ii,1) = P(ii,1)
P(ii,12)=P(ii,1) ! Saving the initial location
P(ii,2) = float(j-1)*dy +ran1(ir)*dy
PInLoc(ii,2) = P(ii,2)
P(ii,13)=P(ii,2)
P(ii,15) = part_tstart + 0.0 !setting insert time to the start time
Z = 0.0d0
do ik = 1, k
Z = Z + dz(ik)
end do
P(ii,3) = Z -dz(k)*ran1(ir)
PInLoc(ii,3) = P(ii,3)
P(ii,14)=P(ii,3) ! Saving the initial location
! assign mass of particle by the volume of the cells
! and the water contained in that cell
P(ii,6) = dx*dy*dz(k)*(Porosity(i,j,k) &
*Saturation(i,j,k))*denh2o*(1.0d0/float(np_ic))
P(ii,7) = 1.0d0
P(ii,8) = 1.0d0
! set up intial concentrations
C(1,i,j,k) = C(1,i,j,k) + P(ii,8)*P(ii,6) / &
(dx*dy*dz(k)*(Porosity(i,j,k) &
*Saturation(i,j,k)))
C(2,i,j,k) = C(2,i,j,k) + P(ii,8)*P(ii,4)*P(ii,6)
C(4,i,j,k) = C(4,i,j,k) + P(ii,8)*P(ii,7)*P(ii,6)
C(3,i,j,k) = C(3,i,j,k) + P(ii,8)*P(ii,6)
end do ! particles per cell
end if ! end-if for PSourceBool
end if ! active domain
else
write(11,*) ' **Warning IC input but no paricles left'
write(11,*) ' **Exiting code gracefully writing restart '
goto 9090
end if
end do ! i
end do ! j
end do ! k
!! if np_ic = -1 then we read a restart file
else if (np_ic == -1) then
write(11,*) 'Reading particle restart File:',trim(runname)//'_particle_restart.bin'
! read in full particle array as binary restart file, should name change?,
! potential overwrite confusion
open(116,file=trim(runname)//'_particle_restart.bin', FORM='unformatted', &
access='stream')
read(116) np_active
read(116) pid
if (np_active < np) then ! check if we have particles left
!do ii = 1, np_active
if (nind > 0) then
read(116) P(1:np_active,1:(nind*2+17))
else
read(116) P(1:np_active,1:17)
endif
!end do !11
close(116)
!pid = np_active
write(11,*) 'RESTART np_active :',np_active
write(11,*) 'RESTART pid:',pid
else
write(11,*) ' **Warning restart IC input but no paricles left'
write(11,*) ' **Exiting code *not* (over)writing restart '
close(116)
stop
end if
end if
!! Define initial particles' locations by surface water
!!
if (np_ic < -1) then
write(11,*) 'River IC activated'
np_active = 0
pid = np_active
PInLoc=0.0d0
!call srand(333)
ir = -3333
do i = 1, nx
do j = 1, ny
k = nz
if (np_active < np) then ! check if we have particles left
!if (saturation(i,j,k) >= 0.95d0) then
do ij = 1, abs(np_ic)
np_active = np_active + 1
pid = pid +1
ii = np_active
P(ii,11) = float(pid)
! assign X, Y, Z locations randomly to each cell
! assign X, Y, Z locations randomly to each cell
P(ii,1) = float(i-1)*dx +ran1(ir)*dx
PInLoc(ii,1) = P(ii,1)
P(ii,12)=P(ii,1) ! Saving the initial location
P(ii,2) = float(j-1)*dy +ran1(ir)*dy
PInLoc(ii,2) = P(ii,2)
P(ii,13)=P(ii,2) ! Saving the initial location
P(ii,15) = part_tstart + 0.0 !setting insert time to the start time
Z = 0.0d0
do ik = 1, k
Z = Z + dz(ik)
end do
P(ii,3) = Z !-dz(k)*ran1(ir)
PInLoc(ii,3) = P(ii,3)
P(ii,14)=P(ii,3) ! Saving the initial location
! assign mass of particle by the volume of the cells
! and the water contained in that cell
P(ii,6) = dx*dy*dz(k)*(Porosity(i,j,k) &
*Saturation(i,j,k))*denh2o*(1.0d0/float(np_ic))
P(ii,7) = 1.0d0
P(ii,8) = 1.0d0
! set up intial concentrations
C(1,i,j,k) = C(1,i,j,k) + P(ii,8)*P(ii,6) / &
(dx*dy*dz(k)*(Porosity(i,j,k) &
*Saturation(i,j,k)))
C(2,i,j,k) = C(2,i,j,k) + P(ii,8)*P(ii,4)*P(ii,6)
C(4,i,j,k) = C(4,i,j,k) + P(ii,8)*P(ii,7)*P(ii,6)
C(3,i,j,k) = C(3,i,j,k) + P(ii,8)*P(ii,6)
end do ! particles per cell
!end if !! saturated at top
else
write(11,*) ' **Warning IC input but no paricles left'
write(11,*) ' **Exiting code gracefully writing restart '
goto 9090
end if
end do ! i
end do ! j
end if !! river IC
! Write out intial concentrations
! normalize ages by mass
where (C(3,:,:,:)>0.0) C(2,:,:,:) = C(2,:,:,:) / C(3,:,:,:)
where (C(3,:,:,:)>0.0) C(4,:,:,:) = C(4,:,:,:) / C(3,:,:,:)
!! Set up output options for VTK grid output
!icwrite = 1
vtk_file=trim(runname)//'_cgrid'
conc_header(1) = 'Concentration'
conc_header(2) = 'Age'
conc_header(3) = 'Mass'
conc_header(4) = 'Comp'
conc_header(5) = 'Delta'
conc_header(6) = 'ET_Npart'
conc_header(7) = 'ET_Mass'
conc_header(8) = 'ET_Age'
conc_header(9) = 'ET_Comp'
if(icwrite > 0) &
call vtk_write(Time_first,C,conc_header,nx,ny,nz,0,n_constituents,Pnts,vtk_file)
!! clear out C arrays
C = 0.0D0
!! write timestep header
write(11,*)
write(11,*)
write(11,*) ' **** Transient Simulation Particle Accounting ****'
write(11,*) ' Timestep PFTimestep OutStep Time Mean_Age Mean_Comp Mean_Mass Total_Mass PrecipIn ETOut &
NP_PrecipIn NP_ETOut &
NP_QOut NP_active_old NP_filtered'
flush(11)
!! open exited particle file and write header
!open(114,file=trim(runname)//'_exited_particles.txt')
open(114,file=trim(runname)//'_exited_particles.bin', FORM='unformatted', &
access='stream')
!write(114,*) 'Time X Y Z PTime Mass Comp ExitStatus'
!flush(114)
if(ipwrite < 0) then
! open/create/write the 3D output file which will write particles out each timestemp, very slowly in parallel
open(214,file=trim(runname)//'_total_particle_trace.3D')
write(214,*) 'X Y Z TIME'
end if !! ipwrite < 0