ROMS
Loading...
Searching...
No Matches
nf_fwrite3d_bry.F
Go to the documentation of this file.
1#include "cppdefs.h"
3
4#ifdef ADJUST_BOUNDARY
5!
6!git $Id$
7!================================================== Hernan G. Arango ===
8! Copyright (c) 2002-2025 The ROMS Group !
9! Licensed under a MIT/X style license !
10! See License_ROMS.md !
11!=======================================================================
12! !
13! This function writes out a generic floating point 3D boundary array !
14! into an output file using either the standard NetCDF library or the !
15! Parallel-IO (PIO) library. !
16! !
17! On Input: !
18! !
19! ng Nested grid number (integer) !
20! model Calling model identifier (integer) !
21! ncname NetCDF output file name (string) !
22! ncvname NetCDF variable name (string) !
23! ncid NetCDF file ID (integer) !
24# if defined PIO_LIB && defined DISTRIBUTE
25!or pioFile PIO file descriptor structure, TYPE(file_desc_t) !
26! pioFile%fh file handler !
27! pioFile%iosystem IO system descriptor (struct) !
28# endif
29! ncvarid NetCDF variable ID (integer) !
30# if defined PIO_LIB && defined DISTRIBUTE
31!or pioVar PIO variable descriptor structure, TYPE(My_VarDesc) !
32! pioVar%vd variable descriptor TYPE(Var_Desc_t)!
33! pioVar%dkind variable data kind !
34! pioVar%gtype variable C-gridtype !
35# endif
36! tindex NetCDF time record index to write (integer) !
37! gtype Grid type (integer) !
38# if defined PIO_LIB && defined DISTRIBUTE
39!or pioDesc IO data decomposition descriptor, TYPE(IO_desc_t) !
40# endif
41! LBij IJ-dimension lower bound (integer) !
42! UBij IJ-dimension upper bound (integer) !
43! LBk K-dimension lower bound (integer) !
44! UBk K-dimension upper bound (integer) !
45! Nrec Number of boundary records (integer) !
46! Ascl Factor to scale field before writing (real) !
47! Abry Boundary field to write out (real) !
48! ExtractField Field extraction flag (integer, OPTIONAL) !
49! ExtractField = 0 no extraction !
50! ExtractField = 1 extraction by interpolation !
51! ExtractField > 1 extraction by decimation !
52! !
53! On Output: !
54! !
55! status Error flag (integer). !
56! MinValue Minimum value (real, OPTIONAL) !
57! MaxValue Maximum value (real, OPTIONAL) !
58! !
59# ifdef POSITIVE_ZERO
60! Starting F95 zero values can be signed (-0 or +0) following the !
61! IEEE 754 floating point standard. This may produce different !
62! output data in serial and parallel applications. Since comparing !
63! serial and parallel output is essential for tracking parallel !
64! partition bugs, "positive zero" is enforced. !
65! !
66# endif
67!=======================================================================
68!
69 USE mod_param
70 USE mod_parallel
71 USE mod_ncparam
72 USE mod_scalars
73!
74# ifdef DISTRIBUTE
76# endif
78!
79 implicit none
80!
82 MODULE PROCEDURE nf90_fwrite3d_bry
83# if defined PIO_LIB && defined DISTRIBUTE
84 MODULE PROCEDURE pio_fwrite3d_bry
85# endif
86 END INTERFACE nf_fwrite3d_bry
87!
88 CONTAINS
89!
90!***********************************************************************
91 FUNCTION nf90_fwrite3d_bry (ng, model, ncname, ncid, &
92 & ncvname, ncvarid, &
93 & tindex, gtype, &
94 & LBij, UBij, LBk, UBk, Nrec, &
95 & Ascl, Abry, &
96 & ExtractField, &
97 & MinValue, MaxValue) RESULT(status)
98!***********************************************************************
99!
100 USE mod_netcdf
101!
102! Imported variable declarations.
103!
104 integer, intent(in) :: ng, model, ncid, ncvarid, tindex, gtype
105 integer, intent(in) :: lbij, ubij, lbk, ubk, nrec
106!
107 integer, intent(in), optional :: extractfield
108!
109 real(dp), intent(in) :: ascl
110!
111 character (len=*), intent(in) :: ncname
112 character (len=*), intent(in) :: ncvname
113!
114# ifdef ASSUMED_SHAPE
115 real(r8), intent(in) :: abry(lbij:,lbk:,:,:)
116# else
117 real(r8), intent(in) :: abry(lbij:ubij,lbk:ubk,4,nrec)
118# endif
119 real(r8), intent(out), optional :: minvalue
120 real(r8), intent(out), optional :: maxvalue
121!
122! Local variable declarations.
123!
124 integer :: iorj, ijklen, klen, npts, i, tile
125 integer :: extract_flag
126 integer :: status
127
128 integer, dimension(5) :: start, total
129!
130 real(r8), parameter :: aspv = 0.0_r8
131
132 real(r8), dimension((UBij-LBij+1)*(UBk-LBk+1)*4*Nrec) :: awrk
133!
134!-----------------------------------------------------------------------
135! Initialize local variables.
136!-----------------------------------------------------------------------
137!
138 status=nf90_noerr
139!
140! Set parallel tile.
141!
142# ifdef DISTRIBUTE
143 tile=myrank
144# else
145 tile=-1
146# endif
147!
148! Set length of boundary data as the value of greater of I- or
149! j-dimension.
150!
151 iorj=iobounds(ng)%IorJ
152 klen=ubk-lbk+1
153 ijklen=iorj*klen
154 npts=ijklen*4*nrec
155!
156! If appropriate, set the field extraction flag to the provided grid
157! geometry through interpolation or decimation.
158!
159 IF (PRESENT(extractfield)) THEN
160 extract_flag=extractfield
161 ELSE
162 extract_flag=0
163 END IF
164!
165! Initialize local array to avoid denormalized numbers. This
166! facilitates processing and debugging.
167!
168 awrk=0.0_r8
169!
170!-----------------------------------------------------------------------
171! Pack 2D boundary data into 1D array.
172!-----------------------------------------------------------------------
173!
174 CALL pack_boundary3d (ng, model, tile, &
175 & gtype, ncvname, tindex, &
176 & extract_flag, &
177 & lbij, ubij, lbk, ubk, nrec, &
178 & ascl, abry, &
179 & start, total, npts, awrk)
180!
181!-----------------------------------------------------------------------
182! If applicable, compute output field minimum and maximum values.
183!-----------------------------------------------------------------------
184!
185 IF (PRESENT(minvalue)) THEN
186 IF (outthread) THEN
187 minvalue=spval
188 maxvalue=-spval
189 DO i=1,npts
190 IF (abs(awrk(i)).lt.spval) THEN
191 minvalue=min(minvalue,awrk(i))
192 maxvalue=max(maxvalue,awrk(i))
193 END IF
194 END DO
195 END IF
196 END IF
197!
198!-----------------------------------------------------------------------
199! Write output buffer into NetCDF file.
200!-----------------------------------------------------------------------
201!
202 status=nf90_noerr
203 IF (outthread) THEN
204 status=nf90_put_var(ncid, ncvarid, awrk, start, total)
205 END IF
206
207# ifdef DISTRIBUTE
208!
209!-----------------------------------------------------------------------
210! Broadcast IO error flag to all nodes.
211!-----------------------------------------------------------------------
212!
213 CALL mp_bcasti (ng, model, status)
214# endif
215!
216 RETURN
217 END FUNCTION nf90_fwrite3d_bry
218
219# if defined PIO_LIB && defined DISTRIBUTE
220!
221!***********************************************************************
222 FUNCTION pio_fwrite3d_bry (ng, model, ncname, pioFile, &
223 & ncvname, pioVar, &
224 & tindex, pioDesc, &
225 & LBij, UBij, LBk, UBk, Nrec, &
226 & Ascl, Abry, &
227 & MinValue, MaxValue) RESULT(status)
228!***********************************************************************
229!
231!
232! Imported variable declarations.
233!
234 integer, intent(in) :: ng, model, tindex
235 integer, intent(in) :: lbij, ubij, lbk, ubk, nrec
236!
237 real(dp), intent(in) :: ascl
238!
239 character (len=*), intent(in) :: ncname
240 character (len=*), intent(in) :: ncvname
241!
242# ifdef ASSUMED_SHAPE
243 real(r8), intent(in) :: abry(lbij:,lbk:,:,:)
244# else
245 real(r8), intent(in) :: abry(lbij:ubij,lbk:ubk,4,nrec)
246# endif
247 real(r8), intent(out), optional :: minvalue
248 real(r8), intent(out), optional :: maxvalue
249!
250 TYPE (file_desc_t), intent(inout) :: piofile
251 TYPE (io_desc_t), intent(inout) :: piodesc
252 TYPE (my_vardesc), intent(inout) :: piovar
253!
254! Local variable declarations.
255!
256 logical, dimension(4) :: bounded
257!
258 integer :: bc, i, ib, ic, ir, j, k, kc, rc
259 integer :: dkind, gtype, tile
260 integer :: iorj, ijklen, imin, imax, jmin, jmax, klen, npts
261 integer :: istr, iend, jstr, jend
262
263 integer, dimension(5) :: start, total
264
265 integer :: status
266!
267 real(r8), parameter :: aspv = 0.0_r8
268
269 real(r8), dimension((UBij-LBij+1)*(UBk-LBk+1)*4*Nrec) :: awrk
270!
271!-----------------------------------------------------------------------
272! Set starting and ending indices to process.
273!-----------------------------------------------------------------------
274!
275 status=pio_noerr
276!
277! Set tile starting and ending bounds.
278!
279 tile=myrank
280 dkind=piovar%dkind
281 gtype=piovar%gtype
282!
283 SELECT CASE (gtype)
284 CASE (p2dvar, p3dvar)
285 imin=bounds(ng)%Istr (tile)
286 imax=bounds(ng)%Iend (tile)
287 jmin=bounds(ng)%Jstr (tile)
288 jmax=bounds(ng)%Jend (tile)
289 CASE (r2dvar, r3dvar)
290 imin=bounds(ng)%IstrR(tile)
291 imax=bounds(ng)%IendR(tile)
292 jmin=bounds(ng)%JstrR(tile)
293 jmax=bounds(ng)%JendR(tile)
294 CASE (u2dvar, u3dvar)
295 imin=bounds(ng)%Istr (tile)
296 imax=bounds(ng)%IendR(tile)
297 jmin=bounds(ng)%JstrR(tile)
298 jmax=bounds(ng)%JendR(tile)
299 CASE (v2dvar, v3dvar)
300 imin=bounds(ng)%IstrR(tile)
301 imax=bounds(ng)%IendR(tile)
302 jmin=bounds(ng)%Jstr (tile)
303 jmax=bounds(ng)%JendR(tile)
304 CASE DEFAULT
305 imin=bounds(ng)%IstrR(tile)
306 imax=bounds(ng)%IendR(tile)
307 jmin=bounds(ng)%JstrR(tile)
308 jmax=bounds(ng)%JendR(tile)
309 END SELECT
310!
311 iorj=iobounds(ng)%IorJ
312 klen=ubk-lbk+1
313 ijklen=iorj*klen
314 npts=ijklen*4*nrec
315!
316! Get tile bounds.
317!
318 istr=bounds(ng)%Istr(tile)
319 iend=bounds(ng)%Iend(tile)
320 jstr=bounds(ng)%Jstr(tile)
321 jend=bounds(ng)%Jend(tile)
322!
323! Set switch to process boundary data by their associated tiles.
324!
325 bounded(iwest )=domain(ng)%Western_Edge (tile)
326 bounded(ieast )=domain(ng)%Eastern_Edge (tile)
327 bounded(isouth)=domain(ng)%Southern_Edge(tile)
328 bounded(inorth)=domain(ng)%Northern_Edge(tile)
329!
330! Set NetCDF dimension counters for processing requested field.
331!
332 start(1)=1
333 total(1)=iorj
334 start(2)=1
335 total(2)=klen
336 start(3)=1
337 total(3)=4
338 start(4)=1
339 total(4)=nrec
340 start(5)=tindex
341 total(5)=1
342!
343!-----------------------------------------------------------------------
344! Pack and scale output data.
345!-----------------------------------------------------------------------
346!
347 awrk=aspv
348!
349 DO ir=1,nrec
350 rc=(ir-1)*ijklen*4
351 DO ib=1,4
352 IF (bounded(ib)) THEN
353 bc=(ib-1)*ijklen+rc
354 IF ((ib.eq.iwest).or.(ib.eq.ieast)) THEN
355 DO k=lbk,ubk
356 kc=(k-lbk)*iorj+bc
357 DO j=jmin,jmax
358 ic=1+(j-lbij)+kc
359 awrk(ic)=abry(j,k,ib,ir)*ascl
360# ifdef POSITIVE_ZERO
361 IF (abs(awrk(ic)).eq.0.0_r8) THEN
362 awrk(ic)=0.0_r8 ! impose positive zero
363 END IF
364# endif
365 END DO
366 END DO
367 ELSE IF ((ib.eq.isouth).or.(ib.eq.inorth)) THEN
368 DO k=lbk,ubk
369 kc=(k-lbk)*iorj+bc
370 DO i=imin,imax
371 ic=1+(i-lbij)+kc
372 awrk(ic)=abry(i,k,ib,ir)*ascl
373# ifdef POSITIVE_ZERO
374 IF (abs(awrk(ic)).eq.0.0_r8) THEN
375 awrk(ic)=0.0_r8 ! impose positive zero
376 END IF
377# endif
378 END DO
379 END DO
380 END IF
381 END IF
382 END DO
383 END DO
384!
385! Collect data from all spawned processes.
386!
387 CALL mp_collect (ng, model, npts, aspv, awrk)
388!
389!-----------------------------------------------------------------------
390! If applicable, compute output field minimum and maximum values.
391!-----------------------------------------------------------------------
392!
393 IF (PRESENT(minvalue)) THEN
394 minvalue=spval
395 maxvalue=-spval
396 DO i=1,npts
397 IF (abs(awrk(i)).lt.spval) THEN
398 minvalue=min(minvalue,awrk(i))
399 maxvalue=max(maxvalue,awrk(i))
400 END IF
401 END DO
402 END IF
403!
404!-----------------------------------------------------------------------
405! Write output buffer into NetCDF file.
406!-----------------------------------------------------------------------
407!
408 status=pio_put_var(piofile, piovar%vd, start, total, awrk)
409!
410 RETURN
411 END FUNCTION pio_fwrite3d_bry
412# endif
413#endif
414 END MODULE nf_fwrite3d_bry_mod
logical outthread
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
integer, parameter r3dvar
Definition mod_param.F:721
type(t_iobounds), dimension(:), allocatable iobounds
Definition mod_param.F:282
integer, parameter u3dvar
Definition mod_param.F:722
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, parameter u2dvar
Definition mod_param.F:718
integer, parameter p2dvar
Definition mod_param.F:716
integer, parameter r2dvar
Definition mod_param.F:717
integer, parameter v2dvar
Definition mod_param.F:719
integer, parameter p3dvar
Definition mod_param.F:720
integer, parameter v3dvar
Definition mod_param.F:723
real(dp), parameter spval
integer, parameter iwest
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth
integer function pio_fwrite3d_bry(ng, model, ncname, piofile, ncvname, piovar, tindex, piodesc, lbij, ubij, lbk, ubk, nrec, ascl, abry, minvalue, maxvalue)
integer function nf90_fwrite3d_bry(ng, model, ncname, ncid, ncvname, ncvarid, tindex, gtype, lbij, ubij, lbk, ubk, nrec, ascl, abry, extractfield, minvalue, maxvalue)
subroutine, public pack_boundary3d(ng, model, tile, gtype, ncvname, tindex, extract_flag, lbij, ubij, lbk, ubk, nrec, bscl, bdat, start, total, npts, bwrk)
Definition pack_field.F:277