ROMS
Loading...
Searching...
No Matches
rp_set_data.F File Reference
#include "cppdefs.h"
#include "tile.h"
#include "set_bounds.h"
Include dependency graph for rp_set_data.F:

Go to the source code of this file.

Functions/Subroutines

subroutine rp_set_data (ng, tile)
 
subroutine rp_set_data_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs)
 

Function/Subroutine Documentation

◆ rp_set_data()

subroutine rp_set_data ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 3 of file rp_set_data.F.

4!
5!git $Id$
6!================================================== Hernan G. Arango ===
7! Copyright (c) 2002-2025 The ROMS Group !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md !
10!=======================================================================
11! !
12! This subroutine processes forcing, boundary, climatology, and !
13! other input data. It time-interpolates between snapshots. !
14! !
15!=======================================================================
16!
17 USE mod_param
18!
19! Imported variable declarations.
20!
21 integer, intent(in) :: ng, tile
22!
23! Local variable declarations.
24!
25 character (len=*), parameter :: MyFile = &
26 & __FILE__
27!
28# include "tile.h"
29!
30# ifdef PROFILE
31 CALL wclock_on (ng, irpm, 4, __line__, myfile)
32# endif
33 CALL rp_set_data_tile (ng, tile, &
34 & lbi, ubi, lbj, ubj, &
35 & imins, imaxs, jmins, jmaxs)
36# ifdef PROFILE
37 CALL wclock_off (ng, irpm, 4, __line__, myfile)
38# endif
39!
40 RETURN
integer, parameter irpm
Definition mod_param.F:664
subroutine rp_set_data_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs)
Definition rp_set_data.F:47
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3

References mod_param::irpm, rp_set_data_tile(), wclock_off(), and wclock_on().

Referenced by rp_main3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ rp_set_data_tile()

subroutine rp_set_data_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs )

Definition at line 44 of file rp_set_data.F.

