-
Notifications
You must be signed in to change notification settings - Fork 26
/
Copy pathmo_Gdata_read.F90
764 lines (607 loc) · 25.6 KB
/
mo_Gdata_read.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
! Copyright (C) 2020 J. Schwinger
!
! This file is part of BLOM/iHAMOCC.
!
! BLOM is free software: you can redistribute it and/or modify it under the
! terms of the GNU Lesser General Public License as published by the Free
! Software Foundation, either version 3 of the License, or (at your option)
! any later version.
!
! BLOM 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 Lesser General Public License for
! more details.
!
! You should have received a copy of the GNU Lesser General Public License
! along with BLOM. If not, see https://www.gnu.org/licenses/.
module mo_Gdata_read
!*************************************************************************************************
! Routines for reading initial condition files for OMIP-BGC, based on WOA 2013 and GLODAPv2
!
! J.Schwinger, *Gfi, Bergen* 2011-05-19
!
! Public routines and variable of this module:
! -subroutine set_Gdata
! Initialise global varibles and read in one data set. Must be
! called before the processing of one data set starts.
! -subroutine clean_Gdata
! Deallocate global fields of this module and reset all global variables.
! Should be called each time, the processing of one data set is finished.
! -subroutine get_profile
! Returns one profile from the currently open data set (opened by a
! previous call to set_Gdata). See header of get profile for details.
! -function get_region
! Returns the index of the region a given point belongs to. If no region
! is found get_region returns 0, which is the index of the 'global region'.
! Note that the regions are defined below in the module header.
!
! Modified
! J.Schwinger, *Uni Climate, BCCR* 2017-07-07
! - adapted this module to read the initial conditions for OMIP-BGC.
! J.Schwinger, *Uni Research, Bergen* 2018-04-12
! - adaptions for reading c-isotope initial values as d13C and d14C
!*************************************************************************************************
use netcdf, only: nf90_noerr,nf90_nowrite,nf90_strerror,nf90_inq_dimid, &
nf90_inquire_dimension,nf90_inq_varid,nf90_get_var, &
nf90_inquire_variable,nf90_get_att,nf90_close,nf90_open
use mod_xc, only: mnproc,xchalt
use mo_control_bgc, only: io_stdo_bgc
implicit none
private
! routines
public :: set_Gdata
public :: get_profile
public :: get_region
public :: clean_Gdata
private :: set_regional_profiles
private :: read_Gdata
private :: calc_mean_profile
! module variables
public :: zlev ! Depth of each z-level [m] in the current data file.
public :: zlev_bnds
public :: nzmax,nz,fillval
public :: inidic,inialk,inipo4,inioxy,inino3,inisil,inid13c,inid14c
! Number of latitudes, longitudes, and z-levels in the WOA and GLODAP data
integer, parameter :: nlon = 360
integer, parameter :: nlat = 180
integer, parameter :: nz_woa = 102 ! Number of z-levels in the WOA data files.
integer, parameter :: nz_glo = 33 ! Number of z-levels in the GLODAP data files.
integer, parameter :: nzmax = nz_woa ! Max nuber of z-levels (=nzwoa)
! Resolution of data in degree
real, parameter :: dres = 1.0
! Max number of gridpoints to select around the center for averaging in
! longitude direction
integer, parameter :: dnmax = 100.0
! Fill value used in this module, original fill values of data files are
! replaced by this fill value during read
real, parameter :: fillval = -1.e+32
! Input file names (incl. full path) set through namelist
character(len=256) :: inidic = ''
character(len=256) :: inialk = ''
character(len=256) :: inipo4 = ''
character(len=256) :: inioxy = ''
character(len=256) :: inino3 = ''
character(len=256) :: inisil = ''
character(len=256) :: inid13c = ''
character(len=256) :: inid14c = ''
character(len=256) :: inic13 = '' ! currently not used
character(len=256) :: inic14 = '' ! currently not used
! Variables set by call to Gdata_set
integer :: nz
real :: cfac, ddeg
real, dimension(:), allocatable :: lon,lat,zlev
real, dimension(:,:), allocatable :: zlev_bnds
real, dimension(:, :, :), allocatable :: rvar,gdata
character(len=16) :: var,ncname
character(len=3) :: dsrc
character(len=256) :: infile
logical :: lset = .false.
!-----------------------------------------
! Definitions for regional mean profiles:
!-----------------------------------------
type region
character(len=64) :: name ! Region name
integer :: idx ! Region index
integer :: npts(nzmax) ! nb of valid data points at each level
real :: clon, clat ! center longitude and latitude
real :: dlon, dlat ! latitude and longitude extent
real :: mprf(nzmax) ! mean profile for region
logical :: global ! global extent T/F
end type region
integer, parameter :: nreg=10
type(region) :: rg(0:nreg)
! Set regions for fall-back profiles
! Global profile;
data rg(0)%idx, rg(0)%name / 0, 'global' /
data rg(0)%clon, rg(0)%clat / 0.0, 0.0 /
data rg(0)%dlon, rg(0)%dlat / 360.0, 180.0 /
data rg(0)%global / .true. /
! Indian Ocean
data rg(1)%idx, rg(1)%name / 1, 'Indian Ocean' /
data rg(1)%clon, rg(1)%clat / 65.0,-10.0 /
data rg(1)%dlon, rg(1)%dlat / 90.0, 80.0 /
data rg(1)%global / .false. /
! North Atlantic
data rg(2)%idx, rg(2)%name / 2, 'North Atlantic' /
data rg(2)%clon, rg(2)%clat / 0.0, 70.0 /
data rg(2)%dlon, rg(2)%dlat / 180.0, 40.0 /
data rg(2)%global / .false. /
! northern subtropical Atlantic
data rg(3)%idx, rg(3)%name / 3, 'Northern subtropical Atlantic' /
data rg(3)%clon, rg(3)%clat / 330.0, 35.0 /
data rg(3)%dlon, rg(3)%dlat / 140.0, 30.0 /
data rg(3)%global / .false. /
! Tropical Atlantic
data rg(4)%idx, rg(4)%name / 4, 'Tropical Atlantic' /
data rg(4)%clon, rg(4)%clat / 335.0, 0.0 /
data rg(4)%dlon, rg(4)%dlat / 90.0, 40.0 /
data rg(4)%global / .false. /
! Southern subtropical Atlantic
data rg(5)%idx, rg(5)%name / 5, 'Southern subtropical Atlantic' /
data rg(5)%clon, rg(5)%clat / 335.0, -35.0 /
data rg(5)%dlon, rg(5)%dlat / 90.0, 30.0 /
data rg(5)%global / .false. /
! North Pacific
data rg(6)%idx, rg(6)%name / 6, 'North Pacific' /
data rg(6)%clon, rg(6)%clat / 180.0, 70.0 /
data rg(6)%dlon, rg(6)%dlat / 180.0, 40.0 /
data rg(6)%global / .false. /
! northern subtropical Pacific
data rg(7)%idx, rg(7)%name / 7, 'Northern subtropical Pacific' /
data rg(7)%clon, rg(7)%clat / 185.0, 35.0 /
data rg(7)%dlon, rg(7)%dlat / 150.0, 30.0 /
data rg(7)%global / .false. /
! Tropical Pacific
data rg(8)%idx, rg(8)%name / 8, 'Tropical Pacific' /
data rg(8)%clon, rg(8)%clat / 200.0, 0.0 /
data rg(8)%dlon, rg(8)%dlat / 180.0, 40.0 /
data rg(8)%global / .false. /
! Southern subtropical Pacific
data rg(9)%idx, rg(9)%name / 9, 'Southern subtropical Pacific' /
data rg(9)%clon, rg(9)%clat / 200.0, -35.0 /
data rg(9)%dlon, rg(9)%dlat / 180.0, 30.0 /
data rg(9)%global / .false. /
! Southern Ocean
data rg(10)%idx, rg(10)%name / 10, 'Southern Ocean' /
data rg(10)%clon, rg(10)%clat / 180.0, -70.0 /
data rg(10)%dlon, rg(10)%dlat / 360.0, 40.0 /
data rg(10)%global / .false. /
contains
subroutine set_Gdata(vname,inddeg)
!--------------------------------------------------------------------------------
! Initialise global varibles and read data set specified by vname. Must be
! called before the first call to any routine of this module.
!--------------------------------------------------------------------------------
! Arguments
character(len=*), intent(in) :: vname ! data set name to read in
real, intent(in) :: inddeg ! extent (in degrees) of region used for averaging
! Valid values of vname are:
! 'pho' - WOA phosphate
! 'nit' - WOA nitrate
! 'sil' - WOA silicate
! 'oxy' - WOA dissolved oxygen
! 'alk' - GLODAP alkalinity
! 'dic' - GLODAP dissolved inorganic carbon
! 'C13' - Dissolved inorganic 13C carbon isotope
! 'd13' - delta13C of dissolved inorganic carbon
! 'C14' - Dissolved inorganic 14C carbon isotope
! 'd14' - delta14C of dissolved inorganic carbon
! Local variables
character(len=*), parameter :: routinestr = 'set_Gdata'
if( allocated(lon) ) deallocate( lon )
if( allocated(lat) ) deallocate( lat )
if( allocated(zlev) ) deallocate( zlev )
if( allocated(zlev_bnds) ) deallocate( zlev_bnds )
if( allocated(rvar) ) deallocate( rvar )
if( allocated(gdata) ) deallocate( gdata )
! Select settings specific to each variable
select case (vname)
case ('pho') ! phosphate
infile = inipo4
ncname = 'po4'
dsrc = 'WOA'
cfac = 1.0e-6 ! data in mumol/L -> kmol/m3
case ('nit') ! nitrate
infile = inino3
ncname = 'no3'
dsrc = 'WOA'
cfac = 1.0e-6 ! data in mumol/L -> kmol/m3
case ('sil') ! silicate
infile = inisil
ncname = 'si'
dsrc = 'WOA'
cfac = 1.0e-6 ! data in mumol/L -> kmol/m3
case ('oxy') ! oxygen
infile = inioxy
ncname = 'o2'
dsrc = 'WOA'
cfac = 44.661*1.0e-6 ! conversion ml/L -> mumol/L -> kmol/m3
case ('alk') ! alkalinity
infile = inialk
ncname = 'At'
dsrc = 'GLO'
cfac = 1.0e-6 ! data in mumol/kg -> mol/kg
case ('dic') ! DIC
infile = inidic
ncname = 'Ct_preind'
dsrc = 'GLO'
cfac = 1.0e-6 ! data in mumol/kg -> mol/kg
case ('C13') ! natural 13C [micromoles/kg]
infile = inic13
ncname = 'C13'
dsrc = 'ISO'
cfac = 1.0e-6 ! data in mumol/kg -> mol/kg
case ('d13') ! natural delta13C [permil]
infile = inid13c
ncname = 'd13C'
dsrc = 'ISO'
cfac = 1.0
case ('C14') ! natural 14C [micromoles/kg]
infile = inic14
ncname = 'C14'
dsrc = 'ISO'
cfac = 1.0e-6 ! data in mumol/kg -> mol/kg
case ('d14') ! natural delta14C [permil]
infile = inid14c
ncname = 'd14C'
dsrc = 'ISO'
cfac = 1.0
case default
call moderr(routinestr,'Invalid vname')
end select
var = vname
ddeg = inddeg
if(mnproc == 1) then
write(io_stdo_bgc,*)
write(io_stdo_bgc,*) '********************************************'
write(io_stdo_bgc,*) 'iHAMOCC: initialising ', trim(vname)
endif
call read_Gdata()
! extend data array by +/-dnmax data points in longitude
allocate( gdata(-dnmax:nlon+dnmax,nlat,nz) )
gdata(:,:,:) = 0.0
gdata( 1:nlon, :,:) = rvar(:,:,:)
gdata(-dnmax:0, :,:) = rvar(nlon-dnmax:nlon,:,:)
gdata(nlon+1:nlon+dnmax,:,:) = rvar(1:dnmax,:,:)
lset = .true.
call set_regional_profiles()
end subroutine set_Gdata
subroutine get_profile(clon,clat,prf)
!--------------------------------------------------------------------------------
! Return a profile suitable for initialisation of HAMCC at point clon/clat.
! A mean profile is calculated by calling calc_mean_profile with the settings
! defined by a previous call to set_Gdata. If no valid data is found for the
! point clon/clat, it is tried to obtain a mean regional profile (e.g. for the
! north atlantic area). These mean profiles are initialised as part of
! set_Gdata.
!--------------------------------------------------------------------------------
! Arguments
real, intent(in) :: clon ! center lon of mean profile
real, intent(in) :: clat ! center lat of mean profile
real, intent(out) :: prf(nzmax) ! mean profile for initialisation
! Local variables
integer :: idx, npts(nzmax)
real :: clon_tmp,clat_tmp
character(len=*), parameter :: routinestr = 'mo_Gdata_read, get_profile'
if( .not. lset ) then
call moderr(routinestr, ' Module not initialised yet')
end if
if( clon < 0 ) then
clon_tmp=clon+360.0
clat_tmp=clat
else
clon_tmp=clon
clat_tmp=clat
endif
! Try to obtain a mean profile for a region centered at clon/clat
call calc_mean_profile(clon_tmp,clat_tmp,ddeg,ddeg,prf,npts)
! Fall back to regional profile if number of valid data points is smaller
! than 3 for the surface layer. A global mean profile is used if
! get_region returns 0.
if( npts(1) < 3 ) then
idx = get_region(clon_tmp,clat_tmp)
prf = rg(idx)%mprf
!write(*,*) 'Region is ', rg(idx)%name, clon, clat
endif
end subroutine get_profile
function get_region(clon,clat)
!--------------------------------------------------------------------------------
! Return index of region the point clon/clat belongs to
! The rectangular regions as defined in the module header (and stored in the
! data type 'rg') are searched. If point clon/clat belongs to region i, the
! index i is the result of this function. If no region is found, get_region
! returns 0, which is the index of the 'global' region defined in the header.
!--------------------------------------------------------------------------------
! ARGUMENTS
real, intent(in) :: clon,clat ! lon/lat of point
! Local variables
integer :: get_region
integer :: i
real :: ll_lon, ur_lon
real :: ll_lat, ur_lat
logical :: boundwithin, found
character(len=*), parameter :: routinestr = 'mo_Gdata_read, get_region'
if( clon < 0 ) call moderr(routinestr, ' clon must be in the range [0,360]')
if( clon > 360.0 ) call moderr(routinestr, ' clon must be in the range [0,360]')
found = .false.
do i=1,nreg
boundwithin = .false.
ll_lon = rg(i)%clon-rg(i)%dlon/2.0
ur_lon = rg(i)%clon+rg(i)%dlon/2.0
ll_lat = rg(i)%clat-rg(i)%dlat/2.0
ur_lat = rg(i)%clat+rg(i)%dlat/2.0
if( ll_lon < 0.0 ) ll_lon = ll_lon+360.0
if( ur_lon > 360.0 ) ur_lon = ur_lon-360.0
if( ll_lon > ur_lon ) boundwithin = .true.
if( clat < ll_lat .or. clat > ur_lat ) cycle
if( boundwithin ) then
if( clon < ll_lon .and. clon > ur_lon ) cycle
else
if( clon < ll_lon .or. clon > ur_lon ) cycle
endif
found = .true.
exit
enddo
if( found ) then
get_region = rg(i)%idx
else
get_region = 0
endif
end function get_region
subroutine set_regional_profiles()
!--------------------------------------------------------------------------------
! Calculate the mean profiles in regions as defined in the module header
!--------------------------------------------------------------------------------
! Local variables
integer :: i
character(len=*), parameter :: routinestr = 'mo_Gdata_read, set_regional_profiles'
if( .not. lset ) call moderr(routinestr, ' Module not initialised yet')
do i=0,nreg
call calc_mean_profile(rg(i)%clon,rg(i)%clat,rg(i)%dlon,rg(i)%dlat, &
rg(i)%mprf,rg(i)%npts,rg(i)%global)
!write(*,*) 'Calculated mean profile for ', rg(i)%name
!write(*,*) '==============='
!write(*,*) rg(i)%mprf
!write(*,*) '==============='
enddo
end subroutine set_regional_profiles
subroutine read_Gdata()
!--------------------------------------------------------------------------------
! Read the WOA or GLODAP data into variables lon/lat/zlev and rvar
!--------------------------------------------------------------------------------
! Local variables
integer :: ncId, vId, dId
integer :: numlon, numlat, numlev
integer :: i, ndim, natts
integer :: dimid(7)
integer :: status
real :: fval
character(len=16) :: lonstr,latstr,depthstr,depthbndsstr,fvalstr
character(len=*), parameter :: routinestr = 'mo_Gdata_read, read_Gdata'
lonstr = 'lon'
latstr = 'lat'
fvalstr = '_FillValue'
select case (dsrc)
case ('WOA')
nz = nz_woa
depthstr='depth'
depthbndsstr='depth_bnds'
case ('GLO')
nz = nz_glo
depthstr='depthz'
depthbndsstr='depthz_bnds'
case ('ISO')
nz = nz_glo
depthstr='depthz'
depthbndsstr='depthz_bnds'
case default
call moderr(routinestr,'Invalid dsrc')
end select
! Open file
if(mnproc == 1) write(io_stdo_bgc,*) 'Reading ', trim(infile)
status = nf90_open(infile,nf90_nowrite,ncid); call ncerr(status)
! Get dimensions
status = nf90_inq_dimid(ncid, trim(lonstr), dId)
call ncerr(status)
status = nf90_inquire_dimension(ncid, dID, len=numlon)
call ncerr(status)
status = nf90_inq_dimid(ncid, trim(latstr), dId)
call ncerr(status)
status = nf90_inquire_dimension(ncid, dID, len=numlat)
call ncerr(status)
status = nf90_inq_dimid(ncid, trim(depthstr), dId)
call ncerr(status)
status = nf90_inquire_dimension(ncid, dId, len=numlev)
call ncerr(status)
if( numlon /= nlon .or. numlat /= nlat .or. numlev /= nz ) &
call moderr(routinestr,'Unexpected nb of elements in data file')
allocate( lon(nlon), lat(nlat), zlev(nz), zlev_bnds(2,nz) )
allocate( rvar(nlon,nlat,nz) )
! Get lon, lat, and lev
status = nf90_inq_varid(ncid, trim(lonstr), vId)
call ncerr(status)
status = nf90_get_var(ncid, vId, lon)
call ncerr(status)
status = nf90_inq_varid(ncid, trim(latstr), vId)
call ncerr(status)
status = nf90_get_var(ncid, vId, lat)
call ncerr(status)
status = nf90_inq_varid(ncid, trim(depthstr), vId)
call ncerr(status)
status = nf90_get_var(ncid, vId, zlev)
call ncerr(status)
status = nf90_inq_varid(ncid, trim(depthbndsstr), vId)
call ncerr(status)
status = nf90_get_var(ncid, vId, zlev_bnds)
call ncerr(status)
! Get varid and fill value
status = nf90_inq_varid(ncid, ncname, vId)
call ncerr(status)
status = nf90_inquire_variable(ncid, vId, ndims=ndim, dimids=dimid, nAtts=natts)
call ncerr(status)
status = nf90_get_att(ncid, vid, trim(fvalstr), fval)
call ncerr(status)
! GetRead the data
status = nf90_get_var(ncid, vId, rvar)
call ncerr(status)
! arrange data to correspond to [0,360] in longitude
select case (dsrc)
case ('WOA')
lon = cshift(lon, -180)
rvar = cshift(rvar,-180,1)
case ('GLO')
lon = cshift(lon, -20)
rvar = cshift(rvar,-20,1)
case ('ISO')
lon = cshift(lon, -180)
rvar = cshift(rvar,-180,1)
end select
do i=1,nlon
if(lon(i)< 0.0) lon(i)=lon(i)+360.0
if(lon(i)>360.0) lon(i)=lon(i)-360.0
enddo
! Fillvalues are assumed to be < 0 currently, otherwise the below code would fail
if(fval > 0.0) call moderr(routinestr,'FillValue > 0 found in data')
where( rvar < fval*0.1 )
! Replace fill values:
rvar = fillval
elsewhere
! unit conversion
rvar = rvar*cfac
end where
! Close data file
status = nf90_close(ncid)
call ncerr(status)
end subroutine read_Gdata
subroutine calc_mean_profile(clon,clat,dlon,dlat,prf,npts,global)
!--------------------------------------------------------------------------------
! Return mean profile around the center point clon/clat.
! The mean profile is calculated from valid data points in the square defined
! by clon+/-dlon/2 clat+/-dlat/2. The number of valid data points per depth
! level is returned in npts. By setting the optional argument global to true,
! all valid data points are used to calculate a global mean profile. clon, clat,
! dlon, and dlat are ignored in this case.
!--------------------------------------------------------------------------------
! Arguments
real, intent(in) :: clon, clat ! center lon/lat of mean profile
real, intent(in) :: dlon, dlat ! lon/lat extent of region to select for averaging
real, intent(out) :: prf(nzmax) ! mean profile calculated from all data in selected region
integer, intent(out) :: npts(nzmax) ! nb of valid data points found for each depth level
logical, optional, intent(in) :: global ! if set to true, calculate mean over the whole data set
! Local variables
integer :: ilonc, ilons, ilone, dnlon
integer :: ilatc, ilats, ilate, dnlat
integer :: l, nelmlon,nelmlat
logical :: gl
character(len=*), parameter :: routinestr = 'mo_Gdata_read, calc_mean_profile'
if( .not. lset ) call moderr(routinestr, ' Module not initialised yet')
if( clon < 0 ) call moderr(routinestr, ' clon must be in the range [0,360]')
if( clon > 360.0 ) call moderr(routinestr, ' clon must be in the range [0,360]')
prf(:) = fillval
npts(:) = 0.0
gl = .false.
if( present(global) ) gl=global
if( gl ) then
ilons=1
ilone=nlon
ilats=1
ilate=nlat
else
! Find index of nearest gridpoint (not exact but okay for this purpose)
do ilonc=1,nlon
if( clon < lon(ilonc) ) exit
enddo
if( ilonc > nlon ) ilonc = nlon
if( lon(ilonc)-clon > dres/2.0 ) ilonc=ilonc-1
if( ilonc < 1 ) ilonc = 1
do ilatc=1,nlat
if( clat < lat(ilatc) ) exit
enddo
if( ilatc > nlat ) ilatc = nlat
if( lat(ilatc)-clat > dres/2.0 ) ilatc=ilatc-1
if( ilatc < 1 ) ilatc = 1
dnlon = int(dlon/2.0*dres) ! Nb of gridpoints to select around the center lon
dnlat = int(dlat/2.0*dres) ! Nb of gridpoints to select around the center lat
nelmlon = 2*dnlon+1
nelmlat = 2*dnlat+1
! Start idices of rectangle:
ilons = ilonc-dnlon
ilats = ilatc-dnlat
! There is no "wrap-around" if southern/northen boundary of rectangle
! goes beyond the pole. Instead rectangle is adjusted such that boundaries
! are aligned with northernmost/sothernmost data gridpoint
if(ilats <= 0 ) ilats=1
if(ilats > nlat-nelmlat+1) ilats= nlat-nelmlat+1
! End indices of rectangle:
ilone = ilons+nelmlon-1
ilate = ilats+nelmlat-1
if( ilons < -dnmax ) call moderr(routinestr,'error: data array too small')
if( ilone > dnmax+nlon ) call moderr(routinestr,'error: data array too small')
endif
! Calculate mean profile:
do l=1,nz
npts(l) = count(gdata(ilons:ilone,ilats:ilate,l) > fillval*0.1)
prf(l) = sum(gdata(ilons:ilone,ilats:ilate,l), mask=gdata(ilons:ilone,ilats:ilate,l) > fillval*0.1)
if( npts(l) > 0) then
prf(l) = prf(l)/npts(l)
else
prf(l) = fillval
endif
enddo
!write(*,*) '================'
!if( gl ) then
! write(*,*) 'global'
!else
! write(*,*) dnlon,dnlat
! write(*,*) ilonc,ilons,ilone,lon(ilonc)
! write(*,*) ilatc,ilats,ilate,lat(ilatc)
!endif
!write(*,*) '================'
!--------------------------------------------------------------------------------
end subroutine calc_mean_profile
subroutine clean_Gdata()
!--------------------------------------------------------------------------------
! Deallocate fields and reset global variables
!--------------------------------------------------------------------------------
if( allocated(lon) ) deallocate( lon )
if( allocated(lat) ) deallocate( lat )
if( allocated(zlev) ) deallocate( zlev )
if( allocated(zlev_bnds) ) deallocate( zlev_bnds )
if( allocated(rvar) ) deallocate( rvar )
if( allocated(gdata) ) deallocate( gdata )
infile = ''
ncname = ''
var = ''
dsrc = ''
cfac = 1.0
ddeg = 0.0
nz = 0
lset = .false.
!--------------------------------------------------------------------------------
end subroutine clean_Gdata
subroutine ncerr(status)
!--------------------------------------------------------------------------------
! Handle netCDF-errors
!--------------------------------------------------------------------------------
integer, intent(in) :: status
if(status == nf90_NoErr) return
write(io_stdo_bgc,*) 'NetCDF error: ',nf90_strerror(status)
write(io_stdo_bgc,*) 'Abort... '
call flush(io_stdo_bgc)
call xchalt('(Module mo_Gdata_read, ncerr)')
stop '(Module mo_Gdata_read, ncerr)'
!--------------------------------------------------------------------------------
end subroutine ncerr
subroutine moderr(routinestr,errstr)
!--------------------------------------------------------------------------------
! Handle errors, which occur in this module
!--------------------------------------------------------------------------------
character(len=*), intent(in) :: routinestr,errstr
write(io_stdo_bgc,'(/3a)') routinestr, ': ', errstr
write(io_stdo_bgc,*) 'Abort... '
call flush(io_stdo_bgc)
call xchalt('(Module mo_Gdata_read)')
stop '(Module mo_Gdata_read)'
!--------------------------------------------------------------------------------
end subroutine moderr
end module mo_Gdata_read