-
Notifications
You must be signed in to change notification settings - Fork 49
/
modstartup.f90
1255 lines (1096 loc) · 47.7 KB
/
modstartup.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
!> \file modstartup.f90
!! Initializes the run
!>
!! Initializes the run.
!>
!! Modstartup reads the namelists and initial data, sets the fields and calls
!! the inits of the other routines where necessary. Reading and writing of the
!! restart files also live in this module.
!! \author Chiel van Heerwaarden, Wageningen U.R.
!! \author Thijs Heus,MPI-M
!! \todo documentation
!! \par Revision list
! This file is part of DALES.
!
! DALES 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.
!
! DALES 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 this program. If not, see <http://www.gnu.org/licenses/>.
!
! Copyright 1993-2009 Delft University of Technology, Wageningen University, Utrecht University, KNMI
!
module modstartup
implicit none
! private
! public :: startup, writerestartfiles,trestart
save
integer (KIND=selected_int_kind(6)) :: irandom= 0 ! * number to seed the randomnizer with
integer :: krand = huge(0), krandumin=1,krandumax=0
real :: randthl= 0.1,randqt=1e-5 ! * thl and qt amplitude of randomnization
real :: randu = 0.5
contains
subroutine startup
!-----------------------------------------------------------------|
! |
! Reads all general options from namoptions |
! |
! Chiel van Heerwaarden 15/06/2007 |
! Thijs Heus 15/06/2007 |
!-----------------------------------------------------------------|
use modglobal, only : initglobal,iexpnr,runtime, dtmax,dtav_glob,timeav_glob,&
lwarmstart,startfile,trestart,&
nsv,itot,jtot,kmax,xsize,ysize,xlat,xlon,xday,xtime,&
lmoist,lcoriol,igrw_damp,geodamptime,lmomsubs,cu, cv,ifnamopt,fname_options,llsadv,&
ibas_prf,lambda_crit,iadv_mom,iadv_tke,iadv_thl,iadv_qt,iadv_sv,courant,peclet,ladaptive,author,lnoclouds,lrigidlid,unudge
use modforces, only : lforce_user
use modsurfdata, only : z0,ustin,wtsurf,wqsurf,wsvsurf,ps,thls,isurf
use modsurface, only : initsurface
use modfields, only : initfields
use modpois, only : initpois
use modradiation, only : initradiation
use modraddata, only : irad,iradiation,&
rad_ls,rad_longw,rad_shortw,rad_smoke,useMcICA,&
timerad,rka,dlwtop,dlwbot,sw0,gc,reff,isvsmoke,lcloudshading
use modtimedep, only : inittimedep,ltimedep
use modboundary, only : initboundary,ksp
use modthermodynamics, only : initthermodynamics,lqlnr, chi_half
use modmicrophysics, only : initmicrophysics
use modsubgrid, only : initsubgrid
use mpi, only : MPI_COMM_WORLD,MPI_INTEGER,MPI_LOGICAL,MPI_CHARACTER
use modmpi, only : initmpi,my_real,myid,nprocx,nprocy,mpierr
implicit none
integer :: ierr
!declare namelists
namelist/RUN/ &
iexpnr,lwarmstart,startfile,runtime,dtmax,dtav_glob,timeav_glob,&
trestart,irandom,randthl,randqt,krand,nsv,courant,peclet,ladaptive,author,&
krandumin, krandumax, randu,&
nprocx,nprocy
namelist/DOMAIN/ &
itot,jtot,kmax,&
xsize,ysize,&
xlat,xlon,xday,xtime,ksp
namelist/PHYSICS/ &
!cstep z0,ustin,wtsurf,wqsurf,wsvsurf,ps,thls,chi_half,lmoist,isurf,lneutraldrag,&
z0,ustin,wtsurf,wqsurf,wsvsurf,ps,thls,lmoist,isurf,chi_half,&
lcoriol,igrw_damp,geodamptime,lmomsubs,ltimedep,irad,timerad,iradiation,rad_ls,rad_longw,rad_shortw,rad_smoke,useMcICA,&
rka,dlwtop,dlwbot,sw0,gc,reff,isvsmoke,lforce_user,lcloudshading,lrigidlid,unudge
namelist/DYNAMICS/ &
llsadv, lqlnr, lambda_crit, cu, cv, ibas_prf, iadv_mom, iadv_tke, iadv_thl, iadv_qt, iadv_sv, lnoclouds
! get myid
call MPI_INIT(mpierr)
call MPI_COMM_RANK( MPI_COMM_WORLD, myid, mpierr )
!read namelists
if(myid==0)then
if (command_argument_count() >=1) then
call get_command_argument(1,fname_options)
end if
write (*,*) fname_options
open(ifnamopt,file=fname_options,status='old',iostat=ierr)
if (ierr /= 0) then
stop 'ERROR:Namoptions does not exist'
end if
read (ifnamopt,RUN,iostat=ierr)
if (ierr > 0) then
print *, 'Problem in namoptions RUN'
print *, 'iostat error: ', ierr
stop 'ERROR: Problem in namoptions RUN'
endif
write(6 ,RUN)
rewind(ifnamopt)
read (ifnamopt,DOMAIN,iostat=ierr)
if (ierr > 0) then
print *, 'Problem in namoptions DOMAIN'
print *, 'iostat error: ', ierr
stop 'ERROR: Problem in namoptions DOMAIN'
endif
write(6 ,DOMAIN)
rewind(ifnamopt)
read (ifnamopt,PHYSICS,iostat=ierr)
if (ierr > 0) then
print *, 'Problem in namoptions PHYSICS'
print *, 'iostat error: ', ierr
stop 'ERROR: Problem in namoptions PHYSICS'
endif
write(6 ,PHYSICS)
rewind(ifnamopt)
read (ifnamopt,DYNAMICS,iostat=ierr)
if (ierr > 0) then
print *, 'Problem in namoptions DYNAMICS'
print *, 'iostat error: ', ierr
stop 'ERROR: Problem in namoptions DYNAMICS'
endif
write(6 ,DYNAMICS)
close(ifnamopt)
end if
!broadcast namelists
call MPI_BCAST(iexpnr ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(lwarmstart ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(startfile ,50,MPI_CHARACTER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(author ,80,MPI_CHARACTER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(runtime ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(trestart ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(dtmax ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(dtav_glob ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(timeav_glob,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(nsv ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(nprocx ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(nprocy ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(itot ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(jtot ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(kmax ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(xsize ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(ysize ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(xlat ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(xlon ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(xday ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(xtime ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(z0 ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(ustin ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
!call MPI_BCAST(lneutraldrag ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(wtsurf ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(wqsurf ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(wsvsurf(1:nsv),nsv,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(ps ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(thls ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(chi_half ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(lmoist ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(lcoriol ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(igrw_damp ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(geodamptime,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(lforce_user,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(lmomsubs ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(ltimedep ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(lrigidlid ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(unudge ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(irad ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(timerad ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(iradiation ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(rad_ls ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(rad_longw ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(rad_shortw ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(rad_smoke ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(useMcIca ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(rka ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(dlwtop ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(dlwbot ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(sw0 ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(gc ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
! CvH call MPI_BCAST(sfc_albedo ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(reff ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(isvsmoke ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(lcloudshading,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(llsadv ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(lqlnr ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(lambda_crit,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(cu ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(cv ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(ksp ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(irandom ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(krand ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(krandumin ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(krandumax ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(randthl ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(randqt ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(randu ,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(ladaptive ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(courant,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(peclet,1,MY_REAL ,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(isurf ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(ibas_prf,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(iadv_mom,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(iadv_tke,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(iadv_thl,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(iadv_qt ,1,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(iadv_sv(1:nsv) ,nsv,MPI_INTEGER,0,MPI_COMM_WORLD,mpierr)
call MPI_BCAST(lnoclouds ,1,MPI_LOGICAL,0,MPI_COMM_WORLD,mpierr)
! Initialize MPI
call initmpi
! Allocate and initialize core modules
call initglobal
call initfields
call initboundary
call initthermodynamics
call initradiation
call initsurface
call initsubgrid
call initpois
call initmicrophysics
call readinitfiles ! moved to obtain the correct btime for the timedependent forcings in case of a warmstart
call inittimedep !depends on modglobal,modfields, modmpi, modsurf, modradiation
call checkinitvalues
end subroutine startup
subroutine checkinitvalues
!-----------------------------------------------------------------|
! |
! Thijs Heus TU Delft 9/2/2006 |
! |
! purpose. |
! -------- |
! |
! checks whether crucial parameters are set correctly |
! |
! interface. |
! ---------- |
! |
! *checkinitvalues* is called from *program*. |
! |
!-----------------------------------------------------------------|
use modsurfdata,only : wtsurf,wqsurf,ustin,thls,isurf,ps,lhetero
use modglobal, only : itot,jtot, ysize,xsize,dtmax,runtime, startfile,lwarmstart,eps1, imax,jmax
use modmpi, only : myid,nprocx,nprocy,mpierr
use modtimedep, only : ltimedep
if(mod(jtot,nprocy) /= 0) then
if(myid==0)then
write(6,*)'STOP ERROR IN NUMBER OF PROCESSORS'
write(6,*)'nprocy must divide jtot!!! '
write(6,*)'nprocy and jtot are: ',nprocy, jtot
end if
call MPI_FINALIZE(mpierr)
stop
else
if(myid==0)then
write(6,*)'jmax = jtot / nprocy = ', jmax
endif
end if
if(mod(itot,nprocx)/=0)then
if(myid==0)then
write(6,*)'STOP ERROR IN NUMBER OF PROCESSORS'
write(6,*)'nprocx must divide itot!!! '
write(6,*)'nprocx and itot are: ',nprocx,itot
end if
call MPI_FINALIZE(mpierr)
stop
else
if(myid==0)then
write(6,*)'imax = itot / nprocx = ', imax
endif
end if
!Check Namroptions
if (runtime < 0)stop 'runtime out of range/not set'
if (dtmax < 0) stop 'dtmax out of range/not set '
if (ps < eps1) stop 'psout of range/not set'
if (thls < eps1) stop 'thls out of range/not set'
if (xsize < 0) stop 'xsize out of range/not set'
if (ysize < 0) stop 'ysize out of range/not set '
if (lwarmstart) then
if (startfile == '') stop 'no restartfile set'
end if
!isurf
if (myid == 0) then
select case (isurf)
case(1)
case(2,10)
case(3:4)
if (wtsurf == -1) stop 'wtsurf not set'
if (wqsurf == -1) stop 'wqsurf not set'
case default
stop 'isurf out of range/not set'
end select
if (isurf ==3) then
if (ustin < 0) stop 'ustin out of range/not set'
end if
end if
if (ltimedep .and. lhetero) then
if (myid == 0) write(6,*)&
'WARNING: You selected to use time dependent (ltimedep) and heterogeneous surface conditions (lhetero) at the same time'
if (myid == 0) write(0,*)&
'WARNING: You selected to use time dependent (ltimedep) and heterogeneous surface conditions (lhetero) at the same time'
endif
end subroutine checkinitvalues
subroutine readinitfiles
use modfields, only : u0,v0,w0,um,vm,wm,thlm,thl0,thl0h,qtm,qt0,qt0h,&
ql0,ql0h,thv0h,sv0,svm,e12m,e120,&
dudxls,dudyls,dvdxls,dvdyls,dthldxls,dthldyls,&
dqtdxls,dqtdyls,dqtdtls,dpdxl,dpdyl,&
wfls,whls,ug,vg,uprof,vprof,thlprof, qtprof,e12prof, svprof,&
v0av,u0av,qt0av,ql0av,thl0av,sv0av,exnf,exnh,presf,presh,rhof,&
thlpcar,thvh,thvf
use modglobal, only : i1,i2,ih,j1,j2,jh,kmax,k1,dtmax,idtmax,dt,rdt,runtime,timeleft,tres,&
rtimee,timee,ntimee,ntrun,btime,dt_lim,nsv,&
zf,dzf,dzh,rv,rd,cp,rlv,pref0,om23_gs,&
ijtot,cu,cv,e12min,dzh,cexpnr,ifinput,lwarmstart,itrestart,&
trestart, ladaptive,llsadv,tnextrestart
use modsubgrid, only : ekm,ekh
use modsurfdata, only : wsvsurf, &
thls,tskin,tskinm,tsoil,tsoilm,phiw,phiwm,Wl,Wlm,thvs,qts,isurf,svs,obl,oblav,&
thvs_patch,lhetero,qskin
use modsurface, only : surface,qtsurf,dthldz
use modboundary, only : boundary
use modmpi, only : slabsum,myid,comm3d,mpierr,my_real
use modthermodynamics, only : thermodynamics,calc_halflev
use moduser, only : initsurf_user
integer i,j,k,n
real, allocatable :: height(:), th0av(:)
real, dimension(2-ih:i1+ih,2-jh:j1+jh,k1) :: thv0
character(80) chmess
allocate (height(k1))
allocate (th0av(k1))
if (.not. lwarmstart) then
!********************************************************************
! 1.0 prepare initial fields from files 'prog.inp' and 'scalar.inp'
! ----------------------------------------------------------------
!--------------------------------------------------------------------
! 1.1 read fields
!-----------------------------------------------------------------
rdt = dtmax / 100.
dt = floor(rdt/tres)
timee = 0
if (myid==0) then
open (ifinput,file='prof.inp.'//cexpnr)
read (ifinput,'(a80)') chmess
write(*, '(a80)') chmess
read (ifinput,'(a80)') chmess
do k=1,kmax
read (ifinput,*) &
height (k), &
thlprof(k), &
qtprof (k), &
uprof (k), &
vprof (k), &
e12prof(k)
end do
close(ifinput)
write(*,*) 'height thl qt u v e12'
do k=kmax,1,-1
write (*,'(f7.1,f8.1,e12.4,3f7.1)') &
height (k), &
thlprof(k), &
qtprof (k), &
uprof (k), &
vprof (k), &
e12prof(k)
end do
if (minval(e12prof(1:kmax)) < e12min) then
write(*,*) 'e12 value is zero (or less) in prof.inp'
do k=1,kmax
e12prof(k) = max(e12prof(k),e12min)
end do
end if
end if ! end if myid==0
! MPI broadcast numbers reading
call MPI_BCAST(thlprof,kmax,MY_REAL ,0,comm3d,mpierr)
call MPI_BCAST(qtprof ,kmax,MY_REAL ,0,comm3d,mpierr)
call MPI_BCAST(uprof ,kmax,MY_REAL ,0,comm3d,mpierr)
call MPI_BCAST(vprof ,kmax,MY_REAL ,0,comm3d,mpierr)
call MPI_BCAST(e12prof,kmax,MY_REAL ,0,comm3d,mpierr)
do k=1,kmax
do j=1,j2
do i=1,i2
thl0(i,j,k) = thlprof(k)
thlm(i,j,k) = thlprof(k)
qt0 (i,j,k) = qtprof (k)
qtm (i,j,k) = qtprof (k)
u0 (i,j,k) = uprof (k) - cu
um (i,j,k) = uprof (k) - cu
v0 (i,j,k) = vprof (k) - cv
vm (i,j,k) = vprof (k) - cv
w0 (i,j,k) = 0.0
wm (i,j,k) = 0.0
e120(i,j,k) = e12prof(k)
e12m(i,j,k) = e12prof(k)
ekm (i,j,k) = 0.0
ekh (i,j,k) = 0.0
end do
end do
end do
!---------------------------------------------------------------
! 1.2 randomnize fields
!---------------------------------------------------------------
krand = min(krand,kmax)
do k = 1,krand
call randomnize(qtm ,k,randqt ,irandom,ih,jh)
call randomnize(qt0 ,k,randqt ,irandom,ih,jh)
call randomnize(thlm,k,randthl,irandom,ih,jh)
call randomnize(thl0,k,randthl,irandom,ih,jh)
end do
do k=krandumin,krandumax
call randomnize(um ,k,randu ,irandom,ih,jh)
call randomnize(u0 ,k,randu ,irandom,ih,jh)
call randomnize(vm ,k,randu ,irandom,ih,jh)
call randomnize(v0 ,k,randu ,irandom,ih,jh)
call randomnize(wm ,k,randu ,irandom,ih,jh)
call randomnize(w0 ,k,randu ,irandom,ih,jh)
end do
svprof = 0.
if(myid==0)then
if (nsv>0) then
open (ifinput,file='scalar.inp.'//cexpnr)
read (ifinput,'(a80)') chmess
read (ifinput,'(a80)') chmess
do k=1,kmax
read (ifinput,*) &
height (k), &
(svprof (k,n),n=1,nsv)
end do
open (ifinput,file='scalar.inp.'//cexpnr)
write (6,*) 'height sv(1) --------- sv(nsv) '
do k=kmax,1,-1
write (6,*) &
height (k), &
(svprof (k,n),n=1,nsv)
end do
end if
end if ! end if myid==0
call MPI_BCAST(wsvsurf,nsv ,MY_REAL ,0,comm3d,mpierr)
call MPI_BCAST(svprof ,k1*nsv,MY_REAL ,0,comm3d,mpierr)
do k=1,kmax
do j=1,j2
do i=1,i2
do n=1,nsv
sv0(i,j,k,n) = svprof(k,n)
svm(i,j,k,n) = svprof(k,n)
end do
end do
end do
end do
!-----------------------------------------------------------------
! 2.2 Initialize surface layer and base profiles
!-----------------------------------------------------------------
select case(isurf)
case(1)
tskin = thls
tskinm = tskin
tsoilm = tsoil
phiwm = phiw
Wlm = Wl
case(2)
tskin = thls
case(3,4)
thls = thlprof(1)
qts = qtprof(1)
tskin = thls
qskin = qts
case(10)
call initsurf_user
end select
! Set initial Obukhov length to -0.1 for iteration
obl = -0.1
oblav = -0.1
call qtsurf
dthldz = (thlprof(1) - thls)/zf(1)
thvs = thls * (1. + (rv/rd - 1.) * qts)
if(lhetero) thvs_patch = thvs !Needed for initialization: thls_patch and qt_patch not yet calculated
u0av(1) = uprof(1)
thl0av(1) = thlprof(1)
svs = svprof(1,:)
call baseprofs ! call baseprofs before thermodynamics
call boundary
call thermodynamics
call surface
! Gradients at the top are now calculated in modboundary, every timestep
! dtheta = (thlprof(kmax)-thlprof(kmax-1)) / dzh(kmax)
! dqt = (qtprof (kmax)-qtprof (kmax-1)) / dzh(kmax)
! do n=1,nsv
! dsv(n) = (svprof(kmax,n)-svprof(kmax-1,n)) / dzh(kmax)
! end do
call boundary
call thermodynamics
else !if lwarmstart
call readrestartfiles
um = u0
vm = v0
wm = w0
thlm = thl0
qtm = qt0
svm = sv0
e12m = e120
call calc_halflev
exnf = (presf/pref0)**(rd/cp)
exnh = (presh/pref0)**(rd/cp)
do j=2,j1
do i=2,i1
do k=2,k1
thv0h(i,j,k) = (thl0h(i,j,k)+rlv*ql0h(i,j,k)/(cp*exnh(k))) &
*(1+(rv/rd-1)*qt0h(i,j,k)-rv/rd*ql0h(i,j,k))
end do
end do
end do
do j=2,j1
do i=2,i1
do k=1,k1
thv0(i,j,k) = (thl0(i,j,k)+rlv*ql0(i,j,k)/(cp*exnf(k))) &
*(1+(rv/rd-1)*qt0(i,j,k)-rv/rd*ql0(i,j,k))
end do
end do
end do
thvh=0.
call slabsum(thvh,1,k1,thv0h,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1) ! redefine halflevel thv using calculated thv
thvh = thvh/ijtot
thvf = 0.0
call slabsum(thvf,1,k1,thv0,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
thvf = thvf/ijtot
u0av = 0.0
v0av = 0.0
thl0av = 0.0
th0av = 0.0
qt0av = 0.0
ql0av = 0.0
sv0av = 0.
call slabsum(u0av ,1,k1,u0 ,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call slabsum(v0av ,1,k1,v0 ,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call slabsum(thl0av,1,k1,thl0,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call slabsum(qt0av ,1,k1,qt0 ,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call slabsum(ql0av ,1,k1,ql0 ,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
do n=1,nsv
call slabsum(sv0av(1,n),1,k1,sv0(1,1,1,n),2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
end do
u0av = u0av /ijtot + cu
v0av = v0av /ijtot + cv
thl0av = thl0av/ijtot
qt0av = qt0av /ijtot
ql0av = ql0av /ijtot
sv0av = sv0av /ijtot
th0av = thl0av + (rlv/cp)*ql0av/exnf
thvh(1) = th0av(1)*(1+(rv/rd-1)*qt0av(1)-rv/rd*ql0av(1)) ! override first level
do k=1,k1
rhof(k) = presf(k)/(rd*thvf(k)*exnf(k))
end do
! CvH - only do this for fixed timestepping. In adaptive dt comes from restartfile
if(ladaptive .eqv. .false.) rdt=dtmax
call baseprofs !call baseprofs
end if
!-----------------------------------------------------------------
! 2.1 read and initialise fields
!-----------------------------------------------------------------
if(myid==0)then
open (ifinput,file='lscale.inp.'//cexpnr)
read (ifinput,'(a80)') chmess
read (ifinput,'(a80)') chmess
write(6,*) ' height u_geo v_geo subs ' &
,' dqtdx dqtdy dqtdtls thl_rad '
do k=1,kmax
read (ifinput,*) &
height (k), &
ug (k), &
vg (k), &
wfls (k), &
dqtdxls(k), &
dqtdyls(k), &
dqtdtls(k), &
thlpcar(k)
end do
close(ifinput)
do k=kmax,1,-1
write (6,'(3f7.1,5e12.4)') &
height (k), &
ug (k), &
vg (k), &
wfls (k), &
dqtdxls(k), &
dqtdyls(k), &
dqtdtls(k), &
thlpcar(k)
end do
end if ! end myid==0
! MPI broadcast variables read in
call MPI_BCAST(ug ,kmax,MY_REAL ,0,comm3d,mpierr)
call MPI_BCAST(vg ,kmax,MY_REAL ,0,comm3d,mpierr)
call MPI_BCAST(wfls ,kmax,MY_REAL ,0,comm3d,mpierr)
call MPI_BCAST(dqtdxls ,kmax,MY_REAL ,0,comm3d,mpierr)
call MPI_BCAST(dqtdyls ,kmax,MY_REAL ,0,comm3d,mpierr)
call MPI_BCAST(dqtdtls ,kmax,MY_REAL ,0,comm3d,mpierr)
call MPI_BCAST(thlpcar ,kmax,MY_REAL ,0,comm3d,mpierr)
!-----------------------------------------------------------------
! 2.3 make large-scale horizontal pressure gradient
!-----------------------------------------------------------------
!******include rho if rho = rho(z) /= 1.0 ***********
do k=1,kmax
dpdxl(k) = om23_gs*vg(k)
dpdyl(k) = -om23_gs*ug(k)
end do
!-----------------------------------------------------------------
! 2.5 make large-scale horizontal gradients
!-----------------------------------------------------------------
whls(1) = 0.0
do k=2,kmax
whls(k) = ( wfls(k)*dzf(k-1) + wfls(k-1)*dzf(k) )/(2*dzh(k))
end do
whls(k1) = (wfls(kmax)+0.5*dzf(kmax)*(wfls(kmax)-wfls(kmax-1)) &
/dzh(kmax))
!******include rho if rho = rho(z) /= 1.0 ***********
if (llsadv) then
if (myid==0) stop 'llsadv should not be used anymore. Large scale gradients were calculated in a non physical way (and lmomsubs had to be set to true to retain conservation of mass)'
end if
dudxls = 0.0
dudyls = 0.0
dvdxls = 0.0
dvdyls = 0.0
dthldxls = 0.0
dthldyls = 0.0
idtmax = floor(dtmax/tres)
btime = timee
timeleft=ceiling(runtime/tres)
dt_lim = timeleft
rdt = real(dt)*tres
ntrun = 0
rtimee = real(timee)*tres
ntimee = nint(timee/dtmax)
itrestart = floor(trestart/tres)
tnextrestart = btime + itrestart
deallocate (height,th0av)
end subroutine readinitfiles
subroutine readrestartfiles
use modsurfdata, only : ustar,thlflux,qtflux,svflux,dthldz,dqtdz,ps,thls,qts,thvs,oblav,&
tsoil,phiw,tskin,Wl,isurf,ksoilmax,Qnet,swdavn,swuavn,lwdavn,lwuavn,nradtime,&
obl,xpatches,ypatches,ps_patch,thls_patch,qts_patch,thvs_patch,oblpatch,lhetero,qskin
use modraddata, only: iradiation, useMcICA
use modfields, only : u0,v0,w0,thl0,qt0,ql0,ql0h,e120,dthvdz,presf,presh,sv0,tmp0,esl,qvsl,qvsi
use modglobal, only : i1,i2,ih,j1,j2,jh,k1,dtheta,dqt,dsv,startfile,timee,&
tres,ifinput,nsv,dt
use modmpi, only : cmyid
use modsubgriddata, only : ekm
character(50) :: name
integer i,j,k,n
!********************************************************************
! 1.0 Read initfiles
!-----------------------------------------------------------------
name = startfile
name(5:5) = 'd'
name(12:19)=cmyid
write(6,*) 'loading ',name
open(unit=ifinput,file=name,form='unformatted', status='old')
read(ifinput) (((u0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
! u0 = u0-cu
read(ifinput) (((v0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
! v0 = v0-cv
read(ifinput) (((w0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
read(ifinput) (((thl0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
read(ifinput) (((qt0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
read(ifinput) (((ql0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
read(ifinput) (((ql0h (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
read(ifinput) (((e120 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
read(ifinput) (((dthvdz(i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
read(ifinput) (((ekm (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
read(ifinput) (((tmp0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
read(ifinput) (((esl (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
read(ifinput) (((qvsl (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
read(ifinput) (((qvsi (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
read(ifinput) ((ustar (i,j ),i=1,i2 ),j=1,j2 )
read(ifinput) ((thlflux (i,j ),i=1,i2 ),j=1,j2 )
read(ifinput) ((qtflux (i,j ),i=1,i2 ),j=1,j2 )
read(ifinput) ((dthldz(i,j ),i=1,i2 ),j=1,j2 )
read(ifinput) ((dqtdz (i,j ),i=1,i2 ),j=1,j2 )
read(ifinput) ( presf ( k) ,k=1,k1)
read(ifinput) ( presh ( k) ,k=1,k1)
read(ifinput) ps,thls,qts,thvs,oblav
read(ifinput) dtheta,dqt,timee,dt,tres
read(ifinput) ((obl (i,j ),i=1,i2 ),j=1,j2 )
read(ifinput) ((tskin(i,j ),i=1,i2 ),j=1,j2 )
read(ifinput) ((qskin(i,j ),i=1,i2 ),j=1,j2 )
if(lhetero) then
read(ifinput) ((ps_patch (i,j),i=1,xpatches),j=1,ypatches)
read(ifinput) ((thls_patch(i,j),i=1,xpatches),j=1,ypatches)
read(ifinput) ((qts_patch (i,j),i=1,xpatches),j=1,ypatches)
read(ifinput) ((thvs_patch(i,j),i=1,xpatches),j=1,ypatches)
read(ifinput) ((oblpatch (i,j),i=1,xpatches),j=1,ypatches)
endif
close(ifinput)
if (nsv>0) then
name(5:5) = 's'
write(6,*) 'loading ',name
open(unit=ifinput,file=name,form='unformatted')
read(ifinput) ((((sv0(i,j,k,n),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1),n=1,nsv)
read(ifinput) (((svflux(i,j,n),i=1,i2),j=1,j2),n=1,nsv)
read(ifinput) (dsv(n),n=1,nsv)
read(ifinput) timee
close(ifinput)
end if
if (isurf == 1) then
name(5:5) = 'l'
write(6,*) 'loading ',name
open(unit=ifinput,file=name,form='unformatted')
read(ifinput) (((tsoil(i,j,k),i=1,i2),j=1,j2),k=1,ksoilmax)
read(ifinput) (((phiw(i,j,k),i=1,i2),j=1,j2),k=1,ksoilmax)
read(ifinput) ((tskin(i,j),i=1,i2),j=1,j2)
read(ifinput) ((Wl(i,j),i=1,i2),j=1,j2)
read(ifinput) ((Qnet(i,j),i=1,i2),j=1,j2)
if(iradiation == 1 .and. useMcICA) then
read(ifinput) (((swdavn(i,j,n),i=1,i2),j=1,j2),n=1,nradtime)
read(ifinput) (((swuavn(i,j,n),i=1,i2),j=1,j2),n=1,nradtime)
read(ifinput) (((lwdavn(i,j,n),i=1,i2),j=1,j2),n=1,nradtime)
read(ifinput) (((lwuavn(i,j,n),i=1,i2),j=1,j2),n=1,nradtime)
end if
read(ifinput) timee
close(ifinput)
end if
end subroutine readrestartfiles
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine writerestartfiles
use modsurfdata,only: ustar,thlflux,qtflux,svflux,dthldz,dqtdz,ps,thls,qts,thvs,oblav,&
tsoil,phiw,tskin,Wl,ksoilmax,isurf,ksoilmax,Qnet,swdavn,swuavn,lwdavn,lwuavn,nradtime,&
obl,xpatches,ypatches,ps_patch,thls_patch,qts_patch,thvs_patch,oblpatch,lhetero,qskin
use modraddata, only: iradiation, useMcICA
use modfields, only : u0,v0,w0,thl0,qt0,ql0,ql0h,e120,dthvdz,presf,presh,sv0,tmp0,esl,qvsl,qvsi
use modglobal, only : i1,i2,ih,j1,j2,jh,k1,dsv,itrestart,tnextrestart,dt_lim,rtimee,timee,tres,cexpnr,&
rtimee,rk3step,ifoutput,nsv,timeleft,dtheta,dqt,dt
use modmpi, only : cmyid,myid
use modsubgriddata, only : ekm
implicit none
logical :: lexitnow = .false.
integer imin,ihour
integer i,j,k,n
character(50) name,linkname
if (timee == 0) return
if (rk3step /=3) return
name = 'exit_now.'//cexpnr
inquire(file=trim(name), EXIST=lexitnow)
if (timee<tnextrestart) dt_lim = min(dt_lim,tnextrestart-timee)
if (timee>=tnextrestart .or. lexitnow) then
tnextrestart = tnextrestart+itrestart
ihour = floor(rtimee/3600)
imin = floor((rtimee-ihour * 3600) /3600. * 60.)
name = 'initd h m .'
write (name(6:7) ,'(i2.2)') ihour
write (name(9:10) ,'(i2.2)') imin
name(12:19)= cmyid
name(21:23)= cexpnr
open (ifoutput,file=name,form='unformatted',status='replace')
write(ifoutput) (((u0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) (((v0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) (((w0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) (((thl0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) (((qt0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) (((ql0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) (((ql0h (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) (((e120 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) (((dthvdz(i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) (((ekm (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) (((tmp0 (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) (((esl (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) (((qvsl (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) (((qvsi (i,j,k),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1)
write(ifoutput) ((ustar (i,j ),i=1,i2 ),j=1,j2 )
write(ifoutput) ((thlflux (i,j ),i=1,i2 ),j=1,j2 )
write(ifoutput) ((qtflux (i,j ),i=1,i2 ),j=1,j2 )
write(ifoutput) ((dthldz(i,j ),i=1,i2 ),j=1,j2 )
write(ifoutput) ((dqtdz (i,j ),i=1,i2 ),j=1,j2 )
write(ifoutput) ( presf ( k) ,k=1,k1)
write(ifoutput) ( presh ( k) ,k=1,k1)
write(ifoutput) ps,thls,qts,thvs,oblav
write(ifoutput) dtheta,dqt,timee, dt,tres
write(ifoutput) ((obl (i,j ),i=1,i2 ),j=1,j2 )
write(ifoutput) ((tskin(i,j ),i=1,i2 ),j=1,j2 )
write(ifoutput) ((qskin(i,j ),i=1,i2 ),j=1,j2 )
if(lhetero) then
write(ifoutput) ((ps_patch (i,j),i=1,xpatches),j=1,ypatches)
write(ifoutput) ((thls_patch(i,j),i=1,xpatches),j=1,ypatches)
write(ifoutput) ((qts_patch (i,j),i=1,xpatches),j=1,ypatches)
write(ifoutput) ((thvs_patch(i,j),i=1,xpatches),j=1,ypatches)
write(ifoutput) ((oblpatch (i,j),i=1,xpatches),j=1,ypatches)
endif
close (ifoutput)
linkname = name
linkname(6:11) = "latest"
call system("ln -sf "//name //" "//linkname)
if (nsv>0) then
name = 'inits h m .'
write (name(6:7) ,'(i2.2)') ihour
write (name(9:10) ,'(i2.2)') imin
name(12:19) = cmyid
name(21:23) = cexpnr
open (ifoutput,file=name,form='unformatted')
write(ifoutput) ((((sv0(i,j,k,n),i=2-ih,i1+ih),j=2-jh,j1+jh),k=1,k1),n=1,nsv)
write(ifoutput) (((svflux(i,j,n),i=1,i2),j=1,j2),n=1,nsv)
write(ifoutput) (dsv(n),n=1,nsv)
write(ifoutput) timee
close (ifoutput)
linkname = name
linkname(6:11) = "latest"
call system("ln -sf "//name //" "//linkname)
end if
if (isurf == 1) then
name = 'initl h m .'
write (name(6:7) ,'(i2.2)') ihour
write (name(9:10) ,'(i2.2)') imin
name(12:19) = cmyid
name(21:23) = cexpnr
open (ifoutput,file=name,form='unformatted')
write(ifoutput) (((tsoil(i,j,k),i=1,i2),j=1,j2),k=1,ksoilmax)
write(ifoutput) (((phiw(i,j,k),i=1,i2),j=1,j2),k=1,ksoilmax)
write(ifoutput) ((tskin(i,j),i=1,i2),j=1,j2)
write(ifoutput) ((Wl(i,j),i=1,i2),j=1,j2)
write(ifoutput) ((Qnet(i,j),i=1,i2),j=1,j2)
if(iradiation == 1 .and. useMcICA) then
write(ifoutput) (((swdavn(i,j,n),i=1,i2),j=1,j2),n=1,nradtime)
write(ifoutput) (((swuavn(i,j,n),i=1,i2),j=1,j2),n=1,nradtime)
write(ifoutput) (((lwdavn(i,j,n),i=1,i2),j=1,j2),n=1,nradtime)
write(ifoutput) (((lwuavn(i,j,n),i=1,i2),j=1,j2),n=1,nradtime)
end if
write(ifoutput) timee
close (ifoutput)
linkname = name
linkname(6:11) = "latest"
call system("ln -sf "//name //" "//linkname)
end if
if (lexitnow) then
timeleft = 0 !jump out of the time loop
end if
if (lexitnow .and. myid == 0 ) then
open(1, file=trim(name), status='old')
close(1,status='delete')
write(*,*) 'Stopped at t=',rtimee
end if
if (myid==0) then
write(*,'(A,F15.7,A,I4)') 'dump at time = ',rtimee,' unit = ',ifoutput
end if
end if
end subroutine writerestartfiles
subroutine exitmodules
use modfields, only : exitfields
use modglobal, only : exitglobal
use modmpi, only : exitmpi
use modboundary, only : exitboundary
use modmicrophysics, only : exitmicrophysics
use modpois, only : exitpois
use modtimedep, only : exittimedep
use modradiation, only : exitradiation
use modsubgrid, only : exitsubgrid
use modsurface, only : exitsurface
use modthermodynamics, only : exitthermodynamics
call exittimedep
call exitthermodynamics
call exitsurface
call exitsubgrid
call exitradiation
call exitpois
call exitmicrophysics
call exitboundary
call exitfields
call exitglobal
call exitmpi
end subroutine exitmodules
!----------------------------------------------------------------
subroutine randomnize(field,klev,ampl,ir,ihl,jhl)
! Adds (pseudo) random noise with given amplitude to the field at level k