ROMS
Loading...
Searching...
No Matches
nf_fwrite2d_bry_mod Module Reference

Data Types

interface  nf_fwrite2d_bry
 

Functions/Subroutines

integer function nf90_fwrite2d_bry (ng, model, ncname, ncid, ncvname, ncvarid, tindex, gtype, lbij, ubij, nrec, ascl, abry, extractfield, minvalue, maxvalue)
 
integer function pio_fwrite2d_bry (ng, model, ncname, piofile, ncvname, piovar, tindex, piodesc, lbij, ubij, nrec, ascl, abry, minvalue, maxvalue)
 

Function/Subroutine Documentation

◆ nf90_fwrite2d_bry()

integer function nf_fwrite2d_bry_mod::nf90_fwrite2d_bry ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
integer, intent(in) ncid,
character (len=*), intent(in) ncvname,
integer, intent(in) ncvarid,
integer, intent(in) tindex,
integer, intent(in) gtype,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) nrec,
real(dp), intent(in) ascl,
real(r8), dimension(lbij:,:,:), intent(in) abry,
integer, intent(in), optional extractfield,
real(r8), intent(out), optional minvalue,
real(r8), intent(out), optional maxvalue )

Definition at line 89 of file nf_fwrite2d_bry.F.

96!***********************************************************************
97!
98 USE mod_netcdf
99!
100! Imported variable declarations.
101!
102 integer, intent(in) :: ng, model, ncid, ncvarid, tindex, gtype
103 integer, intent(in) :: LBij, UBij, Nrec
104!
105 integer, intent(in), optional :: ExtractField
106!
107 real(dp), intent(in) :: Ascl
108!
109 character (len=*), intent(in) :: ncname
110 character (len=*), intent(in) :: ncvname
111!
112# ifdef ASSUMED_SHAPE
113 real(r8), intent(in) :: Abry(LBij:,:,:)
114# else
115 real(r8), intent(in) :: Abry(LBij:UBij,4,Nrec)
116# endif
117 real(r8), intent(out), optional :: MinValue
118 real(r8), intent(out), optional :: MaxValue
119!
120! Local variable declarations.
121!
122 integer :: IorJ, Npts, i, tile
123 integer :: Extract_Flag
124 integer :: status
125!
126 integer, dimension(4) :: start, total
127!
128 real(r8), dimension((UBij-LBij+1)*4*Nrec) :: Awrk
129!
130!-----------------------------------------------------------------------
131! Initialize local variables.
132!-----------------------------------------------------------------------
133!
134 status=nf90_noerr
135!
136! Set parallel tile.
137!
138# ifdef DISTRIBUTE
139 tile=myrank
140# else
141 tile=-1
142# endif
143!
144! Set length of boundary data as the value of greater of I- or
145! j-dimension.
146!
147 iorj=iobounds(ng)%IorJ
148 npts=iorj*4*nrec
149!
150! If appropriate, set the field extraction flag to the provided grid
151! geometry through interpolation or decimation.
152!
153 IF (PRESENT(extractfield)) THEN
154 extract_flag=extractfield
155 ELSE
156 extract_flag=0
157 END IF
158!
159! Initialize local array to avoid denormalized numbers. This
160! facilitates processing and debugging.
161!
162 awrk=0.0_r8
163!
164!-----------------------------------------------------------------------
165! Pack 2D boundary data into 1D array.
166!-----------------------------------------------------------------------
167!
168 CALL pack_boundary2d (ng, model, tile, &
169 & gtype, ncvname, tindex, &
170 & extract_flag, &
171 & lbij, ubij, nrec, &
172 & ascl, abry, &
173 & start, total, npts, awrk)
174!
175!-----------------------------------------------------------------------
176! If applicable, compute output field minimum and maximum values.
177!-----------------------------------------------------------------------
178!
179 IF (PRESENT(minvalue)) THEN
180 IF (outthread) THEN
181 minvalue=spval
182 maxvalue=-spval
183 DO i=1,npts
184 IF (abs(awrk(i)).lt.spval) THEN
185 minvalue=min(minvalue,awrk(i))
186 maxvalue=max(maxvalue,awrk(i))
187 END IF
188 END DO
189 END IF
190 END IF
191!
192!-----------------------------------------------------------------------
193! Write output buffer into NetCDF file.
194!-----------------------------------------------------------------------
195!
196 status=nf90_noerr
197 IF (outthread) THEN
198 status=nf90_put_var(ncid, ncvarid, awrk, start, total)
199 END IF
200
201# ifdef DISTRIBUTE
202!
203!-----------------------------------------------------------------------
204! Broadcast IO error flag to all nodes.
205!-----------------------------------------------------------------------
206!
207 CALL mp_bcasti (ng, model, status)
208# endif
209!
210 RETURN

