forked from CFMIP/COSPv1
-
Notifications
You must be signed in to change notification settings - Fork 0
/
cosp_io.F90
1980 lines (1889 loc) · 80.6 KB
/
cosp_io.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
! (c) British Crown Copyright 2008, the Met Office.
! All rights reserved.
! $Revision: 88 $, $Date: 2013-11-13 07:08:38 -0700 (Wed, 13 Nov 2013) $
! $URL: http://cfmip-obs-sim.googlecode.com/svn/stable/v1.4.0/cosp_io.F90 $
!
! Redistribution and use in source and binary forms, with or without modification, are permitted
! provided that the following conditions are met:
!
! * Redistributions of source code must retain the above copyright notice, this list
! of conditions and the following disclaimer.
! * Redistributions in binary form must reproduce the above copyright notice, this list
! of conditions and the following disclaimer in the documentation and/or other materials
! provided with the distribution.
! * Neither the name of the Met Office nor the names of its contributors may be used
! to endorse or promote products derived from this software without specific prior written
! permission.
!
! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR
! IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND
! FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR
! CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL
! DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER
! IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
! OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
!
! History:
! Jul 2008 - A. Bodas-Salcedo - Initial version
! Oct 2008 - S. Bony - In nc_write_cosp_1d and nc_write_cosp_2d :
! the label of layered cloud fractions was wrong -> corrected
! (before: low was actually mid, mid was high, high was total,
! total was low)
! Sep 2009 - A. Bodas-Salcedo - CMIP5 variable names implemented
!
! Jan 2013 - G. Cesana - Add new phase variables for outputs
#include "cosp_defs.h"
MODULE MOD_COSP_IO
USE MOD_COSP_CONSTANTS
USE MOD_COSP_TYPES
USE cmor_users_functions
USE netcdf
use MOD_COSP_Modis_Simulator
IMPLICIT NONE
! Types to be used as arrays of pointers
TYPE var1d
character(len=16) :: name
character(len=16) :: units
integer :: dimsid(3)
integer :: dimssz(2)
integer :: vid
logical :: lout
real,pointer,dimension(:) :: pntr
END TYPE
TYPE var2d
character(len=16) :: name
character(len=16) :: units
integer :: dimsid(4)
integer :: dimssz(3)
integer :: vid
logical :: lout
real,pointer,dimension(:,:) :: pntr
END TYPE
TYPE var3d
character(len=16) :: name
character(len=16) :: units
integer :: dimsid(5)
integer :: dimssz(4)
integer :: vid
logical :: lout
real,pointer,dimension(:,:,:) :: pntr
END TYPE
CONTAINS
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
!----------------- SUBROUTINE CONSTRUCT_VAR1D --------------------
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
SUBROUTINE CONSTRUCT_VAR1D(name,dimsid,dimssz,pntr,y,units)
! Input arguments
character(len=*),intent(in) :: name
integer,intent(in) :: dimsid(3)
integer,intent(in) :: dimssz(2)
real,dimension(:),target,intent(in) :: pntr
type(var1d),intent(out) :: y
character(len=*),optional,intent(in) :: units
y%name = name
if (present(units)) y%units = units
y%dimsid = dimsid
y%dimssz = dimssz
y%pntr => pntr
END SUBROUTINE CONSTRUCT_VAR1D
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
!----------------- SUBROUTINE CONSTRUCT_VAR2D --------------------
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
SUBROUTINE CONSTRUCT_VAR2D(name,dimsid,dimssz,pntr,y,units)
! Input arguments
character(len=*),intent(in) :: name
integer,intent(in) :: dimsid(4)
integer,intent(in) :: dimssz(3)
real,dimension(:,:),target,intent(in) :: pntr
type(var2d),intent(out) :: y
character(len=*),optional,intent(in) :: units
y%name = name
if (present(units)) y%units = units
y%dimsid = dimsid
y%dimssz = dimssz
y%pntr => pntr
END SUBROUTINE CONSTRUCT_VAR2D
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
!----------------- SUBROUTINE CONSTRUCT_VAR3D --------------------
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
SUBROUTINE CONSTRUCT_VAR3D(name,dimsid,dimssz,pntr,y,units)
! Input arguments
character(len=*),intent(in) :: name
integer,intent(in) :: dimsid(5)
integer,intent(in) :: dimssz(4)
real,dimension(:,:,:),target,intent(in) :: pntr
type(var3d),intent(out) :: y
character(len=*),optional,intent(in) :: units
y%name = name
if (present(units)) y%units = units
y%dimsid = dimsid
y%dimssz = dimssz
y%pntr => pntr
END SUBROUTINE CONSTRUCT_VAR3D
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
!----------------- SUBROUTINE MAP_POINT_TO_LL---------------------
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
SUBROUTINE MAP_POINT_TO_LL(Nx,Ny,geomode,x1,x2,x3,x4,y2,y3,y4,y5)
! Input arguments
integer,intent(in) :: Nx,Ny,geomode
real,intent(in),optional :: x1(:),x2(:,:),x3(:,:,:), &
x4(:,:,:,:)
real,intent(out),optional :: y2(:,:),y3(:,:,:), &
y4(:,:,:,:),y5(:,:,:,:,:)
! Local variables
integer :: Npoints
integer :: px(Nx*Ny),py(Nx*Ny)
integer :: i,j,k,l,m
integer :: Ni,Nj,Nk,Nl
integer :: Mi,Mj,Mk,Ml,Mm
character(len=128) :: proname='MAP_POINT_TO_LL'
Npoints = Nx*Ny
px=0
py=0
! Obtain pointers to do the mapping
if (geomode == 2) then ! (lon,lat) mode
do j=1,Ny
do i=1,Nx
k = (j-1)*Nx+i
px(k) = i
py(k) = j
enddo
enddo
else if (geomode == 3) then ! (lon,lat) mode
do j=1,Nx
do i=1,Ny
k = (j-1)*Ny+i
px(k) = j
py(k) = i
enddo
enddo
else
print *, ' -- '//trim(proname)//': geomode not supported, ',geomode
stop
endif
if (present(x1).and.present(y2)) then
Ni = size(x1,1)
Mi = size(y2,1)
Mj = size(y2,2)
if (Mi*Mj /= Ni) then
print *, ' -- '//trim(proname)//': Nlon*Nlat /= Npoints (opt 1)'
stop
endif
do i=1,Npoints
y2(px(i),py(i)) = x1(i)
enddo
else if (present(x2).and.present(y3)) then
Ni = size(x2,1)
Nj = size(x2,2)
Mi = size(y3,1)
Mj = size(y3,2)
Mk = size(y3,3)
if (Mi*Mj /= Ni) then
print *, ' -- '//trim(proname)//': Nlon*Nlat /= Npoints (opt 2)'
stop
endif
if (Nj /= Mk) then
print *, ' -- '//trim(proname)//': Nj /= Mk (opt 2)'
print *, Ni,Nj,Mi,Mj,Mk
stop
endif
do k=1,Mk
do i=1,Npoints
y3(px(i),py(i),k) = x2(i,k)
enddo
enddo
else if (present(x3).and.present(y4)) then
Ni = size(x3,1)
Nj = size(x3,2)
Nk = size(x3,3)
Mi = size(y4,1)
Mj = size(y4,2)
Mk = size(y4,3)
Ml = size(y4,4)
if (Mi*Mj /= Ni) then
print *, ' -- '//trim(proname)//': Nlon*Nlat /= Npoints (opt 3)'
stop
endif
if (Nj /= Mk) then
print *, ' -- '//trim(proname)//': Nj /= Mk (opt 3)'
stop
endif
if (Nk /= Ml) then
print *, ' -- '//trim(proname)//': Nk /= Ml (opt 3)'
stop
endif
do l=1,Ml
do k=1,Mk
do i=1,Npoints
y4(px(i),py(i),k,l) = x3(i,k,l)
enddo
enddo
enddo
else if (present(x4).and.present(y5)) then
Ni = size(x4,1)
Nj = size(x4,2)
Nk = size(x4,3)
Nl = size(x4,4)
Mi = size(y5,1)
Mj = size(y5,2)
Mk = size(y5,3)
Ml = size(y5,4)
Mm = size(y5,5)
if (Mi*Mj /= Ni) then
print *, ' -- '//trim(proname)//': Nlon*Nlat /= Npoints (opt 4)'
stop
endif
if (Nj /= Mk) then
print *, ' -- '//trim(proname)//': Nj /= Mk (opt 4)'
stop
endif
if (Nk /= Ml) then
print *, ' -- '//trim(proname)//': Nk /= Ml (opt 4)'
stop
endif
if (Nl /= Mm) then
print *, ' -- '//trim(proname)//': Nl /= Mm (opt 4)'
stop
endif
do m=1,Mm
do l=1,Ml
do k=1,Mk
do i=1,Npoints
y5(px(i),py(i),k,l,m) = x4(i,k,l,m)
enddo
enddo
enddo
enddo
else
print *, ' -- '//trim(proname)//': wrong option'
stop
endif
END SUBROUTINE MAP_POINT_TO_LL
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
!----------------- SUBROUTINE MAP_LL_TO_POINT---------------------
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
SUBROUTINE MAP_LL_TO_POINT(Nx,Ny,Np,x2,x3,x4,x5,y1,y2,y3,y4)
! Input arguments
integer,intent(in) :: Nx,Ny,Np
real,intent(in),optional :: x2(:,:),x3(:,:,:), &
x4(:,:,:,:),x5(:,:,:,:,:)
real,intent(out),optional :: y1(:),y2(:,:),y3(:,:,:), &
y4(:,:,:,:)
! Local variables
integer :: px(Nx*Ny),py(Nx*Ny)
integer :: i,j,k,l,m
integer :: Ni,Nj,Nk,Nl,Nm
integer :: Mi,Mj,Mk,Ml
character(len=128) :: proname='MAP_LL_TO_POINT'
px=0
py=0
if (Nx*Ny < Np) then
print *, ' -- '//trim(proname)//': Nx*Ny < Np'
stop
endif
do j=1,Ny
do i=1,Nx
k = (j-1)*Nx+i
px(k) = i
py(k) = j
enddo
enddo
if (present(x2).and.present(y1)) then
Ni = size(x2,1)
Nj = size(x2,2)
Mi = size(y1,1)
if (Ni*Nj < Mi) then
print *, ' -- '//trim(proname)//': Nlon*Nlat < Npoints (opt 1)'
stop
endif
do j=1,Np
y1(j) = x2(px(j),py(j))
enddo
else if (present(x3).and.present(y2)) then
Ni = size(x3,1)
Nj = size(x3,2)
Nk = size(x3,3)
Mi = size(y2,1)
Mj = size(y2,2)
if (Ni*Nj < Mi) then
print *, ' -- '//trim(proname)//': Nlon*Nlat < Npoints (opt 2)'
stop
endif
if (Nk /= Mj) then
print *, ' -- '//trim(proname)//': Nk /= Mj (opt 2)'
stop
endif
do k=1,Nk
do j=1,Np
y2(j,k) = x3(px(j),py(j),k)
enddo
enddo
else if (present(x4).and.present(y3)) then
Ni = size(x4,1)
Nj = size(x4,2)
Nk = size(x4,3)
Nl = size(x4,4)
Mi = size(y3,1)
Mj = size(y3,2)
Mk = size(y3,3)
if (Ni*Nj < Mi) then
print *, ' -- '//trim(proname)//': Nlon*Nlat < Npoints (opt 3)'
stop
endif
if (Nk /= Mj) then
print *, ' -- '//trim(proname)//': Nk /= Mj (opt 3)'
stop
endif
if (Nl /= Mk) then
print *, ' -- '//trim(proname)//': Nl /= Mk (opt 3)'
stop
endif
do l=1,Nl
do k=1,Nk
do j=1,Np
y3(j,k,l) = x4(px(j),py(j),k,l)
enddo
enddo
enddo
else if (present(x5).and.present(y4)) then
Ni = size(x5,1)
Nj = size(x5,2)
Nk = size(x5,3)
Nl = size(x5,4)
Nm = size(x5,5)
Mi = size(y4,1)
Mj = size(y4,2)
Mk = size(y4,3)
Ml = size(y4,4)
if (Ni*Nj < Mi) then
print *, ' -- '//trim(proname)//': Nlon*Nlat < Npoints (opt 4)'
stop
endif
if (Nk /= Mj) then
print *, ' -- '//trim(proname)//': Nk /= Mj (opt 4)'
stop
endif
if (Nl /= Mk) then
print *, ' -- '//trim(proname)//': Nl /= Mk (opt 4)'
stop
endif
if (Nm /= Ml) then
print *, ' -- '//trim(proname)//': Nm /= Ml (opt 4)'
stop
endif
do m=1,Nm
do l=1,Nl
do k=1,Nk
do j=1,Np
y4(j,k,l,m) = x5(px(j),py(j),k,l,m)
enddo
enddo
enddo
enddo
else
print *, ' -- '//trim(proname)//': wrong option'
stop
endif
END SUBROUTINE MAP_LL_TO_POINT
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
!----------------- SUBROUTINE NC_READ_INPUT_FILE -----------------
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
SUBROUTINE NC_READ_INPUT_FILE(fname,Npnts,Nl,Nhydro,lon,lat,p,ph,z,zh,T,qv,rh,tca,cca, &
mr_lsliq,mr_lsice,mr_ccliq,mr_ccice,fl_lsrain,fl_lssnow,fl_lsgrpl, &
fl_ccrain,fl_ccsnow,Reff,dtau_s,dtau_c,dem_s,dem_c,skt,landmask, &
mr_ozone,u_wind,v_wind,sunlit,emsfc_lw,mode,Nlon,Nlat,verbosity)
!Arguments
character(len=512),intent(in) :: fname ! File name
integer,intent(in) :: Npnts,Nl,Nhydro
real,dimension(Npnts),intent(out) :: lon,lat
real,dimension(Npnts,Nl),target,intent(out) :: p,ph,z,zh,T,qv,rh,tca,cca, &
mr_lsliq,mr_lsice,mr_ccliq,mr_ccice,fl_lsrain,fl_lssnow,fl_lsgrpl, &
fl_ccrain,fl_ccsnow,dtau_s,dtau_c,dem_s,dem_c,mr_ozone
real,dimension(Npnts,Nl,Nhydro),intent(out) :: Reff
real,dimension(Npnts),intent(out) :: skt,landmask,u_wind,v_wind,sunlit
real,intent(out) :: emsfc_lw
integer,intent(out) :: mode,Nlon,Nlat
integer,optional :: verbosity
!Local variables
integer :: Npoints,Nlevels,i,j,k
character(len=128) :: vname
integer,parameter :: NMAX_DIM=5
integer :: vrank,vdimid(NMAX_DIM)
character(len=256) :: dimname(NMAX_DIM) ! 256 hardcoded, instead of MAXNCNAM. This works for NetCDF 3 and 4.
integer :: ncid,vid,ndims,nvars,ngatts,recdim,dimsize(NMAX_DIM)
integer :: errst
logical :: Llat,Llon,Lpoint
integer :: Na,Nb,Nc,Nd,Ne
real,dimension(Npnts) :: ll
integer,dimension(:),allocatable :: plon,plat
real,allocatable :: x1(:),x2(:,:),x3(:,:,:),x4(:,:,:,:),x5(:,:,:,:,:) ! Temporary arrays
character(len=64) :: routine_name='NC_READ_INPUT_FILE'
character(len=128) :: errmsg,straux
mode = 0
Nlon = 0
Nlat = 0
Npoints = Npnts
Nlevels = Nl
! Open file
errst = nf90_open(fname, nf90_nowrite, ncid)
if (errst /= 0) then
errmsg="Couldn't open "//trim(fname)
call cosp_error(routine_name,errmsg)
endif
! Get information about dimensions. Curtain mode or lat/lon mode?
Llat =.false.
Llon =.false.
Lpoint=.false.
errst = nf90_inquire(ncid, ndims, nvars, ngatts, recdim)
if (errst /= 0) then
errmsg="Error in nf90_inquire"
call cosp_error(routine_name,errmsg,errcode=errst)
endif
do i = 1,ndims
errst = nf90_Inquire_Dimension(ncid,i,name=dimname(i),len=dimsize(i))
if (errst /= 0) then
write(straux, *) i
errmsg="Error in nf90_Inquire_Dimension, i: "//trim(straux)
call cosp_error(routine_name,errmsg)
endif
if ((trim(dimname(i)).eq.'level').and.(Nlevels > dimsize(i))) then
errmsg='Number of levels selected is greater than in input file '//trim(fname)
call cosp_error(routine_name,errmsg)
endif
if (trim(dimname(i)).eq.'point') then
Lpoint = .true.
if (Npnts > dimsize(i)) then
errmsg='Number of points selected is greater than in input file '//trim(fname)
call cosp_error(routine_name,errmsg)
endif
endif
if (trim(dimname(i)).eq.'lon') then
Llon = .true.
Nlon = dimsize(i)
endif
if (trim(dimname(i)).eq.'lat') then
Llat = .true.
Nlat = dimsize(i)
endif
enddo
! Get lon and lat
if (Llon.and.Llat) then ! 2D mode
if ((Npnts) > Nlon*Nlat) Npoints=Nlon*Nlat
lon = R_UNDEF
lat = R_UNDEF
mode = 2 ! Don't know yet if (lon,lat) or (lat,lon) at this point
else if (Lpoint) then ! 1D mode
Nlon = Npoints
Nlat = Npoints
mode = 1
else
errmsg= trim(fname)//' file contains wrong dimensions'
call cosp_error(routine_name,errmsg)
endif
errst = nf90_inq_varid(ncid, 'lon', vid)
if (errst /= 0) then
errmsg="Error in nf90_inq_varid, var: lon"
call cosp_error(routine_name,errmsg,errcode=errst)
endif
errst = nf90_get_var(ncid, vid, lon, start = (/1/), count = (/Nlon/))
if (errst /= 0) then
errmsg="Error in nf90_get_var, var: lon"
call cosp_error(routine_name,errmsg,errcode=errst)
endif
errst = nf90_inq_varid(ncid, 'lat', vid)
if (errst /= 0) then
errmsg="Error in nf90_inq_varid, var: lat"
call cosp_error(routine_name,errmsg,errcode=errst)
endif
errst = nf90_get_var(ncid, vid, lat, start = (/1/), count = (/Nlat/))
if (errst /= 0) then
errmsg="Error in nf90_get_var, var: lat"
call cosp_error(routine_name,errmsg,errcode=errst)
endif
! Get all variables
do vid = 1,nvars
vdimid=0
errst = nf90_Inquire_Variable(ncid, vid, name=vname, ndims=vrank, dimids=vdimid)
if (errst /= 0) then
write(straux, *) vid
errmsg='Error in nf90_Inquire_Variable, vid '//trim(straux)
call cosp_error(routine_name,errmsg,errcode=errst)
endif
! Read in into temporary array of correct shape
!if (present(verbosity).and.(verbosity == 1)) print *, 'Reading '//trim(vname)//' ...'
if (vrank == 1) then
Na = dimsize(vdimid(1))
allocate(x1(Na))
errst = nf90_get_var(ncid, vid, x1, start=(/1/), count=(/Na/))
endif
if (vrank == 2) then
Na = dimsize(vdimid(1))
Nb = dimsize(vdimid(2))
allocate(x2(Na,Nb))
errst = nf90_get_var(ncid, vid, x2, start=(/1,1/), count=(/Na,Nb/))
endif
if (vrank == 3) then
Na = dimsize(vdimid(1))
Nb = dimsize(vdimid(2))
Nc = dimsize(vdimid(3))
allocate(x3(Na,Nb,Nc))
errst = nf90_get_var(ncid, vid, x3, start=(/1,1,1/), count=(/Na,Nb,Nc/))
if ((mode == 2).or.(mode == 3)) then
if ((Na == Nlon).and.(Nb == Nlat)) then
mode = 2
else if ((Na == Nlat).and.(Nb == Nlon)) then
mode = 3
else
errmsg='Wrong mode for variable '//trim(vname)
call cosp_error(routine_name,errmsg)
endif
endif
endif
if (vrank == 4) then
Na = dimsize(vdimid(1))
Nb = dimsize(vdimid(2))
Nc = dimsize(vdimid(3))
Nd = dimsize(vdimid(4))
allocate(x4(Na,Nb,Nc,Nd))
errst = nf90_get_var(ncid, vid, x4, start=(/1,1,1,1/), count=(/Na,Nb,Nc,Nd/))
endif
if (vrank == 5) then
Na = dimsize(vdimid(1))
Nb = dimsize(vdimid(2))
Nc = dimsize(vdimid(3))
Nd = dimsize(vdimid(4))
Ne = dimsize(vdimid(5))
allocate(x5(Na,Nb,Nc,Nd,Ne))
errst = nf90_get_var(ncid, vid, x5, start=(/1,1,1,1,1/), count=(/Na,Nb,Nc,Nd,Ne/))
endif
if (errst /= 0) then
write(straux, *) vid
errmsg='Error in nf90_get_var, vid '//trim(straux)
call cosp_error(routine_name,errmsg,errcode=errst)
endif
! Map to the right input argument
select case (trim(vname))
case ('pfull')
if (Lpoint) then
p(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=p)
endif
case ('phalf')
if (Lpoint) then
ph(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=ph)
endif
case ('height')
if (Lpoint) then
z(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=z)
endif
case ('height_half')
if (Lpoint) then
zh(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=zh)
endif
case ('T_abs')
if (Lpoint) then
T(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=T)
endif
case ('qv')
if (Lpoint) then
qv(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=qv)
endif
case ('rh')
if (Lpoint) then
rh(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=rh)
endif
case ('tca')
if (Lpoint) then
tca(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=tca)
endif
tca = tca
case ('cca')
if (Lpoint) then
cca(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=cca)
endif
cca = cca
case ('mr_lsliq')
if (Lpoint) then
mr_lsliq(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=mr_lsliq)
endif
case ('mr_lsice')
if (Lpoint) then
mr_lsice(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=mr_lsice)
endif
case ('mr_ccliq')
if (Lpoint) then
mr_ccliq(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=mr_ccliq)
endif
case ('mr_ccice')
if (Lpoint) then
mr_ccice(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=mr_ccice)
endif
case ('fl_lsrain')
if (Lpoint) then
fl_lsrain(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=fl_lsrain)
endif
case ('fl_lssnow')
if (Lpoint) then
fl_lssnow(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=fl_lssnow)
endif
case ('fl_lsgrpl')
if (Lpoint) then
fl_lsgrpl(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=fl_lsgrpl)
endif
case ('fl_ccrain')
if (Lpoint) then
fl_ccrain(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=fl_ccrain)
endif
case ('fl_ccsnow')
if (Lpoint) then
fl_ccsnow(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=fl_ccsnow)
endif
case ('dtau_s')
if (Lpoint) then
dtau_s(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=dtau_s)
endif
case ('dtau_c')
if (Lpoint) then
dtau_c(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=dtau_c)
endif
case ('dem_s')
if (Lpoint) then
dem_s(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=dem_s)
endif
case ('dem_c')
if (Lpoint) then
dem_c(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=dem_c)
endif
case ('Reff')
if (Lpoint) then
Reff(1:Npoints,:,:) = x3(1:Npoints,1:Nlevels,:)
else
call map_ll_to_point(Na,Nb,Npoints,x4=x4,y3=Reff)
endif
case ('skt')
if (Lpoint) then
skt(1:Npoints) = x1(1:Npoints)
else
call map_ll_to_point(Na,Nb,Npoints,x2=x2,y1=skt)
endif
case ('landmask')
if (Lpoint) then
landmask(1:Npoints) = x1(1:Npoints)
else
call map_ll_to_point(Na,Nb,Npoints,x2=x2,y1=landmask)
endif
case ('mr_ozone')
if (Lpoint) then
mr_ozone(1:Npoints,:) = x2(1:Npoints,1:Nlevels)
else
call map_ll_to_point(Na,Nb,Npoints,x3=x3,y2=mr_ozone)
endif
case ('u_wind')
if (Lpoint) then
u_wind(1:Npoints) = x1(1:Npoints)
else
call map_ll_to_point(Na,Nb,Npoints,x2=x2,y1=u_wind)
endif
case ('v_wind')
if (Lpoint) then
v_wind(1:Npoints) = x1(1:Npoints)
else
call map_ll_to_point(Na,Nb,Npoints,x2=x2,y1=v_wind)
endif
case ('sunlit')
if (Lpoint) then
sunlit(1:Npoints) = x1(1:Npoints)
else
call map_ll_to_point(Na,Nb,Npoints,x2=x2,y1=sunlit)
endif
end select
! Free memory
if (vrank == 1) deallocate(x1)
if (vrank == 2) deallocate(x2)
if (vrank == 3) deallocate(x3)
if (vrank == 4) deallocate(x4)
if (vrank == 5) deallocate(x5)
enddo
! SFC emissivity
errst = nf90_inq_varid(ncid, 'emsfc_lw', vid)
if (errst /= 0) then
if (errst == nf90_enotvar) then ! Does not exist, use 1.0
emsfc_lw = 1.0
print *, ' ********* COSP Warning: emsfc_lw does not exist in input file. Set to 1.0.'
else ! Other error, stop
errmsg='Error in nf90_inq_varid, var: emsfc_lw'
call cosp_error(routine_name,errmsg,errcode=errst)
endif
else
errst = nf90_get_var(ncid, vid, emsfc_lw)
if (errst /= 0) then
errmsg='Error in nf90_get_var, var: emsfc_lw'
call cosp_error(routine_name,errmsg,errcode=errst)
endif
endif
! Fill in the lat/lon vectors with the right values for 2D modes
! This might be helpful if the inputs are 2D (gridded) and
! you want outputs in 1D mode
allocate(plon(Npoints),plat(Npoints))
if (mode == 2) then !(lon,lat)
ll = lat
do j=1,Nb
do i=1,Na
k = (j-1)*Na + i
plon(k) = i
plat(k) = j
enddo
enddo
lon(1:Npoints) = lon(plon(1:Npoints))
lat(1:Npoints) = ll(plat(1:Npoints))
else if (mode == 3) then !(lat,lon)
ll = lon
do j=1,Nb
do i=1,Na
k = (j-1)*Na + i
lon(k) = ll(j)
lat(k) = lat(i)
enddo
enddo
lon(1:Npoints) = ll(plon(1:Npoints))
lat(1:Npoints) = lat(plat(1:Npoints))
endif
deallocate(plon,plat)
! Close file
errst = nf90_close(ncid)
if (errst /= 0) then
errmsg='Error in nf90_close'
call cosp_error(routine_name,errmsg,errcode=errst)
endif
END SUBROUTINE NC_READ_INPUT_FILE
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
!--------------- SUBROUTINE NC_CMOR_INIT ---------------------
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
SUBROUTINE NC_CMOR_INIT(cmor_nl,wmode,cfg,vgrid,gb,sg,sgradar,sglidar,&
isccp,misr,modis,rttov,stradar,stlidar,geomode,Nlon,Nlat,N1,N2,N3, &
lon_axid,lat_axid,time_axid,height_axid,height_mlev_axid,grid_id,lonvar_id,latvar_id, &
column_axid,sza_axid,temp_axid,channel_axid,dbze_axid, &
sratio_axid,MISR_CTH_axid,tau_axid,pressure2_axid,reffICE_axid,reffLIQ_axid, &
v1d,v2d,v3d)
! Input arguments
character(len=*),intent(in) :: cmor_nl
character(len=*),intent(in) :: wmode ! Writing mode 'replace' or 'append'
type(cosp_config),intent(in) :: cfg
type(cosp_vgrid),intent(in) :: vgrid
type(cosp_gridbox),intent(in) :: gb
type(cosp_subgrid),intent(in) :: sg
type(cosp_sgradar),intent(in) :: sgradar ! Subgrid radar
type(cosp_sglidar),intent(in) :: sglidar ! Subgrid lidar
type(cosp_isccp),intent(in) :: isccp ! ISCCP outputs
type(cosp_misr),intent(in) :: misr ! MISR outputs
type(cosp_modis),intent(in) :: modis ! MODIS outputs
type(cosp_rttov),intent(in) :: rttov ! RTTOV outputs
type(cosp_radarstats),intent(in) :: stradar ! Summary statistics from radar simulator
type(cosp_lidarstats),intent(in) :: stlidar ! Summary statistics from lidar simulator
integer,intent(in) :: geomode,Nlon,Nlat,N1,N2,N3
integer,intent(out) :: grid_id,latvar_id,lonvar_id,column_axid,height_axid,dbze_axid,height_mlev_axid,sratio_axid, &
tau_axid,pressure2_axid,lon_axid,lat_axid,time_axid,sza_axid,MISR_CTH_axid, &
channel_axid,temp_axid,reffICE_axid,reffLIQ_axid
type(var1d),intent(inout) :: v1d(N1)
type(var2d),intent(inout) :: v2d(N2)
type(var3d),intent(inout) :: v3d(N3)
!--- Local variables ---
integer :: profile_axid
integer :: error_flag,i,j,Npoints,Ncolumns,Nlevels,Nlvgrid,maxtsteps,Nchannels,Dmax
logical :: lfound
real :: lon_ax(Nlon),lat_ax(Nlat)
character(len=512) :: inpath,outpath,start_date,model_id,experiment_id,institution,institute_id,source,calendar, &
contact,history,comment,table,parent_experiment_id,parent_experiment_rip,forcing
character(len=2056) :: references
integer :: realization, nc_action,initialization_method,physics_version
double precision :: branch_time
namelist/CMOR/inpath,outpath,start_date,model_id,experiment_id,branch_time,parent_experiment_id,parent_experiment_rip, &
forcing,institution,institute_id,source,calendar,realization,initialization_method,physics_version, &
contact,history,comment,references,table,maxtsteps
real,dimension(:),allocatable :: profile_ax,column_ax,dbze_ax,channel_ax
real,dimension(:,:),allocatable :: dbze_bounds,vgrid_bounds,sratio_bounds, &
lon_bounds,lat_bounds,mgrid_bounds
integer :: d2(2),d3(3),d4(4),d5(5)
double precision :: tbnds(2,1)
character(len=64) :: pro_name = 'NC_CMOR_INIT'
Npoints = gb%Npoints
Ncolumns = gb%Ncolumns
Nlevels = gb%Nlevels
Nchannels = gb%Nchan
Nlvgrid = vgrid%Nlvgrid
!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
! Safety checks
!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
if (geomode > 1) then
if (Npoints > Nlon*Nlat) then
Npoints = Nlon*Nlat
print *, ' -- '//trim(pro_name)//' Warning: Npoints > Nlon*Nlat'
endif
if (Npoints < Nlon*Nlat) then
print *, ' -- '//trim(pro_name)//': only Npoints >= Nlon*Nlat is supported'
stop
endif
endif
nc_action = CMOR_APPEND_3
if (trim(wmode) == 'replace') nc_action = CMOR_REPLACE_3
!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
! Allocate memory and compute axes and bounds
!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
tbnds(:,1) = gb%time_bnds
allocate(column_ax(Ncolumns),dbze_ax(DBZE_BINS),channel_ax(Nchannels), &
dbze_bounds(2,DBZE_BINS),vgrid_bounds(2,Nlvgrid),mgrid_bounds(2,Nlevels),sratio_bounds(2,SR_BINS), &
lon_bounds(2,Nlon),lat_bounds(2,Nlat))
! Profile
if (geomode == 1) then
allocate(profile_ax(Npoints))
do i=1,Npoints
profile_ax(i) = i
enddo
endif
! Column
do i=1,gb%Ncolumns
column_ax(i) = i
enddo
! Channels
channel_ax = float(gb%ichan)
! Radar Ze
do i=1,DBZE_BINS
dbze_ax(i) = CFAD_ZE_MIN + CFAD_ZE_WIDTH*(i - 0.5)
enddo
do i=1,DBZE_BINS
dbze_bounds(1,i) = CFAD_ZE_MIN + CFAD_ZE_WIDTH*(i - 1)
dbze_bounds(2,i) = CFAD_ZE_MIN + CFAD_ZE_WIDTH*i
enddo
! Height of model levels
do i=1,Nlevels
mgrid_bounds(1,i) = vgrid%mzl(i)
mgrid_bounds(2,i) = vgrid%mzu(i)
enddo
! Height of std grid
do i=1,Nlvgrid
vgrid_bounds(1,i) = vgrid%zl(i)
vgrid_bounds(2,i) = vgrid%zu(i)
enddo
! Lidar scattering ratio bounds (They are output by cosp_cfad_sr->diag_lidar in lmd_ipsl_stats.f90)
sratio_bounds(2,:) = stlidar%srbval(:) ! srbval contains the upper limits from lmd_ipsl_stats.f90
sratio_bounds(1,2:SR_BINS) = stlidar%srbval(1:SR_BINS-1)
sratio_bounds(1,1) = 0.0
sratio_bounds(2,SR_BINS) = 1.e5 ! This matches with Chepfer et al., JGR, 2009. However, it is not consistent
! with the upper limit in lmd_ipsl_stats.f90, which is LIDAR_UNDEF-1=998.999
! Lat lon axes
if (geomode == 2) then
lon_ax = gb%longitude(1:Nlon)
lat_ax = gb%latitude(1:Npoints:Nlon)
else if (geomode == 3) then
lon_ax = gb%longitude(1:Npoints:Nlat)
lat_ax = gb%latitude(1:Nlat)
else if (geomode == 4) then
lon_ax = gb%longitude(1:Nlon)
lat_ax = gb%latitude(1:Nlat)
endif
if (geomode > 1) then
lon_bounds(1,2:Nlon) = (lon_ax(1:Nlon-1) + lon_ax(2:Nlon))/2.0
lon_bounds(1,1) = lon_ax(1) - (lon_bounds(1,2) - lon_ax(1))
lon_bounds(2,1:Nlon-1) = lon_bounds(1,2:Nlon)
lon_bounds(2,Nlon) = lon_ax(Nlon) + (lon_ax(Nlon) - lon_bounds(2,Nlon-1))
lat_bounds(1,2:Nlat) = (lat_ax(1:Nlat-1) + lat_ax(2:Nlat))/2.0
lat_bounds(1,1) = lat_ax(1) - (lat_bounds(1,2) - lat_ax(1))
lat_bounds(2,1:Nlat-1) = lat_bounds(1,2:Nlat)
lat_bounds(2,Nlat) = lat_ax(Nlat) + (lat_ax(Nlat) - lat_bounds(2,Nlat-1))
endif
!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
! Read namelist with information for CMOR output file
!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
open(10,file=cmor_nl,status='old')
read(10,nml=cmor)
close(10)
!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
! Specify path for tables and set up other CMOR options
!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
error_flag = cmor_setup(inpath=trim(inpath),netcdf_file_action=nc_action,create_subdirectories=0)
!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
! Define dataset as output from COSP, and other model details
!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
error_flag = cmor_dataset(outpath=trim(outpath),experiment_id=trim(experiment_id),institution=trim(institution), &
source=trim(source),calendar=trim(calendar),realization=realization,contact=trim(contact), &
history=trim(history),comment=trim(comment),references=trim(references),model_id=trim(model_id), &
branch_time=branch_time,parent_experiment_id=trim(parent_experiment_id),forcing=trim(forcing), &
institute_id=trim(institute_id),parent_experiment_rip=trim(parent_experiment_rip), &
initialization_method=initialization_method,physics_version=physics_version)
error_flag = cmor_set_cur_dataset_attribute('cosp_version',trim(COSP_VERSION))