87
88
90
91
92
93 integer, intent(in) :: ng, model, ncid, ncvarid, tindex, gtype
94 integer, intent(in) :: LBij, UBij, LBk, UBk, Nrec
95
96 integer(i8b), intent(out), optional :: checksum
97
98 real(dp), intent(in) :: Ascl
99 real(r8), intent(out) :: Amin
100 real(r8), intent(out) :: Amax
101
102 character (len=*), intent(in) :: ncname
103 character (len=*), intent(in) :: ncvname
104
105#ifdef ASSUMED_SHAPE
106 real(r8), intent(out) :: Abry(LBij:,:,:,:)
107#else
108 real(r8), intent(out) :: Abry(LBij:UBij,LBk:UBk,4,Nrec)
109#endif
110
111
112
113 logical :: Lchecksum
114 logical, dimension(3) :: foundit
115 logical, dimension(4) :: bounded
116
117 integer :: ghost, i, ib, ij, ir, j, k, tile
118 integer :: IorJ, IJKlen, Imin, Imax, Jmin, Jmax, Klen, Npts
119 integer :: Cgrid, Istr, Iend, Jstr, Jend
120 integer, dimension(5) :: start, total
121
122 integer :: status
123
124 real(r8) :: Afactor, Aoffset, Aspval
125
126 real(r8), allocatable :: Cwrk(:)
127
128 real(r8), dimension(3) :: AttValue
129
130#if !defined PARALLEL_IO && defined DISTRIBUTE
131 real(r8), dimension(3) :: rbuffer
132#endif
133 real(r8), dimension(LBij:UBij,LBk:UBk,4,Nrec) :: wrk
134
135 character (len=12), dimension(3) :: AttName
136
137 character (len=*), parameter :: MyFile = &
138 & __FILE__//", nf90_fread3d_bry"
139
140
141
142
143
144 status=nf90_noerr
145
146
147
148
149
150
151
152
153
154
155 IF (model.eq.iadm) THEN
156 ghost=0
157 ELSE
158 ghost=1
159 END IF
160
161 SELECT CASE (gtype)
162 CASE (p2dvar, p3dvar)
163 cgrid=1
164 CASE (r2dvar, r3dvar)
165 cgrid=2
166 CASE (u2dvar, u3dvar)
167 cgrid=3
168 CASE (v2dvar, v3dvar)
169 cgrid=4
170 CASE DEFAULT
171 cgrid=2
172 END SELECT
173
174#ifdef DISTRIBUTE
175 tile=myrank
176 imin=bounds(ng)%Imin(cgrid,ghost,tile)
177 imax=bounds(ng)%Imax(cgrid,ghost,tile)
178 jmin=bounds(ng)%Jmin(cgrid,ghost,tile)
179 jmax=bounds(ng)%Jmax(cgrid,ghost,tile)
180#else
181 tile=-1
182 imin=lbij
183 imax=ubij
184 jmin=lbij
185 jmax=ubij
186#endif
187 iorj=iobounds(ng)%IorJ
188 klen=ubk-lbk+1
189 ijklen=iorj*klen
190 npts=ijklen*4*nrec
191
192
193
194 istr=bounds(ng)%Istr(tile)
195 iend=bounds(ng)%Iend(tile)
196 jstr=bounds(ng)%Jstr(tile)
197 jend=bounds(ng)%Jend(tile)
198
199
200
201 bounded(iwest )=domain(ng)%Western_Edge (tile)
202 bounded(ieast )=domain(ng)%Eastern_Edge (tile)
203 bounded(isouth)=domain(ng)%Southern_Edge(tile)
204 bounded(inorth)=domain(ng)%Northern_Edge(tile)
205
206
207
208 start(1)=1
209 total(1)=iorj
210 start(2)=1
211 total(2)=klen
212 start(3)=1
213 total(3)=4
214 start(4)=1
215 total(4)=nrec
216 start(5)=tindex
217 total(5)=1
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232 attname(1)='scale_factor'
233 attname(2)='add_offset '
234 attname(3)='_FillValue '
235
237 & attvalue, foundit, &
238 & ncid = ncid)
239 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
240 status=ioerror
241 RETURN
242 END IF
243
244 IF (.not.foundit(1)) THEN
245 afactor=1.0_r8
246 ELSE
247 afactor=attvalue(1)
248 END IF
249
250 IF (.not.foundit(2)) THEN
251 aoffset=0.0_r8
252 ELSE
253 aoffset=attvalue(2)
254 END IF
255
256 IF (.not.foundit(3)) THEN
257 aspval=spval_check
258 ELSE
259 aspval=attvalue(3)
260 END IF
261
262
263
264 IF (PRESENT(checksum)) THEN
265 lchecksum=.true.
266 checksum=0_i8b
267 ELSE
268 lchecksum=.false.
269 END IF
270
271
272
273
274
275 wrk=0.0_r8
276
277 IF (inpthread) THEN
278 status=nf90_get_var(ncid, ncvarid, wrk(lbij:,lbk:,:,:), &
279 & start, total)
280 IF (status.eq.nf90_noerr) THEN
281 amin=spval
282 amax=-spval
283 DO ir=1,nrec
284 DO ib=1,4
285 DO k=lbk,ubk
286 DO ij=lbij,ubij
287 IF (abs(wrk(ij,k,ib,ir)).ge.abs(aspval)) THEN
288 wrk(ij,k,ib,ir)=0.0_r8
289 ELSE
290 wrk(ij,k,ib,ir)=ascl* &
291 & (afactor*wrk(ij,k,ib,ir)+aoffset)
292 amin=min(amin,wrk(ij,k,ib,ir))
293 amax=max(amax,wrk(ij,k,ib,ir))
294 END IF
295 END DO
296 END DO
297 END DO
298 END DO
299 IF ((abs(amin).ge.abs(aspval)).and. &
300 & (abs(amax).ge.abs(aspval))) THEN
301 amin=0.0_r8
302 amax=0.0_r8
303 END IF
304
305 IF (lchecksum) THEN
306 npts=(ubij-lbij+1)*(ubk-lbk+1)*nrec*4
307 IF (.not.allocated(cwrk)) allocate ( cwrk(npts) )
308 cwrk=pack(wrk(lbij:ubij, lbk:ubk, 1:4, 1:nrec), .true.)
309 CALL get_hash (cwrk, npts, checksum)
310 IF (allocated(cwrk)) deallocate (cwrk)
311 END IF
312 END IF
313 END IF
314
315#if !defined PARALLEL_IO && defined DISTRIBUTE
316
317 rbuffer(1)=real(status,r8)
318 rbuffer(2)=amin
319 rbuffer(3)=amax
320 CALL mp_bcastf (ng, model, rbuffer)
321 status=int(rbuffer(1))
322 amin=rbuffer(2)
323 amax=rbuffer(3)
324#endif
325
326 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
327 exit_flag=2
328 ioerror=status
329 RETURN
330 END IF
331
332#if !defined PARALLEL_IO && defined DISTRIBUTE
333
334
335
336 CALL mp_bcastf (ng, model, wrk)
337#endif
338
339
340
341
342
343 abry=0.0_r8
344
345 DO ir=1,nrec
346 DO ib=1,4
347 IF (bounded(ib)) THEN
348 IF ((ib.eq.iwest).or.(ib.eq.ieast)) THEN
349 DO k=lbk,ubk
350 DO j=jmin,jmax
351 abry(j,k,ib,ir)=wrk(j,k,ib,ir)
352 END DO
353 END DO
354 ELSE IF ((ib.eq.isouth).or.(ib.eq.inorth)) THEN
355 DO k=lbk,ubk
356 DO i=imin,imax
357 abry(i,k,ib,ir)=wrk(i,k,ib,ir)
358 END DO
359 END DO
360 END IF
361 END IF
362 END DO
363 END DO
364
365 RETURN