◆ pio_fwrite2d_bry()

integer function nf_fwrite2d_bry_mod::pio_fwrite2d_bry ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
type (file_desc_t), intent(inout) piofile,
character (len=*), intent(in) ncvname,
type (my_vardesc), intent(inout) piovar,
integer, intent(in) tindex,
type (io_desc_t), intent(inout) piodesc,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) nrec,
real(dp), intent(in) ascl,
real(r8), dimension(lbij:,:,:), intent(in) abry,
real(r8), intent(out), optional minvalue,
real(r8), intent(out), optional maxvalue )

Definition at line 216 of file nf_fwrite2d_bry.F.

222!***********************************************************************
223!
225!
226! Imported variable declarations.
227!
228 integer, intent(in) :: ng, model, tindex
229 integer, intent(in) :: LBij, UBij, Nrec
230!
231 real(dp), intent(in) :: Ascl
232!
233 character (len=*), intent(in) :: ncname
234 character (len=*), intent(in) :: ncvname
235!
236# ifdef ASSUMED_SHAPE
237 real(r8), intent(in) :: Abry(LBij:,:,:)
238# else
239 real(r8), intent(in) :: Abry(LBij:UBij,4,Nrec)
240# endif
241 real(r8), intent(out), optional :: MinValue
242 real(r8), intent(out), optional :: MaxValue
243!
244 TYPE (File_desc_t), intent(inout) :: pioFile
245 TYPE (IO_Desc_t), intent(inout) :: pioDesc
246 TYPE (My_VarDesc), intent(inout) :: pioVar
247!
248! Local variable declarations.
249!
250 logical, dimension(4) :: bounded
251!
252 integer :: bc, i, ib, ic, ir, j, rc
253 integer :: dkind, gtype, tile
254 integer :: IorJ, Imin, Imax, Jmin, Jmax, Npts
255 integer :: Istr, Iend, Jstr, Jend
256
257 integer, dimension(4) :: start, total
258
259 integer :: status
260!
261 real(r8), parameter :: Aspv = 0.0_r8
262
263 real(r8), dimension((UBij-LBij+1)*4*Nrec) :: Awrk
264!
265!-----------------------------------------------------------------------
266! Set starting and ending indices to process.
267!-----------------------------------------------------------------------
268!
269 status=pio_noerr
270!
271! Set tile starting and ending bounds.
272!
273 tile=myrank
274 dkind=piovar%dkind
275 gtype=piovar%gtype
276!
277 SELECT CASE (gtype)
278 CASE (p2dvar, p3dvar)
279 imin=bounds(ng)%Istr (tile)
280 imax=bounds(ng)%Iend (tile)
281 jmin=bounds(ng)%Jstr (tile)
282 jmax=bounds(ng)%Jend (tile)
283 CASE (r2dvar, r3dvar)
284 imin=bounds(ng)%IstrR(tile)
285 imax=bounds(ng)%IendR(tile)
286 jmin=bounds(ng)%JstrR(tile)
287 jmax=bounds(ng)%JendR(tile)
288 CASE (u2dvar, u3dvar)
289 imin=bounds(ng)%Istr (tile)
290 imax=bounds(ng)%IendR(tile)
291 jmin=bounds(ng)%JstrR(tile)
292 jmax=bounds(ng)%JendR(tile)
293 CASE (v2dvar, v3dvar)
294 imin=bounds(ng)%IstrR(tile)
295 imax=bounds(ng)%IendR(tile)
296 jmin=bounds(ng)%Jstr (tile)
297 jmax=bounds(ng)%JendR(tile)
298 CASE DEFAULT
299 imin=bounds(ng)%IstrR(tile)
300 imax=bounds(ng)%IendR(tile)
301 jmin=bounds(ng)%JstrR(tile)
302 jmax=bounds(ng)%JendR(tile)
303 END SELECT
304!
305 iorj=iobounds(ng)%IorJ
306 npts=iorj*4*nrec
307!
308! Get tile bounds.
309!
310 istr=bounds(ng)%Istr(tile)
311 iend=bounds(ng)%Iend(tile)
312 jstr=bounds(ng)%Jstr(tile)
313 jend=bounds(ng)%Jend(tile)
314!
315! Set switch to process boundary data by their associated tiles.
316!
317 bounded(iwest )=domain(ng)%Western_Edge (tile)
318 bounded(ieast )=domain(ng)%Eastern_Edge (tile)
319 bounded(isouth)=domain(ng)%Southern_Edge(tile)
320 bounded(inorth)=domain(ng)%Northern_Edge(tile)
321!
322! Set NetCDF dimension counters for processing requested field.
323!
324 start(1)=1
325 total(1)=iorj
326 start(2)=1
327 total(2)=4
328 start(3)=1
329 total(3)=nrec
330 start(4)=tindex
331 total(4)=1
332!
333!-----------------------------------------------------------------------
334! Pack and scale output data.
335!-----------------------------------------------------------------------
336!
337 awrk=aspv
338!
339 DO ir=1,nrec
340 rc=(ir-1)*iorj*4
341 DO ib=1,4
342 IF (bounded(ib)) THEN
343 bc=(ib-1)*iorj+rc
344 IF ((ib.eq.iwest).or.(ib.eq.ieast)) THEN
345 DO j=jmin,jmax
346 ic=1+(j-lbij)+bc
347 awrk(ic)=abry(j,ib,ir)*ascl
348# ifdef POSITIVE_ZERO
349 IF (abs(awrk(ic)).eq.0.0_r8) THEN
350 awrk(ic)=0.0_r8 ! impose positive zero
351 END IF
352# endif
353 END DO
354 ELSE IF ((ib.eq.isouth).or.(ib.eq.inorth)) THEN
355 DO i=imin,imax
356 ic=1+(i-lbij)+bc
357 awrk(ic)=abry(i,ib,ir)*ascl
358# ifdef POSITIVE_ZERO
359 IF (abs(awrk(ic)).eq.0.0_r8) THEN
360 awrk(ic)=0.0_r8 ! impose positive zero
361 END IF
362# endif
363 END DO
364 END IF
365 END IF
366 END DO
367 END DO
368!
369! Collect data from all spawned processes.
370!
371 CALL mp_collect (ng, model, npts, aspv, awrk)
372!
373!-----------------------------------------------------------------------
374! If applicable, compute output field minimum and maximum values.
375!-----------------------------------------------------------------------
376!
377 IF (PRESENT(minvalue)) THEN
378 minvalue=spval
379 maxvalue=-spval
380 DO i=1,npts
381 IF (abs(awrk(i)).lt.spval) THEN
382 minvalue=min(minvalue,awrk(i))
383 maxvalue=max(maxvalue,awrk(i))
384 END IF
385 END DO
386 END IF
387!
388!-----------------------------------------------------------------------
389! Write output buffer into NetCDF file.
390!-----------------------------------------------------------------------
391!
392 status=pio_put_var(piofile, piovar%vd, start, total, awrk)
393!
394 RETURN