47!***********************************************************************
48!
49 USE mod_param
50 USE mod_boundary
51 USE mod_clima
52# if defined FORWARD_READ && defined SOLVE3D
53 USE mod_coupling
54# endif
55 USE mod_forces
56 USE mod_grid
57 USE mod_mixing
58 USE mod_ncparam
59 USE mod_ocean
60 USE mod_stepping
61 USE mod_scalars
62 USE mod_sources
63!
64# ifdef ANALYTICAL
66# endif
69# ifdef SOLVE3D
71# endif
72# ifdef DISTRIBUTE
73# if defined WET_DRY
74 USE distribute_mod, ONLY : mp_boundary
75# endif
77# ifdef SOLVE3D
79# endif
80# endif
81 USE strings_mod, ONLY : founderror
82!
83 implicit none
84!
85! Imported variable declarations.
86!
87 integer, intent(in) :: ng, tile
88 integer, intent(in) :: LBi, UBi, LBj, UBj
89 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
90!
91! Local variable declarations.
92!
93 logical :: SetBC
94 logical :: update = .false.
95# if defined WET_DRY
96 logical :: bry_update
97# endif
98!
99 integer :: ILB, IUB, JLB, JUB
100 integer :: i, ic, itrc, j, k, my_tile
101!
102 real(r8) :: cff, cff1, cff2
103!
104 character (len=*), parameter :: MyFile = &
105 & __FILE__//", rp_set_data_tile"
106
107# include "set_bounds.h"
108!
109! Lower and upper bounds for nontiled (global values) boundary arrays.
110!
111 my_tile=-1 ! for global values
112 ilb=bounds(ng)%LBi(my_tile)
113 iub=bounds(ng)%UBi(my_tile)
114 jlb=bounds(ng)%LBj(my_tile)
115 jub=bounds(ng)%UBj(my_tile)
116
117# ifdef SOLVE3D
118
119# ifdef CLOUDS
120!
121!-----------------------------------------------------------------------
122! Set cloud fraction (nondimensional). Notice that clouds are
123! processed first in case that they are used to adjust shortwave
124! radiation.
125!-----------------------------------------------------------------------
126!
127# ifdef ANA_CLOUD
128 CALL ana_cloud (ng, tile, irpm)
129# else
130 CALL set_2dfld_tile (ng, tile, irpm, idcfra, &
131 & lbi, ubi, lbj, ubj, &
132 & forces(ng)%cloudG, &
133 & forces(ng)%cloud, &
134 & update)
135 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
136# endif
137# endif
138
139# if (defined BULK_FLUXES && !defined FORWARD_FLUXES) || \
140 defined ecosim || \
141 (defined shortwave && defined ana_srflux && defined albedo)
142!
143!-----------------------------------------------------------------------
144! Set surface air temperature (degC).
145!-----------------------------------------------------------------------
146!
147# ifdef ANA_TAIR
148 CALL ana_tair (ng, tile, irpm)
149# else
150 CALL set_2dfld_tile (ng, tile, irpm, idtair, &
151 & lbi, ubi, lbj, ubj, &
152 & forces(ng)%TairG, &
153 & forces(ng)%Tair, &
154 & update)
155 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
156# endif
157# endif
158
159# if (defined BULK_FLUXES && !defined FORWARD_FLUXES) || \
160 defined ecosim || \
161 (defined shortwave && defined ana_srflux && defined albedo)
162!
163!-----------------------------------------------------------------------
164! Set surface air relative or specific humidity.
165!-----------------------------------------------------------------------
166!
167# ifdef ANA_HUMIDITY
168 CALL ana_humid (ng, tile, irpm)
169# else
170 CALL set_2dfld_tile (ng, tile, irpm, idqair, &
171 & lbi, ubi, lbj, ubj, &
172 & forces(ng)%HairG, &
173 & forces(ng)%Hair, &
174 & update)
175 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
176# endif
177# endif
178
179# ifdef SHORTWAVE
180!
181!-----------------------------------------------------------------------
182! Set kinematic surface solar shortwave radiation flux (degC m/s).
183!-----------------------------------------------------------------------
184!
185# ifdef ANA_SRFLUX
186 CALL ana_srflux (ng, tile, irpm)
187# else
188 CALL set_2dfld_tile (ng, tile, irpm, idsrad, &
189 & lbi, ubi, lbj, ubj, &
190 & forces(ng)%srflxG, &
191 & forces(ng)%srflx, &
192 & update)
193 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
194# endif
195
196# if defined DIURNAL_SRFLUX && !defined FORWARD_FLUXES
197!
198! Modulate the averaged shortwave radiation flux by the local diurnal
199! cycle.
200!
201 CALL ana_srflux (ng, tile, irpm)
202# endif
203# endif
204
205# if (defined BULK_FLUXES && !defined FORWARD_FLUXES) && \
206 !defined LONGWAVE && !defined LONGWAVE_OUT
207!
208!-----------------------------------------------------------------------
209! Surface net longwave radiation (degC m/s).
210!-----------------------------------------------------------------------
211!
212 CALL set_2dfld_tile (ng, tile, irpm, idlrad, &
213 & lbi, ubi, lbj, ubj, &
214 & forces(ng)%lrflxG, &
215 & forces(ng)%lrflx, &
216 & update)
217 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
218# endif
219
220# if defined LONGWAVE_OUT && \
221 (defined bulk_fluxes && !defined FORWARD_FLUXES)
222!
223!-----------------------------------------------------------------------
224! Surface downwelling longwave radiation (degC m/s).
225!-----------------------------------------------------------------------
226!
227 CALL set_2dfld_tile (ng, tile, irpm, idldwn, &
228 & lbi, ubi, lbj, ubj, &
229 & forces(ng)%lrflxG, &
230 & forces(ng)%lrflx, &
231 & update)
232 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
233# endif
234
235# if (defined BULK_FLUXES && !defined FORWARD_FLUXES) || \
236 defined ecosim
237!
238!-----------------------------------------------------------------------
239! Set surface winds (m/s).
240!-----------------------------------------------------------------------
241!
242# ifdef ANA_WINDS
243 CALL ana_winds (ng, tile, irpm)
244# else
245 CALL set_2dfld_tile (ng, tile, irpm, iduair, &
246 & lbi, ubi, lbj, ubj, &
247 & forces(ng)%UwindG, &
248 & forces(ng)%Uwind, &
249 & update)
250 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
251
252 CALL set_2dfld_tile (ng, tile, irpm, idvair, &
253 & lbi, ubi, lbj, ubj, &
254 & forces(ng)%VwindG, &
255 & forces(ng)%Vwind, &
256 & update)
257 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
258
259# ifdef CURVGRID
260!
261! If input point surface winds or interpolated from coarse data, rotate
262! to curvilinear grid.
263!
264 IF (.not.linfo(1,iduair,ng).or. &
265 & (iinfo(5,iduair,ng).ne.lm(ng)+2).or. &
266 & (iinfo(6,iduair,ng).ne.mm(ng)+2)) THEN
267 DO j=jstrr,jendr
268 DO i=istrr,iendr
269 cff1=forces(ng)%Uwind(i,j)*grid(ng)%CosAngler(i,j)+ &
270 & forces(ng)%Vwind(i,j)*grid(ng)%SinAngler(i,j)
271 cff2=forces(ng)%Vwind(i,j)*grid(ng)%CosAngler(i,j)- &
272 & forces(ng)%Uwind(i,j)*grid(ng)%SinAngler(i,j)
273 forces(ng)%Uwind(i,j)=cff1
274 forces(ng)%Vwind(i,j)=cff2
275 END DO
276 END DO
277
278 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
279 CALL exchange_r2d_tile (ng, tile, &
280 & lbi, ubi, lbj, ubj, &
281 & forces(ng)%UWind)
282 CALL exchange_r2d_tile (ng, tile, &
283 & lbi, ubi, lbj, ubj, &
284 & forces(ng)%VWind)
285 END IF
286
287# ifdef DISTRIBUTE
288 CALL mp_exchange2d (ng, tile, irpm, 2, &
289 & lbi, ubi, lbj, ubj, &
290 & nghostpoints, &
291 & ewperiodic(ng), nsperiodic(ng), &
292 & forces(ng)%UWind, &
293 & forces(ng)%VWind)
294# endif
295 END IF
296# endif
297# endif
298# endif
299
300# if defined BULK_FLUXES && !defined FORWARD_FLUXES
301!
302!-----------------------------------------------------------------------
303! Set rain fall rate (kg/m2/s).
304!-----------------------------------------------------------------------
305!
306# ifdef ANA_RAIN
307 CALL ana_rain (ng, tile, irpm)
308# else
309 CALL set_2dfld_tile (ng, tile, irpm, idrain, &
310 & lbi, ubi, lbj, ubj, &
311 & forces(ng)%rainG, &
312 & forces(ng)%rain, &
313 & update)
314 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
315# endif
316# endif
317
318# if !defined BULK_FLUXES || defined FORWARD_FLUXES
319!
320!-----------------------------------------------------------------------
321! Set kinematic surface net heat flux (degC m/s).
322!-----------------------------------------------------------------------
323!
324# ifdef ANA_STFLUX
325 CALL ana_stflux (ng, tile, irpm, itemp)
326# else
327 CALL set_2dfld_tile (ng, tile, irpm, idtsur(itemp), &
328 & lbi, ubi, lbj, ubj, &
329 & forces(ng)%stfluxG(:,:,:,itemp), &
330 & forces(ng)%stflux (:,:,itemp), &
331 & update)
332 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
333
334# ifdef TL_IOMS
335 CALL set_2dfld_tile (ng, tile, irpm, idtsur(itemp), &
336 & lbi, ubi, lbj, ubj, &
337 & forces(ng)%stfluxG(:,:,:,itemp), &
338 & forces(ng)%tl_stflux (:,:,itemp), &
339 & update)
340 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
341# endif
342# endif
343# endif
344
345# ifdef QCORRECTION
346!
347!-----------------------------------------------------------------------
348! Set sea surface temperature (SST) and heat flux sensitivity to
349! SST (dQdSST) which are used for surface heat flux correction.
350!-----------------------------------------------------------------------
351!
352# ifdef ANA_SST
353 CALL ana_sst (ng, tile, irpm)
354# else
355 CALL set_2dfld_tile (ng, tile, irpm, idsstc, &
356 & lbi, ubi, lbj, ubj, &
357 & forces(ng)%sstG, &
358 & forces(ng)%sst, &
359 & update)
360 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
361# endif
362!
363# ifdef ANA_DQDSST
364 CALL ana_dqdsst (ng, tile, irpm)
365# else
366 CALL set_2dfld_tile (ng, tile, irpm, iddqdt, &
367 & lbi, ubi, lbj, ubj, &
368 & forces(ng)%dqdtG, &
369 & forces(ng)%dqdt, &
370 & update)
371 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
372# endif
373# endif
374!
375!-----------------------------------------------------------------------
376! Set kinematic bottom net heat flux (degC m/s).
377!-----------------------------------------------------------------------
378!
379# ifdef ANA_BTFLUX
380 CALL ana_btflux (ng, tile, irpm, itemp)
381# else
382 CALL set_2dfld_tile (ng, tile, irpm, idtbot(itemp), &
383 & lbi, ubi, lbj, ubj, &
384 & forces(ng)%btfluxG(:,:,:,itemp), &
385 & forces(ng)%btflux (:,:,itemp), &
386 & update)
387 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
388
389# ifdef TL_IOMS
390 CALL set_2dfld_tile (ng, tile, irpm, idtbot(itemp), &
391 & lbi, ubi, lbj, ubj, &
392 & forces(ng)%btfluxG(:,:,:,itemp), &
393 & forces(ng)%tl_btflux(:,:,itemp), &
394 & update)
395 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
396# endif
397# endif
398
399# ifdef SALINITY
400!
401!-----------------------------------------------------------------------
402! Set kinematic surface freshwater (E-P) flux (m/s).
403!-----------------------------------------------------------------------
404!
405# ifdef ANA_SSFLUX
406 CALL ana_stflux (ng, tile, irpm, isalt)
407# else
408# if !(defined EMINUSP || defined FORWARD_FLUXES || \
409 defined frc_coupling || defined srelaxation)
410 CALL set_2dfld_tile (ng, tile, irpm, idsfwf, &
411 & lbi, ubi, lbj, ubj, &
412 & forces(ng)%stfluxG(:,:,:,isalt), &
413 & forces(ng)%stflux (:,:,isalt), &
414 & update)
415 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
416
417# ifdef TL_IOMS
418 CALL set_2dfld_tile (ng, tile, irpm, idsfwf, &
419 & lbi, ubi, lbj, ubj, &
420 & forces(ng)%stfluxG(:,:,:,isalt), &
421 & forces(ng)%tl_stflux(:,:,isalt), &
422 & update)
423 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
424# endif
425
426# elif (defined EMINUSP || defined FORWARD_FLUXES || \
427 defined frc_coupling)
428 CALL set_2dfld_tile (ng, tile, irpm, idempf, &
429 & lbi, ubi, lbj, ubj, &
430 & forces(ng)%stfluxG(:,:,:,isalt), &
431 & forces(ng)%stflux (:,:,isalt), &
432 & update)
433 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
434
435# ifdef TL_IOMS
436 CALL set_2dfld_tile (ng, tile, irpm, idempf, &
437 & lbi, ubi, lbj, ubj, &
438 & forces(ng)%stfluxG(:,:,:,isalt), &
439 & forces(ng)%tl_stflux (:,:,isalt), &
440 & update)
441 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
442# endif
443# endif
444# endif
445
446# if defined SCORRECTION || defined SRELAXATION
447!
448!-----------------------------------------------------------------------
449! Set surface salinity for freshwater flux correction.
450!-----------------------------------------------------------------------
451!
452# ifdef ANA_SSS
453 CALL ana_sss (ng, tile, irpm)
454# else
455 CALL set_2dfld_tile (ng, tile, irpm, idsssc, &
456 & lbi, ubi, lbj, ubj, &
457 & forces(ng)%sssG, &
458 & forces(ng)%sss, &
459 & update)
460 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
461# endif
462# endif
463!
464!-----------------------------------------------------------------------
465! Set kinematic bottom salt flux (m/s).
466!-----------------------------------------------------------------------
467!
468# ifdef ANA_BSFLUX
469 CALL ana_btflux (ng, tile, irpm, isalt)
470# else
471 CALL set_2dfld_tile (ng, tile, irpm, idtbot(isalt), &
472 & lbi, ubi, lbj, ubj, &
473 & forces(ng)%btfluxG(:,:,:,isalt), &
474 & forces(ng)%btflux (:,:,isalt), &
475 & update)
476 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
477
478# ifdef TL_IOMS
479 CALL set_2dfld_tile (ng, tile, irpm, idtbot(isalt), &
480 & lbi, ubi, lbj, ubj, &
481 & forces(ng)%btfluxG(:,:,:,isalt), &
482 & forces(ng)%tl_btflux(:,:,isalt), &
483 & update)
484 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
485# endif
486# endif
487# endif
488
489# if defined BIOLOGY || defined SEDIMENT || defined T_PASSIVE
490!
491!-----------------------------------------------------------------------
492! Set kinematic surface and bottom passive tracer fluxes (T m/s).
493!-----------------------------------------------------------------------
494!
495 DO itrc=nat+1,nt(ng)
496# ifdef ANA_SPFLUX
497 CALL ana_stflux (ng, tile, irpm, itrc)
498# else
499 CALL set_2dfld_tile (ng, tile, irpm, idtsur(itrc), &
500 & lbi, ubi, lbj, ubj, &
501 & forces(ng)%stfluxG(:,:,:,itrc), &
502 & forces(ng)%stflux (:,:,itrc), &
503 & update)
504 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
505
506# ifdef TL_IOMS
507 CALL set_2dfld_tile (ng, tile, irpm, idtsur(itrc), &
508 & lbi, ubi, lbj, ubj, &
509 & forces(ng)%stfluxG(:,:,:,itrc), &
510 & forces(ng)%tl_stflux (:,:,itrc), &
511 & update)
512 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
513# endif
514# endif
515
516# ifdef ANA_BPFLUX
517 CALL ana_btflux (ng, tile, irpm, itrc)
518# else
519 CALL set_2dfld_tile (ng, tile, irpm, idtbot(itrc), &
520 & lbi, ubi, lbj, ubj, &
521 & forces(ng)%btfluxG(:,:,:,itrc), &
522 & forces(ng)%btflux (:,:,itrc), &
523 & update)
524 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
525
526# ifdef TL_IOMS
527 CALL set_2dfld_tile (ng, tile, irpm, idtbot(itrc), &
528 & lbi, ubi, lbj, ubj, &
529 & forces(ng)%btfluxG(:,:,:,itrc), &
530 & forces(ng)%tl_btflux(:,:,itrc), &
531 & update)
532 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
533# endif
534# endif
535 END DO
536# endif
537# endif
538
539# if !defined BULK_FLUXES || defined FORWARD_FLUXES || \
540 defined frc_coupling
541!
542!-----------------------------------------------------------------------
543! Set kinematic surface momentum flux (m2/s2).
544!-----------------------------------------------------------------------
545!
546# ifdef ANA_SMFLUX
547 CALL ana_smflux (ng, tile, irpm)
548# else
549 CALL set_2dfld_tile (ng, tile, irpm, idusms, &
550 & lbi, ubi, lbj, ubj, &
551 & forces(ng)%sustrG, &
552 & forces(ng)%sustr, &
553 & update)
554 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
555!
556 CALL set_2dfld_tile (ng, tile, irpm, idvsms, &
557 & lbi, ubi, lbj, ubj, &
558 & forces(ng)%svstrG, &
559 & forces(ng)%svstr, &
560 & update)
561 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
562
563# ifdef TL_IOMS
564 CALL set_2dfld_tile (ng, tile, irpm, idusms, &
565 & lbi, ubi, lbj, ubj, &
566 & forces(ng)%sustrG, &
567 & forces(ng)%tl_sustr, &
568 & update)
569 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
570!
571 CALL set_2dfld_tile (ng, tile, irpm, idvsms, &
572 & lbi, ubi, lbj, ubj, &
573 & forces(ng)%svstrG, &
574 & forces(ng)%tl_svstr, &
575 & update)
576 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
577# endif
578
579# ifdef CURVGRID
580!
581! If input point wind stress, rotate to curvilinear grid. Notice
582! that rotation is done at RHO-points. It does not matter.
583!
584 IF (.not.linfo(1,idusms,ng).or. &
585 & (iinfo(5,idusms,ng).ne.lm(ng)+1).or. &
586 & (iinfo(6,idusms,ng).ne.mm(ng)+2)) THEN
587 DO j=jstrr,jendr
588 DO i=istrr,iendr
589 cff1=forces(ng)%sustr(i,j)*grid(ng)%CosAngler(i,j)+ &
590 & forces(ng)%svstr(i,j)*grid(ng)%SinAngler(i,j)
591 cff2=forces(ng)%svstr(i,j)*grid(ng)%CosAngler(i,j)- &
592 & forces(ng)%sustr(i,j)*grid(ng)%SinAngler(i,j)
593 forces(ng)%sustr(i,j)=cff1
594 forces(ng)%svstr(i,j)=cff2
595# ifdef TL_IOMS
596 forces(ng)%tl_sustr(i,j)=cff1
597 forces(ng)%tl_svstr(i,j)=cff2
598# endif
599 END DO
600 END DO
601
602 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
603 CALL exchange_u2d_tile (ng, tile, &
604 & lbi, ubi, lbj, ubj, &
605 & forces(ng)%sustr)
606 CALL exchange_v2d_tile (ng, tile, &
607 & lbi, ubi, lbj, ubj, &
608 & forces(ng)%svstr)
609# ifdef TL_IOMS
610 CALL exchange_u2d_tile (ng, tile, &
611 & lbi, ubi, lbj, ubj, &
612 & forces(ng)%tl_sustr)
613 CALL exchange_v2d_tile (ng, tile, &
614 & lbi, ubi, lbj, ubj, &
615 & forces(ng)%tl_svstr)
616# endif
617 END IF
618
619# ifdef DISTRIBUTE
620 CALL mp_exchange2d (ng, tile, irpm, 2, &
621 & lbi, ubi, lbj, ubj, &
622 & nghostpoints, &
623 & ewperiodic(ng), nsperiodic(ng), &
624 & forces(ng)%sustr, &
625 & forces(ng)%svstr)
626# ifdef TL_IOMS
627 CALL mp_exchange2d (ng, tile, irpm, 2, &
628 & lbi, ubi, lbj, ubj, &
629 & nghostpoints, &
630 & ewperiodic(ng), nsperiodic(ng), &
631 & forces(ng)%tl_sustr, &
632 & forces(ng)%tl_svstr)
633# endif
634# endif
635 END IF
636# endif
637# endif
638# endif
639
640# if (defined BULK_FLUXES && !defined FORWARD_FLUXES) || \
641 defined ecosim || defined atm_press
642!
643!-----------------------------------------------------------------------
644! Set surface air pressure (mb).
645!-----------------------------------------------------------------------
646!
647# ifdef ANA_PAIR
648 CALL ana_pair (ng, tile, irpm)
649# else
650 setbc=.true.
651! SetBC=.FALSE.
652 CALL set_2dfld_tile (ng, tile, irpm, idpair, &
653 & lbi, ubi, lbj, ubj, &
654 & forces(ng)%PairG, &
655 & forces(ng)%Pair, &
656 & update, setbc)
657 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
658# endif
659# endif
660
661# ifdef WAVE_DATA
662!
663!-----------------------------------------------------------------------
664! Set surface wind-induced wave amplitude, direction and period.
665!-----------------------------------------------------------------------
666!
667# ifdef ANA_WWAVE
668 CALL ana_wwave (ng, tile, irpm)
669# else
670# ifdef WAVES_DIR
671 CALL set_2dfld_tile (ng, tile, irpm, idwdir, &
672 & lbi, ubi, lbj, ubj, &
673 & forces(ng)%DwaveG, &
674 & forces(ng)%Dwave, &
675 & update)
676 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
677
678# ifdef CURVGRID
679!
680! If input point-data, rotate direction to curvilinear coordinates.
681!
682 IF (.not.linfo(1,idwdir,ng).or. &
683 & (iinfo(5,idwdir,ng).ne.lm(ng)+2).or. &
684 & (iinfo(6,idwdir,ng).ne.mm(ng)+2)) THEN
685 DO j=jstrr,jendr
686 DO i=istrr,iendr
687 forces(ng)%Dwave(i,j)=forces(ng)%Dwave(i,j)- &
688 & grid(ng)%angler(i,j)
689 END DO
690 END DO
691 END IF
692
693 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
694 CALL exchange_r2d_tile (ng, tile, &
695 & lbi, ubi, lbj, ubj, &
696 & forces(ng)%Dwave)
697 END IF
698
699# ifdef DISTRIBUTE
700 CALL mp_exchange2d (ng, tile, irpm, 1, &
701 & lbi, ubi, lbj, ubj, &
702 & nghostpoints, &
703 & ewperiodic(ng), nsperiodic(ng), &
704 & forces(ng)%Dwave)
705# endif
706# endif
707# endif
708
709# ifdef WAVES_HEIGHT
710!
711 CALL set_2dfld_tile (ng, tile, irpm, idwamp, &
712 & lbi, ubi, lbj, ubj, &
713 & forces(ng)%HwaveG, &
714 & forces(ng)%Hwave, &
715 & update)
716 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
717# endif
718
719# ifdef WAVES_LENGTH
720!
721 CALL set_2dfld_tile (ng, tile, irpm, idwlen, &
722 & lbi, ubi, lbj, ubj, &
723 & forces(ng)%LwaveG, &
724 & forces(ng)%Lwave, &
725 & update)
726 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
727# endif
728
729# ifdef WAVES_TOP_PERIOD
730!
731 CALL set_2dfld_tile (ng, tile, irpm, idwptp, &
732 & lbi, ubi, lbj, ubj, &
733 & forces(ng)%Pwave_topG, &
734 & forces(ng)%Pwave_top, &
735 & update)
736 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
737# endif
738
739# ifdef WAVES_BOT_PERIOD
740!
741 CALL set_2dfld_tile (ng, tile, irpm, idwpbt, &
742 & lbi, ubi, lbj, ubj, &
743 & forces(ng)%Pwave_botG, &
744 & forces(ng)%Pwave_bot, &
745 & update)
746 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
747# endif
748
749# if defined WAVES_UB
750!
751 CALL set_2dfld_tile (ng, tile, irpm, idworb, &
752 & lbi, ubi, lbj, ubj, &
753 & forces(ng)%Ub_swanG, &
754 & forces(ng)%Ub_swan, &
755 & update)
756 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
757# endif
758
759# if defined TKE_WAVEDISS
760!
761 CALL set_2dfld_tile (ng, tile, irpm, idwdis, &
762 & lbi, ubi, lbj, ubj, &
763 & forces(ng)%Wave_dissipG, &
764 & forces(ng)%Wave_dissip, &
765 & update)
766 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
767# endif
768
769# if defined SVENDSEN_ROLLER
770!
771 CALL set_2dfld_tile (ng, tile, irpm, idwbrk, &
772 & lbi, ubi, lbj, ubj, &
773 & forces(ng)%Wave_breakG, &
774 & forces(ng)%Wave_break, &
775 & update)
776 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
777# endif
778# endif
779# endif
780
781# if defined ECOSIM && defined SOLVE3D
782!
783!-----------------------------------------------------------------------
784! Compute spectral irradiance and cosine of average zenith angle of
785! downwelling spectral photons.
786!-----------------------------------------------------------------------
787!
788 CALL ana_specir (ng, tile, irpm)
789# endif
790
791# ifdef ANA_SPINNING
792!
793!-----------------------------------------------------------------------
794! Set time-varying rotation force (centripetal accelerations) for
795! polar coordinate grids.
796!-----------------------------------------------------------------------
797!
798 CALL ana_spinning (ng, tile, irpm)
799# endif
800!
801!-----------------------------------------------------------------------
802! Set point Sources/Sinks (river runoff).
803!-----------------------------------------------------------------------
804!
805# ifdef ANA_PSOURCE
806 IF (luvsrc(ng).or.lwsrc(ng).or.any(ltracersrc(:,ng))) THEN
807 CALL ana_psource (ng, tile, irpm)
808 END IF
809# else
810 IF (domain(ng)%SouthWest_Test(tile)) THEN
811 IF (luvsrc(ng).or.lwsrc(ng)) THEN
812 CALL set_ngfld (ng, irpm, idrtra, 1, nsrc(ng), 1, &
813 & 1, nsrc(ng), 1, &
814 & sources(ng) % QbarG, &
815 & sources(ng) % Qbar, &
816 & update)
817 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
818
819# ifdef SOLVE3D
820 DO k=1,n(ng)
821 DO i=1,nsrc(ng)
822 sources(ng)%Qsrc(i,k)=sources(ng)%Qbar(i)* &
823 & sources(ng)%Qshape(i,k)
824 END DO
825 END DO
826# endif
827 END IF
828
829# ifdef SOLVE3D
830 DO itrc=1,nt(ng)
831 IF (ltracersrc(itrc,ng)) THEN
832 CALL set_ngfld (ng, irpm, idrtrc(itrc), 1, nsrc(ng), n(ng), &
833 & 1, nsrc(ng), n(ng), &
834 & sources(ng) % TsrcG(:,:,:,itrc), &
835 & sources(ng) % Tsrc(:,:,itrc), &
836 & update)
837 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
838 END IF
839 END DO
840# endif
841 END IF
842# endif
843!
844!-----------------------------------------------------------------------
845! Set open boundary conditions fields.
846!-----------------------------------------------------------------------
847!
848! Free-surface.
849!
850 IF (lprocessobc(ng)) THEN
851# ifdef ANA_FSOBC
852 CALL ana_fsobc (ng, tile, irpm)
853# else
854 IF (domain(ng)%SouthWest_Test(tile)) THEN
855 IF (tl_lbc(iwest,isfsur,ng)%acquire) THEN
856 CALL set_ngfld (ng, irpm, idzbry(iwest), jlb, jub, 1, &
857 & 0, mm(ng)+1, 1, &
858 & boundary(ng) % zetaG_west, &
859 & boundary(ng) % zeta_west, &
860 & update)
861 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
862
863 CALL set_ngfld (ng, irpm, idzbry(iwest), jlb, jub, 1, &
864 & 0, mm(ng)+1, 1, &
865 & boundary(ng) % zetaG_west, &
866 & boundary(ng) % tl_zeta_west, &
867 & update)
868 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
869 END IF
870
871 IF (tl_lbc(ieast,isfsur,ng)%acquire) THEN
872 CALL set_ngfld (ng, irpm, idzbry(ieast), jlb, jub, 1, &
873 & 0, mm(ng)+1, 1, &
874 & boundary(ng) % zetaG_east, &
875 & boundary(ng) % zeta_east, &
876 & update)
877 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
878
879 CALL set_ngfld (ng, irpm, idzbry(ieast), jlb, jub, 1, &
880 & 0, mm(ng)+1, 1, &
881 & boundary(ng) % zetaG_east, &
882 & boundary(ng) % tl_zeta_east, &
883 & update)
884 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
885 END IF
886
887 IF (tl_lbc(isouth,isfsur,ng)%acquire) THEN
888 CALL set_ngfld (ng, irpm, idzbry(isouth), ilb, iub, 1, &
889 & 0, lm(ng)+1 ,1, &
890 & boundary(ng) % zetaG_south, &
891 & boundary(ng) % zeta_south, &
892 & update)
893 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
894
895 CALL set_ngfld (ng, irpm, idzbry(isouth), ilb, iub, 1, &
896 & 0, lm(ng)+1 ,1, &
897 & boundary(ng) % zetaG_south, &
898 & boundary(ng) % tl_zeta_south, &
899 & update)
900 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
901 END IF
902
903 IF (tl_lbc(inorth,isfsur,ng)%acquire) THEN
904 CALL set_ngfld (ng, irpm, idzbry(inorth), ilb, iub, 1, &
905 & 0, lm(ng)+1, 1, &
906 & boundary(ng) % zetaG_north, &
907 & boundary(ng) % zeta_north, &
908 & update)
909 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
910
911 CALL set_ngfld (ng, irpm, idzbry(inorth), ilb, iub, 1, &
912 & 0, lm(ng)+1, 1, &
913 & boundary(ng) % zetaG_north, &
914 & boundary(ng) % tl_zeta_north, &
915 & update)
916 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
917 END IF
918 END IF
919# endif
920
921# if defined WET_DRY
922!
923! Ensure that water level on boundary cells is above bed elevation.
924!
925 IF (tl_lbc(iwest,isfsur,ng)%acquire) THEN
926 bry_update=.false.
927 IF (domain(ng)%Western_Edge(tile)) THEN
928 DO j=jstrr,jendr
929 cff=dcrit(ng)-grid(ng)%h(0,j)
930 IF (boundary(ng)%zeta_west(j).le.cff) THEN
931 boundary(ng)%zeta_west(j)=cff
932 boundary(ng)%tl_zeta_west(j)=cff
933 END IF
934 END DO
935 bry_update=.true.
936 END IF
937# ifdef DISTRIBUTE
938 CALL mp_boundary (ng, irpm, jstrr, jendr, jlb, jub, 1, 1, &
939 & bry_update, &
940 & boundary(ng)%zeta_west)
941 CALL mp_boundary (ng, irpm, jstrr, jendr, jlb, jub, 1, 1, &
942 & bry_update, &
943 & boundary(ng)%tl_zeta_west)
944# endif
945 END IF
946
947 IF (tl_lbc(ieast,isfsur,ng)%acquire) THEN
948 bry_update=.false.
949 IF (domain(ng)%Eastern_Edge(tile)) THEN
950 DO j=jstrr,jendr
951 cff=dcrit(ng)-grid(ng)%h(lm(ng)+1,j)
952 IF (boundary(ng)%zeta_east(j).le.cff) THEN
953 boundary(ng)%zeta_east(j)=cff
954 boundary(ng)%tl_zeta_east(j)=cff
955 END IF
956 END DO
957 bry_update=.true.
958 END IF
959# ifdef DISTRIBUTE
960 CALL mp_boundary (ng, irpm, jstrr, jendr, jlb, jub, 1, 1, &
961 & bry_update, &
962 & boundary(ng)%zeta_east)
963 CALL mp_boundary (ng, irpm, jstrr, jendr, jlb, jub, 1, 1, &
964 & bry_update, &
965 & boundary(ng)%tl_zeta_east)
966# endif
967 END IF
968
969 IF (tl_lbc(isouth,isfsur,ng)%acquire) THEN
970 bry_update=.false.
971 IF (domain(ng)%Southern_Edge(tile)) THEN
972 DO i=istrr,iendr
973 cff=dcrit(ng)-grid(ng)%h(i,0)
974 IF (boundary(ng)%zeta_south(i).le.cff) THEN
975 boundary(ng)%zeta_south(i)=cff
976 boundary(ng)%tl_zeta_south(i)=cff
977 END IF
978 END DO
979 bry_update=.true.
980 END IF
981# ifdef DISTRIBUTE
982 CALL mp_boundary (ng, irpm, istrr, iendr, ilb, iub, 1, 1, &
983 & bry_update, &
984 & boundary(ng)%zeta_south)
985 CALL mp_boundary (ng, irpm, istrr, iendr, ilb, iub, 1, 1, &
986 & bry_update, &
987 & boundary(ng)%tl_zeta_south)
988# endif
989 END IF
990
991 IF (tl_lbc(inorth,isfsur,ng)%acquire) THEN
992 bry_update=.false.
993 IF (domain(ng)%Northern_Edge(tile)) THEN
994 DO i=istrr,iendr
995 cff=dcrit(ng)-grid(ng)%h(i,mm(ng)+1)
996 IF (boundary(ng)%zeta_north(i).le.cff) THEN
997 boundary(ng)%zeta_north(i)=cff
998 boundary(ng)%tl_zeta_north(i)=cff
999 END IF
1000 END DO
1001 bry_update=.true.
1002 END IF
1003# ifdef DISTRIBUTE
1004 CALL mp_boundary (ng, irpm, istrr, iendr, ilb, iub, 1, 1, &
1005 & bry_update, &
1006 & boundary(ng)%zeta_north)
1007 CALL mp_boundary (ng, irpm, istrr, iendr, ilb, iub, 1, 1, &
1008 & bry_update, &
1009 & boundary(ng)%tl_zeta_north)
1010# endif
1011 END IF
1012# endif
1013 END IF
1014!
1015! 2D momentum.
1016!
1017 IF (lprocessobc(ng)) THEN
1018# ifdef ANA_M2OBC
1019 CALL ana_m2obc (ng, tile, irpm)
1020# else
1021 IF (domain(ng)%SouthWest_Test(tile)) THEN
1022 IF (tl_lbc(iwest,isubar,ng)%acquire) THEN
1023 CALL set_ngfld (ng, irpm, idu2bc(iwest), jlb, jub, 1, &
1024 & 0, mm(ng)+1, 1, &
1025 & boundary(ng) % ubarG_west, &
1026 & boundary(ng) % ubar_west, &
1027 & update)
1028 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1029
1030 CALL set_ngfld (ng, irpm, idu2bc(iwest), jlb, jub, 1, &
1031 & 0, mm(ng)+1, 1, &
1032 & boundary(ng) % ubarG_west, &
1033 & boundary(ng) % tl_ubar_west, &
1034 & update)
1035 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1036 END IF
1037
1038 IF (tl_lbc(iwest,isvbar,ng)%acquire) THEN
1039 CALL set_ngfld (ng, irpm, idv2bc(iwest), jlb, jub, 1, &
1040 & 1, mm(ng)+1, 1, &
1041 & boundary(ng) % vbarG_west, &
1042 & boundary(ng) % vbar_west, &
1043 & update)
1044 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1045
1046 CALL set_ngfld (ng, irpm, idv2bc(iwest), jlb, jub, 1, &
1047 & 1, mm(ng)+1, 1, &
1048 & boundary(ng) % vbarG_west, &
1049 & boundary(ng) % tl_vbar_west, &
1050 & update)
1051 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1052 END IF
1053
1054 IF (tl_lbc(ieast,isubar,ng)%acquire) THEN
1055 CALL set_ngfld (ng, irpm, idu2bc(ieast), jlb, jub, 1, &
1056 & 0, mm(ng)+1, 1, &
1057 & boundary(ng) % ubarG_east, &
1058 & boundary(ng) % ubar_east, &
1059 & update)
1060 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1061
1062 CALL set_ngfld (ng, irpm, idu2bc(ieast), jlb, jub, 1, &
1063 & 0, mm(ng)+1, 1, &
1064 & boundary(ng) % ubarG_east, &
1065 & boundary(ng) % tl_ubar_east, &
1066 & update)
1067 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1068 END IF
1069
1070 IF (tl_lbc(ieast,isvbar,ng)%acquire) THEN
1071 CALL set_ngfld (ng, irpm, idv2bc(ieast), jlb, jub, 1, &
1072 & 1, mm(ng)+1, 1, &
1073 & boundary(ng) % vbarG_east, &
1074 & boundary(ng) % vbar_east, &
1075 & update)
1076 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1077
1078 CALL set_ngfld (ng, irpm, idv2bc(ieast), jlb, jub, 1, &
1079 & 1, mm(ng)+1, 1, &
1080 & boundary(ng) % vbarG_east, &
1081 & boundary(ng) % tl_vbar_east, &
1082 & update)
1083 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1084 END IF
1085
1086 IF (tl_lbc(isouth,isubar,ng)%acquire) THEN
1087 CALL set_ngfld (ng, irpm, idu2bc(isouth), ilb, iub, 1, &
1088 & 1, lm(ng)+1, 1, &
1089 & boundary(ng) % ubarG_south, &
1090 & boundary(ng) % ubar_south, &
1091 & update)
1092 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1093
1094 CALL set_ngfld (ng, irpm, idu2bc(isouth), ilb, iub, 1, &
1095 & 1, lm(ng)+1, 1, &
1096 & boundary(ng) % ubarG_south, &
1097 & boundary(ng) % tl_ubar_south, &
1098 & update)
1099 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1100 END IF
1101
1102 IF (tl_lbc(isouth,isvbar,ng)%acquire) THEN
1103 CALL set_ngfld (ng, irpm, idv2bc(isouth), ilb, iub, 1, &
1104 & 0, lm(ng)+1, 1, &
1105 & boundary(ng) % vbarG_south, &
1106 & boundary(ng) % vbar_south, &
1107 & update)
1108 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1109
1110 CALL set_ngfld (ng, irpm, idv2bc(isouth), ilb, iub, 1, &
1111 & 0, lm(ng)+1, 1, &
1112 & boundary(ng) % vbarG_south, &
1113 & boundary(ng) % tl_vbar_south, &
1114 & update)
1115 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1116 END IF
1117
1118 IF (tl_lbc(inorth,isubar,ng)%acquire) THEN
1119 CALL set_ngfld (ng, irpm, idu2bc(inorth), ilb, iub, 1, &
1120 & 1, lm(ng)+1, 1, &
1121 & boundary(ng) % ubarG_north, &
1122 & boundary(ng) % ubar_north, &
1123 & update)
1124 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1125
1126 CALL set_ngfld (ng, irpm, idu2bc(inorth), ilb, iub, 1, &
1127 & 1, lm(ng)+1, 1, &
1128 & boundary(ng) % ubarG_north, &
1129 & boundary(ng) % tl_ubar_north, &
1130 & update)
1131 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1132 END IF
1133
1134 IF (tl_lbc(inorth,isvbar,ng)%acquire) THEN
1135 CALL set_ngfld (ng, irpm, idv2bc(inorth), ilb, iub, 1, &
1136 & 0, lm(ng)+1, 1, &
1137 & boundary(ng) % vbarG_north, &
1138 & boundary(ng) % vbar_north, &
1139 & update)
1140 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1141
1142 CALL set_ngfld (ng, irpm, idv2bc(inorth), ilb, iub, 1, &
1143 & 0, lm(ng)+1, 1, &
1144 & boundary(ng) % vbarG_north, &
1145 & boundary(ng) % tl_vbar_north, &
1146 & update)
1147 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1148 END IF
1149 END IF
1150# endif
1151 END IF
1152
1153# ifdef SOLVE3D
1154!
1155! 3D momentum.
1156!
1157 IF (lprocessobc(ng)) THEN
1158# ifdef ANA_M3OBC
1159 CALL ana_m3obc (ng, tile, irpm)
1160# else
1161 IF (domain(ng)%SouthWest_Test(tile)) THEN
1162 IF (tl_lbc(iwest,isuvel,ng)%acquire) THEN
1163 CALL set_ngfld (ng, irpm, idu3bc(iwest), jlb, jub, n(ng), &
1164 & 0, mm(ng)+1, n(ng), &
1165 & boundary(ng) % uG_west, &
1166 & boundary(ng) % u_west, &
1167 & update)
1168 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1169
1170 CALL set_ngfld (ng, irpm, idu3bc(iwest), jlb, jub, n(ng), &
1171 & 0, mm(ng)+1, n(ng), &
1172 & boundary(ng) % uG_west, &
1173 & boundary(ng) % tl_u_west, &
1174 & update)
1175 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1176 END IF
1177
1178 IF (tl_lbc(iwest,isvvel,ng)%acquire) THEN
1179 CALL set_ngfld (ng, irpm, idv3bc(iwest), jlb, jub, n(ng), &
1180 & 1, mm(ng)+1, n(ng), &
1181 & boundary(ng) % vG_west, &
1182 & boundary(ng) % v_west, &
1183 & update)
1184 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1185
1186 CALL set_ngfld (ng, irpm, idv3bc(iwest), jlb, jub, n(ng), &
1187 & 1, mm(ng)+1, n(ng), &
1188 & boundary(ng) % vG_west, &
1189 & boundary(ng) % tl_v_west, &
1190 & update)
1191 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1192 END IF
1193
1194 IF (tl_lbc(ieast,isuvel,ng)%acquire) THEN
1195 CALL set_ngfld (ng, irpm, idu3bc(ieast), jlb, jub, n(ng), &
1196 & 0, mm(ng)+1, n(ng), &
1197 & boundary(ng) % uG_east, &
1198 & boundary(ng) % u_east, &
1199 & update)
1200 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1201
1202 CALL set_ngfld (ng, irpm, idu3bc(ieast), jlb, jub, n(ng), &
1203 & 0, mm(ng)+1, n(ng), &
1204 & boundary(ng) % uG_east, &
1205 & boundary(ng) % tl_u_east, &
1206 & update)
1207 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1208 END IF
1209
1210 IF (tl_lbc(ieast,isvvel,ng)%acquire) THEN
1211 CALL set_ngfld (ng, irpm, idv3bc(ieast), jlb, jub, n(ng), &
1212 & 1, mm(ng)+1, n(ng), &
1213 & boundary(ng) % vG_east, &
1214 & boundary(ng) % v_east, &
1215 & update)
1216 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1217
1218 CALL set_ngfld (ng, irpm, idv3bc(ieast), jlb, jub, n(ng), &
1219 & 1, mm(ng)+1, n(ng), &
1220 & boundary(ng) % vG_east, &
1221 & boundary(ng) % tl_v_east, &
1222 & update)
1223 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1224 END IF
1225
1226 IF (tl_lbc(isouth,isuvel,ng)%acquire) THEN
1227 CALL set_ngfld (ng, irpm, idu3bc(isouth), ilb, iub, n(ng), &
1228 & 1, lm(ng)+1, n(ng), &
1229 & boundary(ng) % uG_south, &
1230 & boundary(ng) % u_south, &
1231 & update)
1232 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1233
1234 CALL set_ngfld (ng, irpm, idu3bc(isouth), ilb, iub, n(ng), &
1235 & 1, lm(ng)+1, n(ng), &
1236 & boundary(ng) % uG_south, &
1237 & boundary(ng) % tl_u_south, &
1238 & update)
1239 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1240 END IF
1241
1242 IF (tl_lbc(isouth,isvvel,ng)%acquire) THEN
1243 CALL set_ngfld (ng, irpm, idv3bc(isouth), ilb, iub, n(ng), &
1244 & 0, lm(ng)+1, n(ng), &
1245 & boundary(ng) % vG_south, &
1246 & boundary(ng) % v_south, &
1247 & update)
1248 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1249
1250 CALL set_ngfld (ng, irpm, idv3bc(isouth), ilb, iub, n(ng), &
1251 & 0, lm(ng)+1, n(ng), &
1252 & boundary(ng) % vG_south, &
1253 & boundary(ng) % tl_v_south, &
1254 & update)
1255 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1256 END IF
1257
1258 IF (tl_lbc(inorth,isuvel,ng)%acquire) THEN
1259 CALL set_ngfld (ng, irpm, idu3bc(inorth), ilb, iub, n(ng), &
1260 & 1, lm(ng)+1, n(ng), &
1261 & boundary(ng) % uG_north, &
1262 & boundary(ng) % u_north, &
1263 & update)
1264 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1265
1266 CALL set_ngfld (ng, irpm, idu3bc(inorth), ilb, iub, n(ng), &
1267 & 1, lm(ng)+1, n(ng), &
1268 & boundary(ng) % uG_north, &
1269 & boundary(ng) % tl_u_north, &
1270 & update)
1271 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1272 END IF
1273
1274 IF (tl_lbc(inorth,isvvel,ng)%acquire) THEN
1275 CALL set_ngfld (ng, irpm, idv3bc(inorth), ilb, iub, n(ng), &
1276 & 0, lm(ng)+1, n(ng), &
1277 & boundary(ng) % vG_north, &
1278 & boundary(ng) % v_north, &
1279 & update)
1280 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1281
1282 CALL set_ngfld (ng, irpm, idv3bc(inorth), ilb, iub, n(ng), &
1283 & 0, lm(ng)+1, n(ng), &
1284 & boundary(ng) % vG_north, &
1285 & boundary(ng) % tl_v_north, &
1286 & update)
1287 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1288 END IF
1289 END IF
1290# endif
1291 END IF
1292!
1293! Tracers.
1294!
1295 IF (lprocessobc(ng)) THEN
1296# ifdef ANA_TOBC
1297 CALL ana_tobc (ng, tile, irpm)
1298# else
1299 IF (domain(ng)%SouthWest_Test(tile)) THEN
1300 DO itrc=1,nt(ng)
1301 IF (tl_lbc(iwest,istvar(itrc),ng)%acquire) THEN
1302 CALL set_ngfld (ng, irpm, idtbry(iwest,itrc), &
1303 & jlb, jub, n(ng), 0, mm(ng)+1, n(ng), &
1304 & boundary(ng) % tG_west(:,:,:,itrc), &
1305 & boundary(ng) % t_west(:,:,itrc), &
1306 & update)
1308 & __line__, myfile)) RETURN
1309
1310 CALL set_ngfld (ng, irpm, idtbry(iwest,itrc), &
1311 & jlb, jub, n(ng), 0, mm(ng)+1, n(ng), &
1312 & boundary(ng) % tG_west(:,:,:,itrc), &
1313 & boundary(ng) % tl_t_west(:,:,itrc), &
1314 & update)
1316 & __line__, myfile)) RETURN
1317 END IF
1318
1319 IF (tl_lbc(ieast,istvar(itrc),ng)%acquire) THEN
1320 CALL set_ngfld (ng, irpm, idtbry(ieast,itrc), &
1321 & jlb, jub, n(ng), 0, mm(ng)+1, n(ng), &
1322 & boundary(ng) % tG_east(:,:,:,itrc), &
1323 & boundary(ng) % t_east(:,:,itrc), &
1324 & update)
1326 & __line__, myfile)) RETURN
1327
1328 CALL set_ngfld (ng, irpm, idtbry(ieast,itrc), &
1329 & jlb, jub, n(ng), 0, mm(ng)+1, n(ng), &
1330 & boundary(ng) % tG_east(:,:,:,itrc), &
1331 & boundary(ng) % tl_t_east(:,:,itrc), &
1332 & update)
1334 & __line__, myfile)) RETURN
1335 END IF
1336
1337 IF (tl_lbc(isouth,istvar(itrc),ng)%acquire) THEN
1338 CALL set_ngfld (ng, irpm, idtbry(isouth,itrc), &
1339 & ilb, iub, n(ng), 0, lm(ng)+1, n(ng), &
1340 & boundary(ng) % tG_south(:,:,:,itrc), &
1341 & boundary(ng) % t_south(:,:,itrc), &
1342 & update)
1344 & __line__, myfile)) RETURN
1345
1346 CALL set_ngfld (ng, irpm, idtbry(isouth,itrc), &
1347 & ilb, iub, n(ng), 0, lm(ng)+1, n(ng), &
1348 & boundary(ng) % tG_south(:,:,:,itrc), &
1349 & boundary(ng) % tl_t_south(:,:,itrc), &
1350 & update)
1352 & __line__, myfile)) RETURN
1353 END IF
1354
1355 IF (tl_lbc(inorth,istvar(itrc),ng)%acquire) THEN
1356 CALL set_ngfld (ng, irpm, idtbry(inorth,itrc), &
1357 & ilb, iub, n(ng), 0, lm(ng)+1, n(ng), &
1358 & boundary(ng) % tG_north(:,:,:,itrc), &
1359 & boundary(ng) % t_north(:,:,itrc), &
1360 & update)
1362 & __line__, myfile)) RETURN
1363
1364 CALL set_ngfld (ng, irpm, idtbry(inorth,itrc), &
1365 & ilb, iub, n(ng), 0, lm(ng)+1, n(ng), &
1366 & boundary(ng) % tG_north(:,:,:,itrc), &
1367 & boundary(ng) % tl_t_north(:,:,itrc), &
1368 & update)
1370 & __line__, myfile)) RETURN
1371 END IF
1372 END DO
1373 END IF
1374# endif
1375 END IF
1376# endif
1377!
1378!-----------------------------------------------------------------------
1379! Set sea surface height climatology (m).
1380!-----------------------------------------------------------------------
1381!
1382 IF (lsshclm(ng)) THEN
1383# ifdef ANA_SSH
1384 CALL ana_ssh (ng, tile, irpm)
1385# else
1386 CALL set_2dfld_tile (ng, tile, irpm, idsshc, &
1387 & lbi, ubi, lbj, ubj, &
1388 & clima(ng)%sshG, &
1389 & clima(ng)%ssh, &
1390 & update)
1391 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1392# endif
1393 END IF
1394!
1395!-----------------------------------------------------------------------
1396! Set 2D momentum climatology (m/s).
1397!-----------------------------------------------------------------------
1398!
1399 IF (lm2clm(ng)) THEN
1400# ifdef ANA_M2CLIMA
1401 CALL ana_m2clima (ng, tile, irpm)
1402# else
1403 CALL set_2dfld_tile (ng, tile, irpm, idubcl, &
1404 & lbi, ubi, lbj, ubj, &
1405 & clima(ng)%ubarclmG, &
1406 & clima(ng)%ubarclm, &
1407 & update)
1408 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1409!
1410 CALL set_2dfld_tile (ng, tile, irpm, idvbcl, &
1411 & lbi, ubi, lbj, ubj, &
1412 & clima(ng)%vbarclmG, &
1413 & clima(ng)%vbarclm, &
1414 & update)
1415 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1416# endif
1417 END IF
1418
1419# ifdef SOLVE3D
1420!
1421!-----------------------------------------------------------------------
1422! Set 3D momentum climatology (m/s).
1423!-----------------------------------------------------------------------
1424!
1425 IF (lm3clm(ng)) THEN
1426# ifdef ANA_M3CLIMA
1427 CALL ana_m3clima (ng, tile, irpm)
1428# else
1429 CALL set_3dfld_tile (ng, tile, irpm, iduclm, &
1430 & lbi, ubi, lbj, ubj, 1, n(ng), &
1431 & clima(ng)%uclmG, &
1432 & clima(ng)%uclm, &
1433 & update)
1434 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1435!
1436 CALL set_3dfld_tile (ng, tile, irpm, idvclm, &
1437 & lbi, ubi, lbj, ubj, 1, n(ng), &
1438 & clima(ng)%vclmG, &
1439 & clima(ng)%vclm, &
1440 & update)
1441 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1442# endif
1443 END IF
1444!
1445!-----------------------------------------------------------------------
1446! Set tracer climatology.
1447!-----------------------------------------------------------------------
1448!
1449# ifdef ANA_TCLIMA
1450 IF (any(ltracerclm(:,ng))) THEN
1451 CALL ana_tclima (ng, tile, irpm)
1452 END IF
1453# else
1454 ic=0
1455 DO itrc=1,nt(ng)
1456 IF (ltracerclm(itrc,ng)) THEN
1457 ic=ic+1
1458 CALL set_3dfld_tile (ng, tile, irpm, idtclm(itrc), &
1459 & lbi, ubi, lbj, ubj, 1, n(ng), &
1460 & clima(ng)%tclmG(:,:,:,:,ic), &
1461 & clima(ng)%tclm (:,:,:,ic), &
1462 & update)
1463 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1464 END IF
1465 END DO
1466# endif
1467# endif
1468
1469# ifdef FORWARD_READ
1470!
1471!-----------------------------------------------------------------------
1472! Set forward solution needed by Tangent/Adjoint models.
1473!-----------------------------------------------------------------------
1474!
1475! Set forward free-surface.
1476!
1477 DO k=1,3
1478 CALL set_2dfld_tile (ng, tile, irpm, idfsur, &
1479 & lbi, ubi, lbj, ubj, &
1480 & ocean(ng)%zetaG, &
1481 & ocean(ng)%zeta(:,:,k), &
1482 & update)
1483 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1484 END DO
1485
1486# ifdef SOLVE3D
1487 DO j=jstrr,jendr
1488 DO i=istrr,iendr
1489 coupling(ng)%Zt_avg1(i,j)=ocean(ng)%zeta(i,j,1)
1490 END DO
1491 END DO
1492# endif
1493!
1494! Set forward 2D momentum.
1495!
1496 DO k=1,3
1497 CALL set_2dfld_tile (ng, tile, irpm, idubar, &
1498 & lbi, ubi, lbj, ubj, &
1499 & ocean(ng)%ubarG, &
1500 & ocean(ng)%ubar(:,:,k), &
1501 & update)
1502 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1503
1504 CALL set_2dfld_tile (ng, tile, irpm, idvbar, &
1505 & lbi, ubi, lbj, ubj, &
1506 & ocean(ng)%vbarG, &
1507 & ocean(ng)%vbar(:,:,k), &
1508 & update)
1509 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1510 END DO
1511
1512# ifdef FORWARD_RHS
1513!
1514! Set forward variables associated with 2D right-hand-side terms.
1515!
1516 DO k=1,2
1517 CALL set_2dfld_tile (ng, tile, irpm, idrzet, &
1518 & lbi, ubi, lbj, ubj, &
1519 & ocean(ng)%rzetaG, &
1520 & ocean(ng)%rzeta(:,:,k), &
1521 & update)
1522 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1523!
1524 CALL set_2dfld_tile (ng, tile, irpm, idru2d, &
1525 & lbi, ubi, lbj, ubj, &
1526 & ocean(ng)%rubarG, &
1527 & ocean(ng)%rubar(:,:,k), &
1528 & update)
1529 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1530!
1531 CALL set_2dfld_tile (ng, tile, irpm, idrv2d, &
1532 & lbi, ubi, lbj, ubj, &
1533 & ocean(ng)%rvbarG, &
1534 & ocean(ng)%rvbar(:,:,k), &
1535 & update)
1536 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1537 END DO
1538# endif
1539
1540# ifdef SOLVE3D
1541!
1542! Set forward time-averaged barotropic fluxes.
1543!
1544 CALL set_2dfld_tile (ng, tile, irpm, idufx1, &
1545 & lbi, ubi, lbj, ubj, &
1546 & coupling(ng)%DU_avg1G, &
1547 & coupling(ng)%DU_avg1, &
1548 & update)
1549 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1550!
1551 CALL set_2dfld_tile (ng, tile, irpm, idufx2, &
1552 & lbi, ubi, lbj, ubj, &
1553 & coupling(ng)%DU_avg2G, &
1554 & coupling(ng)%DU_avg2, &
1555 & update)
1556 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1557!
1558 CALL set_2dfld_tile (ng, tile, irpm, idvfx1, &
1559 & lbi, ubi, lbj, ubj, &
1560 & coupling(ng)%DV_avg1G, &
1561 & coupling(ng)%DV_avg1, &
1562 & update)
1563 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1564!
1565 CALL set_2dfld_tile (ng, tile, irpm, idvfx2, &
1566 & lbi, ubi, lbj, ubj, &
1567 & coupling(ng)%DV_avg2G, &
1568 & coupling(ng)%DV_avg2, &
1569 & update)
1570 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1571!
1572! Set 3D momentum.
1573!
1574 DO k=1,2
1575 CALL set_3dfld_tile (ng, tile, irpm, iduvel, &
1576 & lbi, ubi, lbj, ubj, 1, n(ng), &
1577 & ocean(ng)%uG, &
1578 & ocean(ng)%u(:,:,:,k), &
1579 & update)
1580 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1581!
1582 CALL set_3dfld_tile (ng, tile, irpm, idvvel, &
1583 & lbi, ubi, lbj, ubj, 1, n(ng), &
1584 & ocean(ng)%vG, &
1585 & ocean(ng)%v(:,:,:,k), &
1586 & update)
1587 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1588 END DO
1589
1590# ifdef FORWARD_RHS
1591!
1592! Set variables associated with 3D momentum right-hand-side terms.
1593!
1594 DO k=1,2
1595 CALL set_3dfld_tile (ng, tile, irpm, idru3d, &
1596 & lbi, ubi, lbj, ubj, 1, n(ng), &
1597 & ocean(ng)%ruG, &
1598 & ocean(ng)%ru(:,:,:,k), &
1599 & update)
1600 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1601!
1602 CALL set_3dfld_tile (ng, tile, irpm, idrv3d, &
1603 & lbi, ubi, lbj, ubj, 1, n(ng), &
1604 & ocean(ng)%rvG, &
1605 & ocean(ng)%rv(:,:,:,k), &
1606 & update)
1607 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1608 END DO
1609# endif
1610!
1611! Set 3D tracers.
1612!
1613 DO itrc=1,nt(ng)
1614 DO k=1,3
1615 CALL set_3dfld_tile (ng, tile, irpm, idtvar(itrc), &
1616 & lbi, ubi, lbj, ubj, 1, n(ng), &
1617 & ocean(ng)%tG(:,:,:,:,itrc), &
1618 & ocean(ng)%t(:,:,:,k,itrc), &
1619 & update)
1620 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1621 END DO
1622 END DO
1623
1624# ifdef FORWARD_MIXING
1625!
1626! Set forward vertical mixing variables.
1627!
1628 DO itrc=1,nat
1629 CALL set_3dfld_tile (ng, tile, irpm, iddiff(itrc), &
1630 & lbi, ubi, lbj, ubj, 0, n(ng), &
1631 & mixing(ng)%AktG(:,:,:,:,itrc), &
1632 & mixing(ng)%Akt(:,:,:,itrc), &
1633 & update)
1634 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1635 END DO
1636!
1637 CALL set_3dfld_tile (ng, tile, irpm, idvvis, &
1638 & lbi, ubi, lbj, ubj, 0, n(ng), &
1639 & mixing(ng)%AkvG, &
1640 & mixing(ng)%Akv, &
1641 & update)
1642 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1643# endif
1644
1645# if defined MY25_MIXING_NOT_YET || defined GLS_MIXING_NOT_YET
1646!
1647! Set forward turbulent kinetic energy.
1648!
1649 DO k=1,3
1650 CALL set_3dfld_tile (ng, tile, irpm, idmtke, &
1651 & lbi, ubi, lbj, ubj, 0, n(ng), &
1652 & mixing(ng)%tkeG, &
1653 & mixing(ng)%tke(:,:,:,k), &
1654 & update)
1655 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1656 END DO
1657!
1658! Set forward turbulent kinetic energy times length scale.
1659!
1660 DO k=1,3
1661 CALL set_3dfld_tile (ng, tile, irpm, idmtls, &
1662 & lbi, ubi, lbj, ubj, 0, n(ng), &
1663 & mixing(ng)%glsG, &
1664 & mixing(ng)%gls(:,:,:,k), &
1665 & update)
1666 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1667 END DO
1668!
1669! Set forward vertical mixing length scale.
1670!
1671 CALL set_3dfld_tile (ng, tile, irpm, idvmls, &
1672 & lbi, ubi, lbj, ubj, 0, n(ng), &
1673 & mixing(ng)%LscaleG, &
1674 & mixing(ng)%Lscale, &
1675 & update)
1676 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1677!
1678! Set forward vertical mixing coefficient for turbulent kinetic energy.
1679!
1680 CALL set_3dfld_tile (ng, tile, irpm, idvmkk, &
1681 & lbi, ubi, lbj, ubj, 0, n(ng), &
1682 & mixing(ng)%AkkG, &
1683 & mixing(ng)%Akk, &
1684 & update)
1685 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1686
1687# ifdef GLS_MIXING_NOT_YET
1688!
1689! Set forward vertical mixing coefficient for turbulent length scale.
1690!
1691 CALL set_3dfld_tile (ng, tile, irpm, idvmkp, &
1692 & lbi, ubi, lbj, ubj, 0, n(ng), &
1693 & mixing(ng)%AkpG, &
1694 & mixing(ng)%Akp, &
1695 & update)
1696 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1697# endif
1698# endif
1699
1700# ifdef LMD_MIXING_NOT_YET
1701!
1702! Set forward depth of surface oceanic boundary layer.
1703!
1704 CALL set_2dfld_tile (ng, tile, irpm, idhsbl, &
1705 & lbi, ubi, lbj, ubj, &
1706 & mixing(ng)%hsblG, &
1707 & mixing(ng)%hsbl, &
1708 & update)
1709 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1710# endif
1711
1712# ifdef LMD_BKPP_NOT_YET
1713!
1714! Set forward depth of bottom oceanic boundary layer.
1715!
1716 CALL set_2dfld_tile (ng, tile, irpm, idhbbl, &
1717 & lbi, ubi, lbj, ubj, &
1718 & mixing(ng)%hbblG, &
1719 & mixing(ng)%hbbl, &
1720 & update)
1721 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1722# endif
1723
1724# ifdef LMD_NONLOCAL_NOT_YET
1725!
1726! Set forward boundary layer nonlocal transport.
1727!
1728 DO itrc=1,nat
1729 CALL set_3dfld_tile (ng, tile, irpm, idghat(itrc), &
1730 & lbi, ubi, lbj, ubj, 0, n(ng), &
1731 & mixing(ng)%ghatsG(:,:,:,:,itrc), &
1732 & mixing(ng)%ghat(:,:,:,itrc), &
1733 & update)
1734 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1735 END DO
1736# endif
1737# endif
1738
1739# if defined R4DVAR || defined R4DVAR_ANA_SENSITIVITY || \
1740 defined tl_r4dvar
1741!
1742!-----------------------------------------------------------------------
1743! Set weak contraint frequent impulse forcing.
1744!-----------------------------------------------------------------------
1745!
1746 IF (frequentimpulse(ng)) THEN
1747!
1748! Set free-surface forcing.
1749!
1750 CALL set_2dfld_tile (ng, tile, irpm, idztlf, &
1751 & lbi, ubi, lbj, ubj, &
1752 & ocean(ng)%f_zetaG, &
1753 & ocean(ng)%f_zeta, &
1754 & update)
1755 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1756
1757# ifndef SOLVE3D
1758!
1759! Set 2D momentum forcing.
1760!
1761 CALL set_2dfld_tile (ng, tile, irpm, idubtf, &
1762 & lbi, ubi, lbj, ubj, &
1763 & ocean(ng)%f_ubarG, &
1764 & ocean(ng)%f_ubar, &
1765 & update)
1766 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1767
1768 CALL set_2dfld_tile (ng, tile, irpm, idvbtf, &
1769 & lbi, ubi, lbj, ubj, &
1770 & ocean(ng)%f_vbarG, &
1771 & ocean(ng)%f_vbar, &
1772 & update)
1773 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1774# else
1775!
1776! Set 3D momentum forcing.
1777!
1778 CALL set_3dfld_tile (ng, tile, irpm, idutlf, &
1779 & lbi, ubi, lbj, ubj, 1, n(ng), &
1780 & ocean(ng)%f_uG, &
1781 & ocean(ng)%f_u, &
1782 & update)
1783 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1784!
1785 CALL set_3dfld_tile (ng, tile, irpm, idvtlf, &
1786 & lbi, ubi, lbj, ubj, 1, n(ng), &
1787 & ocean(ng)%f_vG, &
1788 & ocean(ng)%f_v, &
1789 & update)
1790 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1791!
1792! Set 3D tracers forcing.
1793!
1794 DO itrc=1,nt(ng)
1795 CALL set_3dfld_tile (ng, tile, irpm, idttlf(itrc), &
1796 & lbi, ubi, lbj, ubj, 1, n(ng), &
1797 & ocean(ng)%f_tG(:,:,:,:,itrc), &
1798 & ocean(ng)%f_t(:,:,:,itrc), &
1799 & update)
1800 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1801 END DO
1802# endif
1803 END IF
1804# endif
1805# endif
1806!
1807 RETURN
subroutine ana_spinning(ng, tile, model)
Definition ana_spinning.h:3
subroutine ana_fsobc(ng, tile, model)
Definition ana_fsobc.h:3
subroutine ana_wwave(ng, tile, model)
Definition ana_wwave.h:3
subroutine ana_sst(ng, tile, model)
Definition ana_sst.h:3
subroutine ana_m3obc(ng, tile, model)
Definition ana_m3obc.h:3
subroutine ana_dqdsst(ng, tile, model)
Definition ana_dqdsst.h:3
subroutine ana_tclima(ng, tile, model)
Definition ana_tclima.h:3
subroutine ana_btflux(ng, tile, model, itrc)
Definition ana_btflux.h:3
subroutine ana_ssh(ng, tile, model)
Definition ana_ssh.h:3
subroutine ana_m2obc(ng, tile, model)
Definition ana_m2obc.h:3
subroutine ana_tobc(ng, tile, model)
Definition ana_tobc.h:3
subroutine ana_sss(ng, tile, model)
Definition ana_sss.h:3
subroutine ana_psource(ng, tile, model)
Definition ana_psource.h:3
subroutine ana_winds(ng, tile, model)
Definition ana_winds.h:3
subroutine ana_srflux(ng, tile, model)
Definition ana_srflux.h:3
subroutine ana_smflux(ng, tile, model)
Definition ana_smflux.h:3
subroutine ana_pair(ng, tile, model)
Definition ana_pair.h:3
subroutine ana_m2clima(ng, tile, model)
Definition ana_m2clima.h:3
subroutine ana_tair(ng, tile, model)
Definition ana_tair.h:3
subroutine ana_m3clima(ng, tile, model)
Definition ana_m3clima.h:3
subroutine ana_stflux(ng, tile, model, itrc)
Definition ana_stflux.h:3
subroutine ana_rain(ng, tile, model)
Definition ana_rain.h:3
subroutine ana_specir(ng, tile, model)
Definition ana_specir.h:3
subroutine ana_humid(ng, tile, model)
Definition ana_humid.h:3
subroutine ana_cloud(ng, tile, model)
Definition ana_cloud.h:3
subroutine mp_boundary(ng, model, imin, imax, lbi, ubi, lbk, ubk, update, a)
subroutine exchange_r2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine exchange_u2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine exchange_v2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
type(t_boundary), dimension(:), allocatable boundary
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
type(t_coupling), dimension(:), allocatable coupling
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_mixing), dimension(:), allocatable mixing
Definition mod_mixing.F:399
integer iddqdt
integer idvair
integer idvmls
integer idubtf
integer, dimension(2) iddiff
integer, dimension(4) idu3bc
integer idrtra
integer idwdis
integer idrv3d
integer idubar
integer, dimension(4) idzbry
integer idvvel
integer idhsbl
integer, dimension(4) idu2bc
integer idvsms
logical, dimension(:,:,:), allocatable linfo
integer idcfra
integer idwlen
integer isvvel
integer isvbar
integer idvclm
integer, dimension(:), allocatable idrtrc
integer idpair
integer, dimension(:), allocatable idtbot
integer idvfx2
integer, dimension(:), allocatable idtsur
integer idru2d
integer idvmkp
integer idempf
integer, dimension(:), allocatable idtclm
integer idfsur
integer idwbrk
integer, dimension(:), allocatable idtvar
integer, dimension(:), allocatable istvar
integer idhbbl
integer idvfx1
integer idldwn
integer isuvel
integer idufx2
integer idsssc
integer isfsur
integer iduclm
integer idztlf
integer idsfwf
integer idvbtf
integer, dimension(4) idv3bc
integer iduair
integer, dimension(:), allocatable idttlf
integer idmtke
integer iduvel
integer, dimension(2) idghat
integer idqair
integer idutlf
integer isubar
integer, dimension(:,:,:), allocatable iinfo
integer idubcl
integer, dimension(4) idv2bc
integer idlrad
integer idvtlf
integer idru3d
integer idusms
integer idvmkk
integer idwamp
integer idwdir
integer idvvis
integer idwptp
integer idrzet
integer idufx1
integer idrain
integer idvbcl
integer idsrad
integer idmtls
integer idsshc
integer idwpbt
integer idworb
integer idsstc
integer, dimension(:,:), allocatable idtbry
integer idtair
integer idrv2d
integer idvbar
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer nat
Definition mod_param.F:499
integer, dimension(:), allocatable n
Definition mod_param.F:479
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
type(t_lbc), dimension(:,:,:), allocatable tl_lbc
Definition mod_param.F:379
integer nghostpoints
Definition mod_param.F:710
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, dimension(:), allocatable lm
Definition mod_param.F:455
integer, dimension(:), allocatable nt
Definition mod_param.F:489
integer, dimension(:), allocatable mm
Definition mod_param.F:456
logical, dimension(:), allocatable luvsrc
logical, dimension(:,:), allocatable ltracersrc
real(r8), dimension(:), allocatable dcrit
logical, dimension(:), allocatable lprocessobc
logical, dimension(:), allocatable ewperiodic
integer, parameter iwest
logical, dimension(:), allocatable nsperiodic
logical, dimension(:), allocatable lm3clm
logical, dimension(:), allocatable lsshclm
logical, dimension(:), allocatable frequentimpulse
logical, dimension(:), allocatable lwsrc
integer exit_flag
integer isalt
integer itemp
integer, parameter isouth
logical, dimension(:), allocatable lm2clm
integer, parameter ieast
logical, dimension(:,:), allocatable ltracerclm
integer, parameter inorth
integer noerror
type(t_sources), dimension(:), allocatable sources
Definition mod_sources.F:90
integer, dimension(:), allocatable nsrc
Definition mod_sources.F:97
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)
subroutine mp_exchange3d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, nghost, ew_periodic, ns_periodic, a, b, c, d)
subroutine set_2dfld_tile(ng, tile, model, ifield, lbi, ubi, lbj, ubj, finp, fout, update, setbc)
Definition set_2dfld.F:25
subroutine set_3dfld_tile(ng, tile, model, ifield, lbi, ubi, lbj, ubj, lbk, ubk, finp, fout, update, setbc)
Definition set_3dfld.F:26
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52
subroutine set_ngfld(ng, model, ifield, lbi, ubi, ubj, istr, iend, jrec, finp, fout, update)
Definition set_ngfld.F:5

