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