-
Notifications
You must be signed in to change notification settings - Fork 5
/
Copy pathmetdata_ecmwf_uioformat.f90
1854 lines (1624 loc) · 66.5 KB
/
metdata_ecmwf_uioformat.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
!//=========================================================================
!// Oslo CTM3
!//=========================================================================
!// Based on UCI CTM core p-7.1 (1/2013).
!//
!// Ole Amund Sovde, April 2015
!//=========================================================================
!// Routine for reading meteorological data from ECMWF, in UiO EC-files.
!//=========================================================================
module metdata_ecmwf
!//-----------------------------------------------------------------------
!// MODULE: metdata_ecmwf
!// DESCRIPTION: Routine for reading meteorological data from ECMWF,
!// stored on UiO EC-file format.
!//
!// Contains
!// subroutine update_metdata
!// subroutine CFRMIN
!// subroutine CIWMIN
!// subroutine fluxfilter2
!// subroutine r4data2mpblocks
!//-----------------------------------------------------------------------
use cmn_size, only: LPAR
!//-----------------------------------------------------------------------
implicit none
!//-----------------------------------------------------------------------
integer :: LMAP(LPAR+1)
!//-----------------------------------------------------------------------
character(len=*), parameter, private :: f90file = 'metdata_ecmwf.f90'
!//-----------------------------------------------------------------------
private
public update_metdata
save
!//-----------------------------------------------------------------------
contains
!//-----------------------------------------------------------------------
subroutine update_metdata(ILOOP, DTMET, NMET)
!//---------------------------------------------------------------------
!// Read input data. ECMWF data has three files, one for spectral
!// data (.b01), one for 2-d gridpoint data (.b03) and one for
!// 3-d gridpoint data (.b02).
!// Spectral transformation is done in this routine, and sumation of
!// physical fields from Txx ---> T42, T21 or other resolution.
!// Allow several different formats for 2-d gridpoint files by labelling
!// the record numbers we need (we need N15R records)
!//
!// This routine assumes meteorology is to be updated each NMET. E.g.
!// ERA-40 had 6-hour fields, and was tweaked to be read every other
!// NMET, but with the double time step. This has been removed for now.
!//
!// Ole Amund Sovde, April 2015
!//---------------------------------------------------------------------
use cmn_precision, only: r8, r4
use cmn_size, only: IPARW, JPARW, LPARW, IPAR, JPAR, IDGRD, JDGRD, &
LPAR, LWEPAR, LWDPAR, &
FFTPTS, NRFFT, NRNM, NMMAX, NTRUNW, LOSLOCHEM, LDUST
use cmn_ctm, only: JYEAR, JMON, JDAY, GMTAU, LMMAP, XLMMAP, TMET, &
ALP, ALPV, TRIG, IFAX, LDEG, ZDEGI, ZDEGJ, IMAP, JMAP, &
XYZA, XYZB, IMEPZ, DD, SS, ETAAW, ETABW, ETAA, ETAB, &
WGLYE, WGLYG, YDGRD, AREAXYW, AREAXY, PLAND
use cmn_met, only: P, U, V, T, Q, CWETN, CWETE, CWETD, CENTU, CENTD, &
PRECCNV, PRECLS, ZOFLE, &
CLDFR, CLDLWC, CLDIWC, &
SLH, SHF, SMF, SFT, SFQ, SFU, SFV, BLH, SA, &
PVU, UMS, VMS, CI, SD, SMLT, ES, SNFL, PhotActRad, MSLP, USTR, &
PMEAN, PMEANW, MPATH1, MPATH2, MFILE3, MYEAR, &
metTYPE, metCYCLE, metREVNR, HnativeRES, VRESW
use cmn_parameters, only: M_AIR, R_UNIV, R_AIR, R_H2O, A0, CPI
use cloudjx, only: SWSTORE, CLDSTORE, TYPSTORE, RAN4, IRAN0, &
LCLDAVG, LCLDQMD, LCLDQMN, LCLDRANA, LCLDRANQ
use regridding, only: TRUNG4, TRUNG8
use utilities, only: ctmExitC, CFRMIN, CIWMIN
use dust_oslo, only: dust_set_ssrd, dust_set_strd, dust_set_SWVL1
use physics_oslo, only: get_pvu, ijlw2lij
!//---------------------------------------------------------------------
implicit none
!//---------------------------------------------------------------------
!// Input
integer, intent(in) :: ILOOP, & ! <=0, init, > 0 ordinary run
NMET ! the meteorological timestep
real(r8), intent(in) :: DTMET ! meteorological time step [s]
!//---------------------------------------------------------------------
logical, parameter :: VERBOSE = .false.
real(r8), parameter :: EPS = 0.01_r8 ! minimum cloud fraction
!// Limits for convective fluxes
real(r8), parameter :: mincwetelim = 3.8e-4_r8
real(r8), parameter :: maxcwetdlim = -2.0e-5_r8
real(r8), parameter :: mincdetulim = 1.e-4_r8
real(r8), parameter :: mincdetdlim = 1.e-6_r8
!// Local parameters
real(r8) :: &
VEDGE, BAND, DETA, DETB, DELP, &
ZCOS, ZDT, SFTD, LV, ESAT, SDEN, &
PT, PB, &
DELZ, LWC, IWC, TFACT, QMIN, &
DMASS, ZAREA, PSRF, &
POFLE(LPAR+1)
!// Indices
integer :: I,J,II,JJ,L,LL
integer :: NF, NrOf3Dfields, NrOf2Dfields, ioerr
!// To be read from file
integer :: FLD_ID, ISEC(16)
logical :: fex
logical :: LUV !// Flag if U or V are calculated by SPE2GP
logical :: LPVU !// Need to calculate PVU?
character(len=120) :: FNAME1, FNAME2, FNAME3
character(len=4) :: CYYx
!// To check if RAIN is split into CONVRAIN/LSRAIN or not
!// If 217 and 218 are found instead of 216 we have CONVRAIN/LSRAIN:
logical :: found216 !// field 216 exist; only total rain
real(r8) :: CNVFRAC, CONVTOT, RAINTOT
!//---------------------------------------------------------------------
!// Allocatable arrays - single precision
real(r4), allocatable :: &
GRDATA2(:,:), GRDATA2HI(:,:), SPWK5(:,:), SPWK6(:,:)
real(r4), allocatable :: &
GRDATA3(:,:,:), GRDATA3HI(:,:,:)
real(r4), allocatable :: &
PSPEHI(:)
!// Allocatable arrays - double precision
real(r8), allocatable :: &
FFTAR(:,:), WORK(:), &
VTMP(:,:), UTMP(:)
real(r8), dimension(:,:), allocatable :: &
PW, EWSS, NSSS, LSPREC, CNVPREC, R8XY
real(r8), dimension(:,:,:), allocatable :: &
TW, QW, CLDFRW, CLDIWCW, CLDLWCW, CDETU, CDETD, &
ZOFLEW, WK3DT, WK3DU, WK3DV, UMSW, VMSW, &
TMPCRAIN, TMPLRAIN
!//---------------------------------------------------------------------
character(len=*), parameter :: subr = 'update_metdata'
!//---------------------------------------------------------------------
!// Allocate spectral arrays
allocate( FFTAR(FFTPTS,NRFFT), WORK(FFTPTS*NRFFT), &
SPWK5(NRNM,LPARW), SPWK6(NRNM,LPARW), PSPEHI(NRNM) )
!// Allocate 3D arrays - native resolution
allocate( VTMP(IPARW,JPARW), UTMP(IPARW), &
WK3DU(IPARW,JPARW,LPARW), WK3DV(IPARW,JPARW,LPARW), &
WK3DT(IPARW,JPARW,LPARW), &
ZOFLEW(LPAR+1,IPARW,JPARW), QW(IPARW,JPARW,LPAR), &
TW(IPARW,JPARW,LPAR), CLDFRW(IPARW,JPARW,LWEPAR), &
CLDIWCW(IPARW,JPARW,LWEPAR), CLDLWCW(IPARW,JPARW,LWEPAR), &
GRDATA3HI(IPARW,JPARW,LPARW), &
UMSW(IPARW,JPARW,LPARW), VMSW(IPARW,JPARW,LPARW) )
!// Allocate 3D arrays - window resolution (IPAR/JPAR)
allocate( TMPCRAIN(IPAR,JPAR,LPARW), TMPLRAIN(IPAR,JPAR,LPARW), &
CDETU(IPAR,JPAR,LWEPAR), CDETD(IPAR,JPAR,LWDPAR), &
GRDATA3(IPAR,JPAR,LPARW) )
!// Allocate 2D arrays - native resolution
allocate( PW(IPARW,JPARW), GRDATA2HI(IPARW,JPARW) )
!// Allocate 2D arrays - window resolution (IPAR/JPAR)
allocate( LSPREC(IPAR,JPAR),CNVPREC(IPAR,JPAR), &
EWSS(IPAR,JPAR), NSSS(IPAR,JPAR), &
GRDATA2(IPAR,JPAR), R8XY(IPAR,JPAR) )
!//---------------------------------------------------------------------
!// Initialize
found216 = .false.
LUV = .false.
!// Time step for meteorological data is DTMET
ZDT = 1._r8 / DTMET
!//Define minimum humidity QMIN in kg/kg
QMIN = 3.e-6_r8 * 18._r8 / M_AIR
!//UPDATE
!// locate the position of random number sequence based on year/day/hour
IRAN0 = 1 + 3*(JYEAR - 1900) + 7*JDAY + 11*nint(GMTAU)
if (metTYPE(1:10) .eq. 'ECMWF_oIFS') then
!// ECMWF openIFS EC-format
write(MPATH2(1:4),'(i4.4)') MYEAR
write(MFILE3(7:13),'(i4.4,A3)') MYEAR,TMET
if (metREVNR .eq. 1) then
write(CYYx(1:4),'(a1,i2.2,a1)') 'C',metCYCLE,'a'
else
write(6,'(a)') '*** Unknown metdata revision number!'
write(6,'(a,2i4)') ' metCYCLE / metREVNR: ',metCYCLE,metREVNR
stop
end if
FNAME3 = trim(MFILE3)
FNAME1 = trim(MPATH1)//trim(MPATH2)//trim(FNAME3)//'.b01'
FNAME2 = trim(MPATH1)//trim(MPATH2)//trim(FNAME3)//'.b02'
FNAME3 = trim(MPATH1)//trim(MPATH2)//trim(FNAME3)//'.b03'
else if (trim(metTYPE) .eq. 'ECMWF_IFS') then
!// Old style metdata names
write(MPATH2(1:4),'(i4)') MYEAR
write(MFILE3(3:5),'(A3)') TMET
FNAME3 = trim(MFILE3)
FNAME1 = trim(MPATH1)//trim(MPATH2)//trim(FNAME3)//'.b01'
FNAME2 = trim(MPATH1)//trim(MPATH2)//trim(FNAME3)//'.b02'
FNAME3 = trim(MPATH1)//trim(MPATH2)//trim(FNAME3)//'.b03'
else
write(6,'(a)') f90file//':'//subr// &
': Not set up for metTYPE: '//trim(metTYPE)
stop
end if
!//---------------------------------------------------------------------
!// Initial step - setup dry airmass (P() and Q())
!//---------------------------------------------------------------------
if (ILOOP .le. 0) then
!// Set up level weightings if vertical resolution degraded
do LL = LPARW+1, 1, -1
L = LMMAP(LL)
LMAP(L) = LL
end do
!//print out info, open files and read in data, transform and sum up
write(6,'(a)') 'Initializing meteorological data'
write(6,'(2x,a)') trim(FNAME1)
write(6,'(2x,a)') trim(FNAME2)
inquire(FILE=trim(FNAME1), exist=fex)
if (.not. fex) call ctmExitC('update_metdata: No such file: '//trim(FNAME1))
close(13)
open(13,FILE=trim(FNAME1),FORM='UNFORMATTED',STATUS='OLD',iostat=ioerr)
if (ioerr .ne. 0) call ctmExitC('update_metdata: ERROR on .b01 file')
inquire(FILE=trim(FNAME2), exist=fex)
if (.not. fex) call ctmExitC('update_metdata: No such file: '//trim(FNAME2))
close(14)
open(14,FILE=trim(FNAME2),FORM='UNFORMATTED',STATUS='OLD',iostat=ioerr)
if (ioerr .ne. 0) call ctmExitC('update_metdata: ERROR on .b02 file')
!// Pressure field
read(13) FLD_ID
if (FLD_ID .ne. 152) &
call ctmExitC('update_metdata: ERROR, input Data PSPEHI')
read(13) PSPEHI
!// Pressure; transform from spectral to grid point data
call SPE2GP(PSPEHI, NRNM, LUV, PW, IPARW, JPARW, 1, ALP, NMMAX, &
FFTAR, WORK, TRIG, IFAX, (NTRUNW+1), FFTPTS, NRFFT)
!// save high-res pressure Log(Ps) --> hPa
do J = 1, JPARW
do I = 1, IPARW
PW(I,J) = EXP(PW(I,J)) * 1.e-2_r8
end do
end do
if (ldeg) then
!// Truncate to window resolution
call TRUNG8(PW, P, ZDEGI, ZDEGJ, IMAP, JMAP, IDGRD, &
JDGRD, IPARW, JPARW, IPAR, JPAR, 1, 1)
else
do J = 1,JPAR
do I = 1,IPAR
P(I,J) = PW(I,J)
end do
end do
end if
!// Temperature (originally not in the UCI code)
read(13) FLD_ID
if (FLD_ID .ne. 130) &
call ctmExitC('update_metdata: ERROR, input Data T')
read(13) SPWK5
!// Temperature; transform from spectral to grid point data
call SPE2GP(SPWK5, NRNM, LUV, WK3DT, IPARW, JPARW, LPARW, ALP, &
NMMAX, FFTAR, WORK, TRIG, IFAX, (NTRUNW+1), FFTPTS, NRFFT)
TW(:,:,:) = 0._r8
do L = 1,LPAR
do LL = LMAP(L),LMAP(L+1)-1
do J = 1,JPARW
do I = 1,IPARW
TW(I,J,L) = TW(I,J,L) + WK3DT(I,J,LL) * XLMMAP(LL)
end do
end do
end do
end do
if (ldeg) then
call TRUNG8(TW, T, ZDEGI, ZDEGJ, IMAP, JMAP, IDGRD, &
JDGRD, IPARW, JPARW, IPAR, JPAR, LPAR, 1)
else
do L = 1,LPAR
do J = 1,JPAR
do I = 1,IPAR
T(I,J,L) = TW(I,J,L)
end do
end do
end do
end if
!// Test to check temperature field
if (minval(WK3DT(:,:,1)) .le. 0.) then
print*,'update_metdata: Temperature is wrong'
print*,'MIN T',minval(WK3DT(:,:,1))
print*,'Maybe ALP=0? min/max ALP',minval(alp), maxval(alp)
stop
end if
!// Close file
close(13)
!// Water vapour
read(14) FLD_ID
if (FLD_ID .ne. 133) &
call ctmExitC('update_metdata: ERROR, input Data Q')
read(14) GRDATA3HI
!// Collapse layers
QW(:,:,:) = 0._r8
do L = 1,LPAR
do LL = LMAP(L),LMAP(L+1)-1
do J = 1,JPARW
do I = 1,IPARW
QW(I,J,L) = QW(I,J,L) + GRDATA3HI(I,J,LL) * XLMMAP(LL)
end do
end do
end do
end do
!// There may be some entries of small negative numbers
if (minval(QW) .lt. 0._r8) then
print*,'update_metdata: min QW:',minval(QW),', setting to',QMIN
QW(:,:,:) = max(QW(:,:,:), QMIN)
end if
if (ldeg) then
call TRUNG8(QW, Q, ZDEGI, ZDEGJ, IMAP, JMAP, IDGRD, &
JDGRD, IPARW, JPARW, IPAR, JPAR, LPAR,1)
else
Q(:,:,:) = QW(:,:,:)
end if
!// Close file
close(14)
!// Polar cap filtering
call EPZ_TQ(T, Q, XYZA, XYZB, P, IMEPZ,IPAR,JPAR,LPAR,IPAR,JPAR,LPAR)
call EPZ_P(P, IMEPZ, IPAR,JPAR, IPAR, JPAR)
!// Initialization is done - these steps are for all time steps
return
end if !// if (ILOOP .le. 0) then
!//---------------------------------------------------------------------
!// All time steps - update all meteorological fields
!//---------------------------------------------------------------------
if (verbose) write(6,'(a,i5)') f90file//':'//subr// &
': Reading new metdata JDAY: '//TMET//', NMET:',NMET
!// Clear arrays
U(:,:,:) = 0._r8
V(:,:,:) = 0._r8
T(:,:,:) = 0._r8
Q(:,:,:) = 0._r8
CWETN(:,:,:) = 0._r8
CWETE(:,:,:) = 0._r8
CWETD(:,:,:) = 0._r8
CENTU(:,:,:) = 0._r8
CENTD(:,:,:) = 0._r8
CDETU(:,:,:) = 0._r8
CDETD(:,:,:) = 0._r8
PRECCNV(:,:,:) = 0._r8
PRECLS(:,:,:) = 0._r8
CLDFR(:,:,:) = 0._r8
CLDLWC(:,:,:) = 0._r8
CLDIWC(:,:,:) = 0._r8
!// No need to clear 2D arrays; they are fully updated below
SLH(:,:) = 0._r8
SHF(:,:) = 0._r8
SMF(:,:) = 0._r8
SFT(:,:) = 0._r8
SFQ(:,:) = 0._r8
SFU(:,:) = 0._r8
SFV(:,:) = 0._r8
BLH(:,:) = 0._r8
MSLP(:,:) = 0._r8
SA(:,:) = 0._r8
!// New day? If so, open new files
if (NMET .eq. 1) then
inquire(FILE=trim(FNAME1), exist=fex)
if (.not. fex) call ctmExitC('update_metdata: No such file: '//trim(FNAME1))
close(13)
open(13,FILE=trim(FNAME1),FORM='UNFORMATTED',STATUS='OLD',iostat=ioerr)
if (ioerr .ne. 0) call ctmExitC('update_metdata: ERROR on .b01 file')
inquire(FILE=trim(FNAME2), exist=fex)
if (.not. fex) call ctmExitC('update_metdata: No such file: '//trim(FNAME2))
close(14)
open(14,FILE=trim(FNAME2),FORM='UNFORMATTED',STATUS='OLD',iostat=ioerr)
if (ioerr .ne. 0) call ctmExitC('update_metdata: ERROR on .b02 file')
inquire(FILE=trim(FNAME3), exist=fex)
if (.not. fex) call ctmExitC('update_metdata: No such file: '//trim(FNAME3))
close(15)
open(15,FILE=trim(FNAME3),FORM='UNFORMATTED',STATUS='OLD',iostat=ioerr)
if (ioerr .ne. 0) call ctmExitC('update_metdata: ERROR on .b03 file')
end if
!//---------------------------------------------------------------------
!// SPECTRAL DATA
!//---------------------------------------------------------------------
!// Read in spectral data - error check only on Pressure since it
!// is assumed that the ordering of fields is fine.
!// Pressure field
read(13) FLD_ID
if (FLD_ID .ne. 152) &
call ctmExitC('update_metdata: ERROR, input Data PSPEHI')
read(13) PSPEHI
!// Pressure; transform from spectral to grid point data
call SPE2GP(PSPEHI, NRNM, LUV, PW, IPARW, JPARW, 1, ALP, &
NMMAX, FFTAR, WORK, TRIG, IFAX, (NTRUNW+1), FFTPTS, NRFFT)
!// save high-res pressure Log(Ps) --> hPa
do J = 1,JPARW
do I = 1,IPARW
PW(I,J) = EXP(PW(I,J))*1.e-2_r8
end do
end do
if (ldeg) then
Call TRUNG8(PW, P, ZDEGI, ZDEGJ, IMAP, JMAP, IDGRD, &
JDGRD, IPARW, JPARW, IPAR, JPAR, 1, 1)
else
do J = 1, JPAR
do I = 1, IPAR
P(I,J) = PW(I,J)
end do
end do
end if
!// Temperature
read(13) FLD_ID
If (FLD_ID .ne. 130) &
call ctmExitC('update_metdata: ERROR, input Data PSPEHI')
read(13) SPWK5
!// Temperature; transform from spectral to grid point data
call SPE2GP(SPWK5, NRNM, LUV, WK3DT, IPARW, JPARW, LPARW, ALP, &
NMMAX, FFTAR, WORK, TRIG, IFAX, (NTRUNW+1), FFTPTS, NRFFT)
TW(:,:,:) = 0._r8 !// Clear array
do L = 1,LPAR
do LL = LMAP(L),LMAP(L+1)-1
do J = 1,JPARW
do I = 1,IPARW
TW(I,J,L) = TW(I,J,L) + WK3DT(I,J,LL) * XLMMAP(LL)
end do
end do
end do
end do
if (ldeg) then
call TRUNG8(TW, T, ZDEGI, ZDEGJ, IMAP, JMAP, IDGRD, &
JDGRD, IPARW, JPARW, IPAR, JPAR, LPAR, 1)
else
do L = 1, LPAR
do J = 1, JPAR
do I = 1, IPAR
T(I,J,L) = TW(I,J,L)
end do
end do
end do
end if
!// Vorticity and Divergence - find U and V
read(13) FLD_ID
if (FLD_ID .ne. 138) &
call ctmExitC('update_metdata: ERROR, input Data 138')
read(13) SPWK5
read(13) FLD_ID
if (FLD_ID .ne. 155) &
call ctmExitC('update_metdata: ERROR, input Data 155')
read(13) SPWK6
!// Vorticity and Divergence - find U and V,
!// then do spectral to GP transform
call ZD2UV(SPWK5, SPWK6, WK3DU, WK3DV, DD, SS, IPARW, JPARW, LPARW, &
ALP, ALPV, NMMAX, NRNM, FFTAR, WORK, TRIG, IFAX, &
NTRUNW, NRFFT, FFTPTS)
!// scale from U=u*cos() to u
do L = 1, LPARW
DETA = ETAAW(L) - ETAAW(L+1)
DETB = ETABW(L) - ETABW(L+1)
do J = 1, JPARW
BAND = A0 * (WGLYE(J+1) - WGLYE(J))
ZCOS = 1._r8 / COS(WGLYG(J))
do I = 1, IPARW
UTMP(I) = WK3DU(I,J,L) * ZCOS * BAND
!// Save center values for m/s
UMSW(I,J,L) = WK3DU(I,J,L)*ZCOS !// Save U(m/s)
end do
!// Get edge values (unit conversion needs pressure on edge)
do I = 2, IPARW
DELP = DETA + DETB * (PW(I-1,J) + PW(I,J)) * 0.5_r8
WK3DU(I,J,L) = (UTMP(I-1) + UTMP(I)) * 0.5_r8 * DELP
end do
DELP = DETA + DETB * (PW(1,J) + PW(IPARW,J)) * 0.5_r8
WK3DU(1,J,L) = (UTMP(IPARW) + UTMP(1)) * 0.5_r8 * DELP
end do
end do
!// scale from V=v*cos() to v (remember V is on edge already)
do L = 1, LPARW
DETA = ETAAW(L) - ETAAW(L+1)
DETB = ETABW(L) - ETABW(L+1)
VTMP(:,:) = 0._r8 !// Needed for getting VMSW
!// WK3DV does not start at SP, but at the J=2. It means
!// WK3DV is flux out of box J, not into J.
!// Given all edge points, the number of boxes should be
!// JPARW+1, covering 1:JPARW+1. But both the poles should
!// have zero wind, so we only need JPARW-1 boxes.
!//
!// But in addition, the values at JPARW/2 and JPARW/2+1
!// are duplicated, both representing wind across Equator.
!// So WK3DV has the size JPARW instead of JPARW-1.
!//
!// A bit confusing, the CTM uses V as flux into J-box instead of
!// out of it. This will be taken care of at the end.
!// SH flux out of J-box:
do J = 1, JPARW/2
ZCOS = 1._r8 / COS(WGLYE(J+1))
do I = 1, IPARW
VTMP(I,J) = WK3DV(I,J,L) * ZCOS
end do
end do
!// Now VTMP covers the V up to Equator, starting at the
!// models index J=2, but stored in VTMP at J=1.
!// NH flux out of J-box (JPARW/2 and JPARW/2+1 are duplicated):
do J = JPARW/2 + 2, JPARW
ZCOS = 1._r8 / COS(WGLYE(J))
do I = 1, IPARW
VTMP(I,J-1) = WK3DV(I,J,L) * ZCOS
end do
end do
!// Thus, at J=JPARW, VTMP (with its indicing) represents NP
!// and is therefore zero.
!//change unit, m/s ==> Kg/s (100./G0 is done in pdyn0.f)
!// Save center values for m/s J=1 (assume V=0 at SP)
VMSW(:,1,L) = VTMP(:,1) * 0.5_r8
!// Map the VTMP back to correct model edge grid (from J to J+1),
!// i.e. converting from flux out of J-1 to flux into J.
do J = 2, JPARW
VEDGE = 2._r8 * CPI * A0 * COS(WGLYE(J))/real(IPARW, r8)
do I = 1, IPARW
DELP = DETA + DETB * (PW(I,J-1) + PW(I,J)) * 0.5_r8
!// Save center values for m/s J>1
VMSW(I,J,L) = (VTMP(I,J) + VTMP(I,J-1)) * 0.5_r8
WK3DV(I,J,L) = VTMP(I,J-1) * VEDGE * DELP
end do
end do
end do
if (ldeg) then
do L = 1, LPAR
do LL = LMAP(L),LMAP(L+1)-1
do I = 1, IPAR
II = IMAP(1,I)
do J = 1, JPAR
do JJ = 1, JDGRD
U(I,J,L) = U(I,J,L) + WK3DU(II,JMAP(JJ,J),LL)
end do
end do
end do
do J = 1, JPAR
JJ = JMAP(1,J)
do I = 1, IPAR
do II = 1, IDGRD
V(I,J,L) = V(I,J,L) + WK3DV(IMAP(II,I),JJ,LL)
end do
end do
end do
end do
end do
else
do L = 1, LPAR
do LL = LMAP(L),LMAP(L+1)-1
do J = 1, JPAR
do I = 1, IPAR
U(I,J,L) = U(I,J,L) + WK3DU(I,J,LL)
V(I,J,L) = V(I,J,L) + WK3DV(I,J,LL)
end do
end do
end do
end do
end if
!// VERTICAL VELOCITY [Pa/s] phony read in if present
!// Some of the oldest metdata also contain W (FLD_ID 135).
!// W is not used, so I have included a dummy read-in for that.
!read(13, iostat=ioerr) FLD_ID
!if (ioerr .eq. 0) then
! if (FLD_ID .ne. 135) then
! read(13) RSUM
! else
! backspace(13)
! end if
!end if !// if (ioerr .eq. 0) then
!//---------------------------------------------------------------------
!// 3-d GRID POINT DATA
!//---------------------------------------------------------------------
!// Read in 3-d Gridpoint data, check on field type
!// WATER VAPOUR
read(14) FLD_ID
if (FLD_ID .ne. 133) &
call ctmExitC('update_metdata: ERROR, input Data Q')
read(14) GRDATA3HI
!// Collapse layers
QW(:,:,:) = 0._r8
do L = 1, LPAR
do LL = LMAP(L), LMAP(L+1)-1
do J = 1, JPARW
do I = 1, IPARW
QW(I,J,L) = QW(I,J,L) + GRDATA3HI(I,J,LL) * XLMMAP(LL)
end do
end do
end do
end do
if (minval(QW) .lt. 0._r8) then
print*,'update_metdata: min QW:',minval(QW),', setting to',QMIN
QW(:,:,:) = max(QW(:,:,:), QMIN)
end if
if (ldeg) then
call TRUNG8(QW, Q, ZDEGI, ZDEGJ, IMAP, JMAP, IDGRD, &
JDGRD, IPARW, JPARW, IPAR, JPAR, LPAR, 1)
else
do L = 1, LPAR
do J = 1, JPAR
do I = 1, IPAR
Q(I,J,L) = QW(I,J,L)
end do
end do
end do
end if
!// Altitudes
!// R = 287.*(1-Q) + 461.5*Q -- assume 0.5% bl w.v.==> R = 288.
!// delta-z (m) = dln(P) * R * T / g where R/g = 288/9.81 = 29.36
!// CTM3: We use Q directly instead of assuming 0.5%;
!// Using Rd=287 and Tv=T*(1 + 0.6*q); 287/9.80665=29.26586
do J = 1, JPAR
do I = 1, IPAR
PSRF = P(I,J)
!// Surface pressure; must be set for DELZ to be calculated
POFLE(1) = PSRF
!// Height of box bottom, above sea level (acts as topography)
ZOFLE(1,I,J) = 16.e3_r8 * log10(1013.25_r8 / PMEAN(I,J))
if (ZOFLE(1,I,J) .ne. ZOFLE(1,I,J)) then
print*,'update_metdata: ZOFLE a',1,i,j,PMEAN(I,J),PMEANW(I,J)
print*,pmean(:,j)
stop
end if
do L = 2, LPAR + 1
!// Pressure of box bottom
POFLE(L) = ETAA(L) + ETAB(L) * PSRF
!DELZ = -29.36_r8 * T(I,J,L-1) * log(POFLE(L)/POFLE(L-1))
!// Thickness of layer (L-1) (remember ZOFLE starts at topography)
DELZ = -29.26586_r8 * T(I,J,L-1) * (1._r8 + 0.6_r8 * Q(I,J,L-1)) &
* log(POFLE(L)/POFLE(L-1))
!// Add DELZ of (L-1) to get box bottom height of L
ZOFLE(L,I,J) = ZOFLE(L-1,I,J) + DELZ
if (ZOFLE(L,I,J) .ne. ZOFLE(L,I,J)) then
print*,'update_metdata: ZOFLE b',l,i,j,delz,q(i,j,l-1),&
T(I,J,L-1),psrf
stop
end if
end do
end do
end do
!// Also set up ZOFLEW for metdata grid
do J = 1, JPARW
do I = 1, IPARW
PSRF = PW(I,J)
POFLE(1) = PSRF
ZOFLEW(1,I,J) = 16e3_r8 * log10(1013.25_r8 / PMEANW(I,J))
do L = 2, LPAR + 1
POFLE(L) = ETAA(L) + ETAB(L)*PSRF
!DELZ = -29.36_r8 * TW(I,J,L-1) * log(POFLE(L)/POFLE(L-1))
!// Thickness of layer below
DELZ = -29.26586_r8 * TW(I,J,L-1)*(1._r8 + 0.6_r8*QW(I,J,L-1)) &
* log(POFLE(L)/POFLE(L-1))
ZOFLEW(L,I,J) = ZOFLEW(L-1,I,J) + DELZ
end do
end do
end do
if (LPARW .eq. 40) then
if (MYEAR .lt. 2001) then
NrOf3Dfields = 8
else
! Splitting rain in large scale and convective
NrOf3Dfields = 9
end if
else if (LPARW .eq. 60) then
NrOf3Dfields = 9
!// Note that ERA-40 has 8 3D fields
end if
do NF = 1, NrOf3Dfields
read(14) FLD_ID, ISEC
!// Read native resolution
GRDATA3HI(:,:,:) = 0._r4
read(14) (((GRDATA3HI(I,J,L),I=1,IPARW),J=1,JPARW),L=ISEC(13),ISEC(14))
!// Degrade or not
GRDATA3(:,:,:) = 0._r4
if (ldeg) then
call TRUNG4(GRDATA3HI, GRDATA3, ZDEGI, ZDEGJ, IMAP, JMAP, IDGRD, &
JDGRD,IPARW,JPARW,IPAR,JPAR,LPARW,1)
else
do L = 1, LPAR
do J = 1, JPAR
do I = 1, IPAR
GRDATA3(I,J,L) = GRDATA3HI(I,J,L)
end do
end do
end do
end if
select case(FLD_ID)
case (212)
!// Mass Flux updrafts [accumulated kg/(m^2*s)]
!// - Mass flux is through BOTTOM EDGE of grid box (i.e. model
!// half layers).
!// - EC-data are stored from half layer 2, since the flux into
!// layer 1 is always zero (no flux in through surface).
!// - Scale with area and dt --> [kg/s]
CWETE(:,:,:) = 0._r8 !// No flux up from surface
!// Start retrieval at layer 2
do L = 2, LWEPAR
!// LMAP is full-level, but for since mass fluxes are not stored
!// for surface, and starts at layer 2, we need to subtract 1.
!// What comes into L=2, comes from LMAP(2)-1.
!// For L40, L=2 and LMAP(2)-1 = 1, while
!// for L37, L=2 and LMAP(2)-1 = 3
LL = LMAP(L) - 1
do J = 1, JPAR
do I = 1, IPAR
!// Filter values. File data has minimum values less than zero.
!// Hence we treat values less that that as as zero
if (GRDATA3(I,J,LL) .gt. mincwetelim) then
CWETE(I,J,L) = GRDATA3(I,J,LL) * AREAXY(I,J) * ZDT
else
CWETE(I,J,L) = 0._r8
end if
end do
end do
end do
case (213)
!// Mass Flux downdrafts [accumulated kg/(m^2*s)]
!// - Mass flux is through BOTTOM EDGE of grid box (i.e. model
!// half layers).
!// This means that the flux for layer 1 is zero, and that
!// downward flux is NEGATIVE, going out at the bottom of the box.
!// - EC-data are stored from half layer 2, since the flux into
!// layer 1 is always zero (no flux in through surface).
!// - For a grid box in layer L, -FD(L) goes out at bottom
!// and -FD(L+1) comes in from above (FD is negative). The
!// balance with entrainment E and detrainment D is:
!// FD(L) - FD(L+1) = E - D
!// - Scale with area and dt --> [kg/s]
CWETD(:,:,:) = 0._r8 !// No flux down to surface
!// Start retrieval at layer 2
do L = 2, LWDPAR
!// See Case(212) for comment on LMAP
LL = LMAP(L) - 1
do J = 1, JPAR
do I = 1, IPAR
if (GRDATA3(I,J,LL) .lt. maxcwetdlim) then
CWETD(I,J,L) = GRDATA3(I,J,LL) * AREAXY(I,J) * ZDT
else
CWETD(I,J,L) = 0._r8
end if
end do
end do
end do
case (214)
!// Updrafts detrainment rate [accumulated kg/(m3*s); kg/(m2*s)
!// per gridbox height]
!// - Entrainment has to be built from detrainment.
!// - The rate is per height, so we have to multiply with dZ and area.
!// - Detrainment/entrainment rates are given at GRID CENTER
!// (i.e. model full layers).
!// - Detrainment must be summed up when collapsing layers!
!// - Scale with box height, area and dt --> [kg/s]
CDETU(:,:,:) = 0._r8
do J = 1, JPAR
do I = 1, IPAR
do L = 1, LWEPAR
do LL = LMAP(L), LMAP(L+1) - 1
!// Detrainment must be summed up when collapsing layers
if (GRDATA3(I,J,LL) .gt. mincdetulim) then
CDETU(I,J,L) = CDETU(I,J,L) &
+ GRDATA3(I,J,LL) * XLMMAP(LL) &
* (ZOFLE(L+1,I,J) - ZOFLE(L,I,J)) &
* AREAXY(I,J) * ZDT
end if
end do
end do
end do
end Do
case (215)
!// Downdrafts detrainment rate [acc. kg/(m3*s); kg/(m2*s) per
!// gridbox height]
!// - Entrainment has to be built from detrainment.
!// - The rate is per height, so we have to multiply with dZ and
!// area.
!// - Detrainment/entrainment rates are given at GRID CENTER
!// (i.e. model full layers).
!// - Detrainment must be summed up when collapsing layers!
!// - Scale with box height and dt --> [kg/s]
CDETD(:,:,:) = 0._r8
do J = 1, JPAR
do I = 1, IPAR
do L = 1, LWDPAR
do LL = LMAP(L),LMAP(L+1)-1
if (GRDATA3(I,J,LL) .gt. mincdetdlim) then
CDETD(I,J,L) = CDETD(I,J,L) &
+ GRDATA3(I,J,LL) * XLMMAP(LL) &
* (ZOFLE(L+1,I,J) - ZOFLE(L,I,J)) &
* AREAXY(I,J) * ZDT
end if
end do
end do
end do
end do
case (216)
!/ Rainfall [kg/(m^2)] (accumulated kg/(m^2*s)) - edge value
!// scale with area and dt --> [kg/s]
do L = 1, LWEPAR
LL = LMAP(L)
do J = 1, JPAR
do I = 1, IPAR
PRECLS(I,J,L) = max(0._r8, GRDATA3(I,J,LL) * AREAXY(I,J) * ZDT)
end do
end do
end do
found216 = .TRUE.
case (217)
!// CONVECTIVE RAINFALL [kg/(m^2)] (accumulated kg/(m^2*s))
!// EDGE VALUE
!// Scale with area and dt or re-initialize --> [kg/s]
TMPCRAIN(:,:,:) = 0._r8
do L = 1, LWEPAR
LL = LMAP(L)
do J = 1, JPAR
do I = 1, IPAR
TMPCRAIN(I,J,L) = max(0._r8, GRDATA3(I,J,LL) * AREAXY(I,J) * ZDT)
end do
end do
end do
PRECCNV(:,:,:) = TMPCRAIN(:,:,:)
case (218)
!// LARGE SCALE RAINFALL 2 [kg/(m^2)] (accumulated kg/(m^2*s))
!// EDGE VALUE
!// Scale with area and dt or re-initialize --> [kg/s]
TMPLRAIN(:,:,:) = 0._r8
do L = 1, LWEPAR
LL = LMAP(L)
do J = 1, JPAR
do I = 1, IPAR
TMPLRAIN(I,J,L) = max(0._r8, GRDATA3(I,J,LL) * AREAXY(I,J) * ZDT)
end do
end do
end do
case (246)
!// Cloud Liquid Water Content [Kg/Kg] - midpoint value
!// Cloud routine needs CLDLWCW in any case, so we calculate
!// that first. Could repeat with GRDATA3, but rather do
!// degradation of CLDLWCW.
CLDLWCW(:,:,:) = 0._r8
do L = 1, LWEPAR
do LL = LMAP(L), LMAP(L+1)-1
do J = 1, JPARW
do I = 1, IPARW
CLDLWCW(I,J,L) = CLDLWCW(I,J,L) &
+ GRDATA3HI(I,J,LL) * XLMMAP(LL)
end do
end do
end do
end do
if (ldeg) then
call TRUNG8(CLDLWCW, CLDLWC, ZDEGI, ZDEGJ, IMAP, JMAP, IDGRD, &
JDGRD, IPARW, JPARW, IPAR, JPAR, LWEPAR, 1)
else
CLDLWC(:,:,:) = 0._r8
do L = 1, LWEPAR
do J = 1,JPAR
do I = 1,IPAR
CLDLWC(I,J,L) = CLDLWCW(I,J,L)
end do
end do
end do
end if
case (247)
!// Cloud Ice Water Content [Kg/Kg] - midpoint value
!// Cloud routine needs CLDIWCW in any case, so we calculate
!// that first.
CLDIWCW(:,:,:) = 0._r8
do L = 1, LWEPAR
do LL = LMAP(L), LMAP(L+1)-1
do J = 1, JPARW
do I = 1, IPARW
CLDIWCW(I,J,L) = CLDIWCW(I,J,L) &
+ GRDATA3HI(I,J,LL) * XLMMAP(LL)
end do
end do
end do
end do
if (ldeg) then
call TRUNG8(CLDIWCW, CLDIWC, ZDEGI, ZDEGJ, IMAP, JMAP, IDGRD, &
JDGRD, IPARW, JPARW, IPAR, JPAR, LWEPAR, 1)
else
CLDIWC(:,:,:) = 0._r8
do L = 1, LWEPAR
do J = 1, JPAR
do I = 1, IPAR
CLDIWC(I,J,L) = CLDIWCW(I,J,L)
end do
end do
end do
end if
case (248)
!// Cloud Fraction [0,1] - midpoint value
!// Cloud routine needs CLDFRW in any case, so we calculate
!// that first.
CLDFRW(:,:,:) = 0._r8
do L = 1, LWEPAR
do LL = LMAP(L), LMAP(L+1)-1
do J = 1, JPARW
do I = 1, IPARW
CLDFRW(I,J,L) = CLDFRW(I,J,L) &
+ XLMMAP(LL) * max(min(GRDATA3HI(I,J,LL),1._r4),0._r4)
end do
end do
end do
end do
if (ldeg) then
call TRUNG8(CLDFRW, CLDFR, ZDEGI, ZDEGJ, IMAP, JMAP, IDGRD, &
JDGRD, IPARW, JPARW, IPAR, JPAR, LWEPAR, 1)
else
do L = 1, LWEPAR
do J = 1, JPAR
do I = 1, IPAR
CLDFR(I,J,L) = CLDFRW(I,J,L)
end do
end do
end do
end if
!// Cloud cover routines