-
Notifications
You must be signed in to change notification settings - Fork 3
/
get_obs_arg.f90
374 lines (318 loc) · 15 KB
/
get_obs_arg.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
subroutine get_obs_arg
!---------------------------------------------------------------------------
! !
! Copyright 2018 Anna Teruzzi, OGS, Trieste !
! !
! This file is part of 3DVarBio.
! 3DVarBio is based on OceanVar (Dobricic, 2006) !
! !
! 3DVarBio 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. !
! !
! 3DVarBio 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 OceanVar. If not, see <http://www.gnu.org/licenses/>. !
! !
!---------------------------------------------------------------------------
!-----------------------------------------------------------------------
! !
! Load ARGO observations !
! !
! Version 1: A. Teruzzi 2018 !
!-----------------------------------------------------------------------
use set_knd
use drv_str
use grd_str
use obs_str
use mpi_str
use filenames
implicit none
INTEGER(i4) :: k
INTEGER(i4) :: i1, kk, i
REAL(r8), ALLOCATABLE, DIMENSION(:) :: TmpFlc, TmpPar, TmpLon, TmpLat
REAL(r8), ALLOCATABLE, DIMENSION(:) :: TmpDpt, TmpTim, TmpRes, TmpErr, TmpIno
INTEGER(i4) :: GlobalArgNum, Counter, ierr
character(len=1024) :: filename
arg%no = 0
arg%nc = 0
! ---
! Allocate memory for observations
if(MyId .eq. 0) then
open(511,file=trim(ARGO_FILE))
read(511,'(I4)') GlobalArgNum
write(drv%dia,*)'Number of ARGO observations: ', GlobalArgNum
endif
call MPI_Bcast(GlobalArgNum, 1, MPI_INT, 0, Var3DCommunicator, ierr)
if(GlobalArgNum .eq. 0)then
if(MyId .eq. 0) &
close(511)
return
endif
ALLOCATE( TmpFlc(GlobalArgNum), TmpPar(GlobalArgNum))
ALLOCATE( TmpLon(GlobalArgNum), TmpLat(GlobalArgNum))
ALLOCATE( TmpDpt(GlobalArgNum), TmpTim(GlobalArgNum))
ALLOCATE( TmpRes(GlobalArgNum), TmpErr(GlobalArgNum))
ALLOCATE( TmpIno(GlobalArgNum))
if(MyId .eq. 0) then
! process 0 reads all the argo observations
do k=1,GlobalArgNum
read (511,*) &
TmpFlc(k), TmpPar(k), &
TmpLon(k), TmpLat(k), &
TmpDpt(k), TmpTim(k), &
TmpRes(k), TmpErr(k), TmpIno(k)
end do
close (511)
endif
call MPI_Bcast(TmpFlc, GlobalArgNum, MPI_REAL8, 0, Var3DCommunicator, ierr)
call MPI_Bcast(TmpPar, GlobalArgNum, MPI_REAL8, 0, Var3DCommunicator, ierr)
call MPI_Bcast(TmpLon, GlobalArgNum, MPI_REAL8, 0, Var3DCommunicator, ierr)
call MPI_Bcast(TmpLat, GlobalArgNum, MPI_REAL8, 0, Var3DCommunicator, ierr)
call MPI_Bcast(TmpDpt, GlobalArgNum, MPI_REAL8, 0, Var3DCommunicator, ierr)
call MPI_Bcast(TmpTim, GlobalArgNum, MPI_REAL8, 0, Var3DCommunicator, ierr)
call MPI_Bcast(TmpRes, GlobalArgNum, MPI_REAL8, 0, Var3DCommunicator, ierr)
call MPI_Bcast(TmpErr, GlobalArgNum, MPI_REAL8, 0, Var3DCommunicator, ierr)
call MPI_Bcast(TmpIno, GlobalArgNum, MPI_REAL8, 0, Var3DCommunicator, ierr)
! Counting the number of observations that falls in the domain
Counter = 0
do k=1,GlobalArgNum
if( TmpLon(k) .ge. grd%lon(1,1) .and. TmpLon(k) .lt. grd%NextLongitude .and. &
TmpLat(k) .ge. grd%lat(1,1) .and. TmpLat(k) .lt. grd%lat(grd%im,grd%jm) ) then
Counter = Counter + 1
endif
enddo
if(drv%Verbose .eq. 1) &
print*, "MyId", MyId, "has",Counter,"ARGO observations"
arg%no = Counter
ALLOCATE ( arg%ino(arg%no), arg%flg(arg%no), arg%flc(arg%no), arg%par(arg%no))
ALLOCATE ( arg%lon(arg%no), arg%lat(arg%no), arg%dpt(arg%no), arg%tim(arg%no))
ALLOCATE ( arg%inc(arg%no))
ALLOCATE ( arg%err(arg%no))
ALLOCATE ( arg%res(arg%no))
ALLOCATE ( arg%ib(arg%no), arg%jb(arg%no), arg%kb(arg%no))
ALLOCATE ( arg%pb(arg%no), arg%qb(arg%no), arg%rb(arg%no))
ALLOCATE ( arg%pq1(arg%no), arg%pq2(arg%no), arg%pq3(arg%no), arg%pq4(arg%no))
ALLOCATE ( arg%pq5(arg%no), arg%pq6(arg%no), arg%pq7(arg%no), arg%pq8(arg%no))
Counter = 0
do k=1,GlobalArgNum
if( TmpLon(k) .ge. grd%lon(1,1) .and. TmpLon(k) .lt. grd%NextLongitude .and. &
TmpLat(k) .ge. grd%lat(1,1) .and. TmpLat(k) .lt. grd%lat(grd%im,grd%jm) ) then
Counter = Counter + 1
arg%flc(Counter) = TmpFlc(k)
arg%par(Counter) = TmpPar(k)
arg%lon(Counter) = TmpLon(k)
arg%lat(Counter) = TmpLat(k)
arg%dpt(Counter) = TmpDpt(k)
arg%res(Counter) = TmpRes(k)
arg%err(Counter) = TmpErr(k)
arg%ino(Counter) = TmpIno(k)
endif
enddo
! DECOMMENT FOLLOWING TWO LINES TO MAKE FILTER TEST
! arg%res(:) = 1
! arg%err(:) = 1.d-2
! ---
! Initialise quality flag
arg%flg(:) = 1
arg%rb(:) = 0
! ---
! Vertical interpolation parameters
do k = 1,arg%no
if(arg%flg(k).eq.1)then
arg%kb(k) = grd%km-1
do kk = 1,grd%km-1
if( arg%dpt(k).ge.grd%dep(kk) .and. arg%dpt(k).lt.grd%dep(kk+1) ) then
arg%kb(k) = kk
arg%rb(k) = (arg%dpt(k) - grd%dep(kk)) / (grd%dep(kk+1) - grd%dep(kk))
else if ( arg%dpt(k).ge.0 .and. arg%dpt(k).lt.grd%dep(1)) then
arg%kb(k) = 1
endif
enddo
endif
enddo
! ---
! Count good observations
arg%nc = 0
do k=1,arg%no
if(arg%flg(k).eq.1)then
arg%nc = arg%nc + 1
else
arg%res(k) = 0.
arg%inc(k) = 0.
arg%pq1(k) = 0.
arg%pq2(k) = 0.
arg%pq3(k) = 0.
arg%pq4(k) = 0.
arg%pq5(k) = 0.
arg%pq6(k) = 0.
arg%pq7(k) = 0.
arg%pq8(k) = 0.
endif
enddo
arg%flc(:) = arg%flg(:)
DEALLOCATE( TmpFlc, TmpPar)
DEALLOCATE( TmpLon, TmpLat)
DEALLOCATE( TmpDpt, TmpTim)
DEALLOCATE( TmpRes, TmpErr)
DEALLOCATE( TmpIno)
end subroutine get_obs_arg
subroutine int_par_arg
!-----------------------------------------------------------------------
! !
! Get interpolation parameters for a grid !
! !
! Version 1: A. Teruzzi 2018 !
!-----------------------------------------------------------------------
use set_knd
use drv_str
use grd_str
use obs_str
use mpi_str
implicit none
INTEGER(i4) :: i, j, k, ierr
INTEGER(i4) :: i1, j1, k1, idep
REAL(r8) :: p1, q1, r1
REAL(r8) :: msk4, div_x, div_y
LOGICAL :: ins
ins(i,i1) = i.ge.1 .and. i.le.i1
if(arg%no.gt.0) then
arg%flc(:) = arg%flg(:)
! ---
! Horizontal interpolation parameters
do k = 1,arg%no
do j=1,grd%jm-1
do i=1,grd%im-1
if( grd%lat(i,j).le.arg%lat(k) .and. grd%lat(i,j+1).gt.arg%lat(k) .and. &
grd%lon(i,j).le.arg%lon(k) .and. grd%lon(i+1,j).gt.arg%lon(k) ) then
j1 = j
i1 = i
q1 = j1 + (arg%lat(k) - grd%lat(i,j)) / (grd%lat(i,j+1) - grd%lat(i,j))
p1 = i1 + (arg%lon(k) - grd%lon(i,j)) / (grd%lon(i+1,j) - grd%lon(i,j))
else if( i .eq. grd%im-1 .and. grd%lat(i,j).le.arg%lat(k) .and. grd%lat(i,j+1).gt.arg%lat(k) .and. &
grd%lon(grd%im,j).le.arg%lon(k) .and. grd%NextLongitude.gt.arg%lon(k) ) then
j1 = j
i1 = grd%im
q1 = j1 + (arg%lat(k) - grd%lat(i,j)) / (grd%lat(i,j+1) - grd%lat(i,j))
p1 = i1 + (arg%lon(k) - grd%lon(grd%im,j)) / (grd%NextLongitude - grd%lon(grd%im,j))
endif
enddo
enddo
! q1 = (arg%lat(k) - grd%lat(1,1)) / grd%dlt + 1.0
! j1 = int(q1)
! p1 = (arg%lon(k) - grd%lon(1,1)) / grd%dln + 1.0
! i1 = int(p1)
if(ins(j1,grd%jm) .and. ins(i1,grd%im)) then
arg%ib(k) = i1
arg%jb(k) = j1
arg%pb(k) = (p1-i1)
arg%qb(k) = (q1-j1)
else
arg%flc(k) = 0
endif
enddo
! ---
! Undefine masked for multigrid
do k = 1,arg%no
if(arg%flc(k).eq.1)then
i1 = arg%ib(k)
j1 = arg%jb(k)
idep = arg%kb(k)+1
msk4 = grd%global_msk(GlobalRowOffset+i1,j1,idep) + grd%global_msk(GlobalRowOffset+i1+1,j1,idep) + &
grd%global_msk(GlobalRowOffset+i1,j1+1,idep) + grd%global_msk(GlobalRowOffset+i1+1,j1+1,idep)
if(msk4.lt.1.) arg%flc(k) = 0
endif
enddo
! ---
! Horizontal interpolation parameters for each masked grid
do k = 1,arg%no
if(arg%flc(k) .eq. 1) then
i1=arg%ib(k)
p1=arg%pb(k)
j1=arg%jb(k)
q1=arg%qb(k)
k1=arg%kb(k)
div_y = (1.-q1) * max(grd%global_msk(GlobalRowOffset+i1,j1 ,k1),grd%global_msk(GlobalRowOffset+i1+1,j1 ,k1)) &
+ q1 * max(grd%global_msk(GlobalRowOffset+i1,j1+1,k1),grd%global_msk(GlobalRowOffset+i1+1,j1+1,k1))
div_x = (1.-p1) * grd%global_msk(GlobalRowOffset+i1 ,j1,k1) + p1 * grd%global_msk(GlobalRowOffset+i1+1,j1,k1)
arg%pq1(k) = grd%global_msk(GlobalRowOffset+i1,j1,k1) &
* max(grd%global_msk(GlobalRowOffset+i1,j1,k1),grd%global_msk(GlobalRowOffset+i1+1,j1,k1)) &
* (1.-p1) * (1.-q1) &
/( div_x * div_y + 1.e-16 )
arg%pq2(k) = grd%global_msk(GlobalRowOffset+i1+1,j1,k1) &
* max(grd%global_msk(GlobalRowOffset+i1,j1,k1),grd%global_msk(GlobalRowOffset+i1+1,j1,k1)) &
* p1 * (1.-q1) &
/( div_x * div_y + 1.e-16 )
div_x = (1.-p1) * grd%global_msk(GlobalRowOffset+i1 ,j1+1,k1) + p1 * grd%global_msk(GlobalRowOffset+i1+1,j1+1,k1)
arg%pq3(k) = grd%global_msk(GlobalRowOffset+i1,j1+1,k1) &
* max(grd%global_msk(GlobalRowOffset+i1,j1+1,k1),grd%global_msk(GlobalRowOffset+i1+1,j1+1,k1)) &
* (1.-p1) * q1 &
/( div_x * div_y + 1.e-16 )
arg%pq4(k) = grd%global_msk(GlobalRowOffset+i1+1,j1+1,k1) &
* max(grd%global_msk(GlobalRowOffset+i1,j1+1,k1),grd%global_msk(GlobalRowOffset+i1+1,j1+1,k1)) &
* p1 * q1 &
/( div_x * div_y + 1.e-16 )
k1=arg%kb(k) + 1
div_y = (1.-q1) * max(grd%global_msk(GlobalRowOffset+i1,j1 ,k1),grd%global_msk(GlobalRowOffset+i1+1,j1 ,k1)) &
+ q1 * max(grd%global_msk(GlobalRowOffset+i1,j1+1,k1),grd%global_msk(GlobalRowOffset+i1+1,j1+1,k1))
div_x = (1.-p1) * grd%global_msk(GlobalRowOffset+i1 ,j1,k1) + p1 * grd%global_msk(GlobalRowOffset+i1+1,j1,k1)
arg%pq5(k) = grd%global_msk(GlobalRowOffset+i1,j1,k1) &
* max(grd%global_msk(GlobalRowOffset+i1,j1,k1),grd%global_msk(GlobalRowOffset+i1+1,j1,k1)) &
* (1.-p1) * (1.-q1) &
/( div_x * div_y + 1.e-16 )
arg%pq6(k) = grd%global_msk(GlobalRowOffset+i1+1,j1,k1) &
* max(grd%global_msk(GlobalRowOffset+i1,j1,k1),grd%global_msk(GlobalRowOffset+i1+1,j1,k1)) &
* p1 * (1.-q1) &
/( div_x * div_y + 1.e-16 )
div_x = (1.-p1) * grd%global_msk(GlobalRowOffset+i1 ,j1+1,k1) + p1 * grd%global_msk(GlobalRowOffset+i1+1,j1+1,k1)
arg%pq7(k) = grd%global_msk(GlobalRowOffset+i1,j1+1,k1) &
* max(grd%global_msk(GlobalRowOffset+i1,j1+1,k1),grd%global_msk(GlobalRowOffset+i1+1,j1+1,k1)) &
* (1.-p1) * q1 &
/( div_x * div_y + 1.e-16 )
arg%pq8(k) = grd%global_msk(GlobalRowOffset+i1+1,j1+1,k1) &
* max(grd%global_msk(GlobalRowOffset+i1,j1+1,k1),grd%global_msk(GlobalRowOffset+i1+1,j1+1,k1)) &
* p1 * q1 &
/( div_x * div_y + 1.e-16 )
r1=arg%rb(k)
arg%pq1(k) = (1.-r1) * arg%pq1(k)
arg%pq2(k) = (1.-r1) * arg%pq2(k)
arg%pq3(k) = (1.-r1) * arg%pq3(k)
arg%pq4(k) = (1.-r1) * arg%pq4(k)
arg%pq5(k) = r1 * arg%pq5(k)
arg%pq6(k) = r1 * arg%pq6(k)
arg%pq7(k) = r1 * arg%pq7(k)
arg%pq8(k) = r1 * arg%pq8(k)
if(arg%pq1(k) .lt. 1.E-16) arg%pq1(k) = dble(0)
if(arg%pq2(k) .lt. 1.E-16) arg%pq2(k) = dble(0)
if(arg%pq3(k) .lt. 1.E-16) arg%pq3(k) = dble(0)
if(arg%pq4(k) .lt. 1.E-16) arg%pq4(k) = dble(0)
if(arg%pq5(k) .lt. 1.E-16) arg%pq5(k) = dble(0)
if(arg%pq6(k) .lt. 1.E-16) arg%pq6(k) = dble(0)
if(arg%pq7(k) .lt. 1.E-16) arg%pq7(k) = dble(0)
if(arg%pq8(k) .lt. 1.E-16) arg%pq8(k) = dble(0)
endif
enddo
! ---
! Count good observations
arg%nc = 0
do k=1,arg%no
if(arg%flc(k).eq.1)then
arg%nc = arg%nc + 1
endif
enddo
endif
arg%nc_global = 0
call MPI_Allreduce(arg%nc, arg%nc_global, 1, MPI_INT, MPI_SUM, Var3DCommunicator, ierr)
if(MyId .eq. 0) then
write(drv%dia,*)'Real number of ARGO observations: ',arg%nc_global
print*,'Good argo observations: ',arg%nc_global
end if
DEALLOCATE ( arg%ino, arg%flg, arg%par)
DEALLOCATE ( arg%lon, arg%lat, arg%dpt, arg%tim)
DEALLOCATE ( arg%pb, arg%qb, arg%rb)
end subroutine int_par_arg