References analytical_mod::ana_btflux(), analytical_mod::ana_cloud(), analytical_mod::ana_dqdsst(), analytical_mod::ana_fsobc(), analytical_mod::ana_humid(), analytical_mod::ana_m2clima(), analytical_mod::ana_m2obc(), analytical_mod::ana_m3clima(), analytical_mod::ana_m3obc(), analytical_mod::ana_pair(), analytical_mod::ana_psource(), analytical_mod::ana_rain(), analytical_mod::ana_smflux(), analytical_mod::ana_specir(), analytical_mod::ana_spinning(), analytical_mod::ana_srflux(), analytical_mod::ana_ssh(), analytical_mod::ana_sss(), analytical_mod::ana_sst(), analytical_mod::ana_stflux(), analytical_mod::ana_tair(), analytical_mod::ana_tclima(), analytical_mod::ana_tobc(), analytical_mod::ana_winds(), analytical_mod::ana_wwave(), mod_boundary::boundary, mod_param::bounds, mod_clima::clima, mod_coupling::coupling, mod_scalars::dcrit, mod_param::domain, mod_scalars::ewperiodic, exchange_2d_mod::exchange_r2d_tile(), exchange_2d_mod::exchange_u2d_tile(), exchange_2d_mod::exchange_v2d_tile(), mod_scalars::exit_flag, mod_forces::forces, strings_mod::founderror(), mod_scalars::frequentimpulse, mod_grid::grid, mod_ncparam::idcfra, mod_ncparam::iddiff, mod_ncparam::iddqdt, mod_ncparam::idempf, mod_ncparam::idfsur, mod_ncparam::idghat, mod_ncparam::idhbbl, mod_ncparam::idhsbl, mod_ncparam::idldwn, mod_ncparam::idlrad, mod_ncparam::idmtke, mod_ncparam::idmtls, mod_ncparam::idpair, mod_ncparam::idqair, mod_ncparam::idrain, mod_ncparam::idrtra, mod_ncparam::idrtrc, mod_ncparam::idru2d, mod_ncparam::idru3d, mod_ncparam::idrv2d, mod_ncparam::idrv3d, mod_ncparam::idrzet, mod_ncparam::idsfwf, mod_ncparam::idsrad, mod_ncparam::idsshc, mod_ncparam::idsssc, mod_ncparam::idsstc, mod_ncparam::idtair, mod_ncparam::idtbot, mod_ncparam::idtbry, mod_ncparam::idtclm, mod_ncparam::idtsur, mod_ncparam::idttlf, mod_ncparam::idtvar, mod_ncparam::idu2bc, mod_ncparam::idu3bc, mod_ncparam::iduair, mod_ncparam::idubar, mod_ncparam::idubcl, mod_ncparam::idubtf, mod_ncparam::iduclm, mod_ncparam::idufx1, mod_ncparam::idufx2, mod_ncparam::idusms, mod_ncparam::idutlf, mod_ncparam::iduvel, mod_ncparam::idv2bc, mod_ncparam::idv3bc, mod_ncparam::idvair, mod_ncparam::idvbar, mod_ncparam::idvbcl, mod_ncparam::idvbtf, mod_ncparam::idvclm, mod_ncparam::idvfx1, mod_ncparam::idvfx2, mod_ncparam::idvmkk, mod_ncparam::idvmkp, mod_ncparam::idvmls, mod_ncparam::idvsms, mod_ncparam::idvtlf, mod_ncparam::idvvel, mod_ncparam::idvvis, mod_ncparam::idwamp, mod_ncparam::idwbrk, mod_ncparam::idwdir, mod_ncparam::idwdis, mod_ncparam::idwlen, mod_ncparam::idworb, mod_ncparam::idwpbt, mod_ncparam::idwptp, mod_ncparam::idzbry, mod_ncparam::idztlf, mod_scalars::ieast, mod_ncparam::iinfo, mod_scalars::inorth, mod_param::irpm, mod_scalars::isalt, mod_ncparam::isfsur, mod_scalars::isouth, mod_ncparam::istvar, mod_ncparam::isubar, mod_ncparam::isuvel, mod_ncparam::isvbar, mod_ncparam::isvvel, mod_scalars::itemp, mod_scalars::iwest, mod_ncparam::linfo, mod_param::lm, mod_scalars::lm2clm, mod_scalars::lm3clm, mod_scalars::lprocessobc, mod_scalars::lsshclm, mod_scalars::ltracerclm, mod_scalars::ltracersrc, mod_scalars::luvsrc, mod_scalars::lwsrc, mod_mixing::mixing, mod_param::mm, distribute_mod::mp_boundary(), mp_exchange_mod::mp_exchange2d(), mp_exchange_mod::mp_exchange3d(), mod_param::n, mod_param::nat, mod_param::nghostpoints, mod_scalars::noerror, mod_scalars::nsperiodic, mod_sources::nsrc, mod_param::nt, mod_ocean::ocean, set_2dfld_mod::set_2dfld_tile(), set_3dfld_mod::set_3dfld_tile(), set_ngfld(), mod_sources::sources, and mod_param::tl_lbc.

Referenced by rp_set_data().

Here is the call graph for this function:
Here is the caller graph for this function: