New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
trcsms_fabm.F90 in branches/UKMO/AMM15_v3_6_STABLE_package_collate_utils366_fabmv1/NEMOGCM/NEMO/TOP_SRC/FABM – NEMO

source: branches/UKMO/AMM15_v3_6_STABLE_package_collate_utils366_fabmv1/NEMOGCM/NEMO/TOP_SRC/FABM/trcsms_fabm.F90 @ 13241

Last change on this file since 13241 was 13241, checked in by dford, 4 years ago

Update NEMO-FABM coupler for compatability with FABM v1.0.

File size: 24.7 KB
Line 
1MODULE trcsms_fabm
2   !!======================================================================
3   !!                         ***  MODULE trcsms_fabm  ***
4   !! TOP :   Main module of the FABM tracers
5   !!======================================================================
6   !! History :   1.0  !  2015-04  (PML) Original code
7   !! History :   1.1  !  2020-06  (PML) Update to FABM 1.0, improved performance
8   !!----------------------------------------------------------------------
9#if defined key_fabm
10   !!----------------------------------------------------------------------
11   !!   'key_fabm'                                               FABM tracers
12   !!----------------------------------------------------------------------
13   !! trc_sms_fabm       : FABM model main routine
14   !! trc_sms_fabm_alloc : allocate arrays specific to FABM sms
15   !!----------------------------------------------------------------------
16   USE par_trc         ! TOP parameters
17   USE oce_trc         ! Ocean variables
18   USE trc             ! TOP variables
19   USE trcbc
20   USE trd_oce
21   USE trdtrc
22   USE traqsr,only: ln_qsr_spec
23#if defined key_trdtrc && defined key_iomput
24   USE trdtrc_oce
25#endif
26
27   USE oce, only: tsn  ! Needed?
28   USE sbc_oce, only: lk_oasis,fr_i
29   USE dom_oce
30   USE zdf_oce
31   USE iom
32   USE xios
33   USE cpl_oasis3
34   USE st2D_fabm
35   USE inputs_fabm
36   USE vertical_movement_fabm
37   USE zdfmxl
38   USE asmbgc, ONLY: mld_choice_bgc
39   USE lbclnk
40   USE fabm_types, ONLY: standard_variables
41
42   !USE fldread         !  time interpolation
43
44   IMPLICIT NONE
45
46#  include "domzgr_substitute.h90"
47#  include "vectopt_loop_substitute.h90"
48
49   PRIVATE
50
51   PUBLIC   trc_sms_fabm            ! called by trcsms.F90 module
52   PUBLIC   trc_sms_fabm_alloc      ! called by trcini_fabm.F90 module
53   PUBLIC   trc_sms_fabm_check_mass ! called by trcwri_fabm.F90
54
55   ! Work arrays
56   REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: flux          ! Cross-interface flux of pelagic variables (# m-2 s-1)
57   REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: current_total ! Totals of conserved quantities
58
59   ! Arrays for environmental variables that are computed by the coupler
60   REAL(wp), PUBLIC, TARGET, ALLOCATABLE, DIMENSION(:,:,:) :: prn,rho
61   REAL(wp), PUBLIC, TARGET, ALLOCATABLE, DIMENSION(:,:) :: taubot
62   REAL(wp), PUBLIC, TARGET :: daynumber_in_year
63   REAL(wp), PUBLIC, POINTER :: swrad(:,:,:)  ! pointer to ERSEM spectral heating term
64
65   ! State repair counters
66   INTEGER, SAVE :: repair_interior_count = 0
67   INTEGER, SAVE :: repair_surface_count  = 0
68   INTEGER, SAVE :: repair_bottom_count   = 0
69
70   ! Coupler parameters
71   INTEGER, PUBLIC :: nn_adv  ! Vertical advection scheme for sinking/floating/movement
72                              ! (1: 1st order upwind, 3: 3rd order TVD)
73
74   ! Flag indicating whether model%start has been called (will be done on-demand)
75   LOGICAL, SAVE :: started = .false.
76
77   !!----------------------------------------------------------------------
78   !! NEMO/TOP 3.3 , NEMO Consortium (2010)
79   !! $Id$
80   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
81   !!----------------------------------------------------------------------
82CONTAINS
83
84   SUBROUTINE trc_sms_fabm( kt )
85      !!----------------------------------------------------------------------
86      !!                     ***  trc_sms_fabm  ***
87      !!
88      !! ** Purpose :   main routine of FABM model
89      !!
90      !! ** Method  : -
91      !!----------------------------------------------------------------------
92      !
93      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
94      INTEGER :: jn, jk
95      REAL(wp), POINTER, DIMENSION(:,:,:) :: ztrfabm, pdat
96      REAL(wp), DIMENSION(jpi,jpj)    :: vint
97
98!!----------------------------------------------------------------------
99      !
100      IF( nn_timing == 1 )  CALL timing_start('trc_sms_fabm')
101      !
102      IF(lwp) WRITE(numout,*)
103      IF(lwp) WRITE(numout,'(a,i0,a,i4.4,a,i2.2,a,i2.2,a,i5,a)') &
104          ' trc_sms_fabm:  FABM model, iteration ',kt,' ', &
105          nyear,'-',nmonth,'-',nday,' ',nsec_day," secs"
106      IF(lwp) WRITE(numout,*) ' ~~~~~~~~~~~~~~'
107
108      IF (.NOT. started) CALL nemo_fabm_start
109
110      CALL update_inputs( kt )
111
112      CALL compute_fabm( kt )
113
114      CALL compute_vertical_movement( kt, nn_adv )
115
116      CALL st2d_fabm_nxt( kt )
117     
118      CALL asmdiags_fabm( kt )
119
120      IF( l_trdtrc )  CALL wrk_alloc( jpi, jpj, jpk, ztrfabm )
121
122      CALL trc_bc_read  ( kt )       ! tracers: surface and lateral Boundary Conditions
123      CALL trc_rnf_fabm ( kt ) ! River forcings
124
125      ! Send 3D diagnostics to output (these apply to time "n")
126      DO jn = 1, size(model%interior_diagnostic_variables)
127         IF (model%interior_diagnostic_variables(jn)%save) THEN
128            ! Save 3D field
129            pdat => model%get_interior_diagnostic_data(jn)
130            CALL iom_put(model%interior_diagnostic_variables(jn)%name, pdat)
131
132            ! Save depth integral if selected for output in XIOS
133            IF (iom_use(TRIM(model%interior_diagnostic_variables(jn)%name)//'_VINT')) THEN
134               vint = 0._wp
135               DO jk = 1, jpkm1
136                  vint = vint + pdat(:,:,jk) * fse3t(:,:,jk) * tmask(:,:,jk)
137               END DO
138               CALL iom_put(TRIM(model%interior_diagnostic_variables(jn)%name)//'_VINT', vint)
139            END IF
140         END IF
141      END DO
142
143      ! Send 2D diagnostics to output (these apply to time "n")
144      DO jn = 1, size(model%horizontal_diagnostic_variables)
145         IF (model%horizontal_diagnostic_variables(jn)%save) &
146             CALL iom_put( model%horizontal_diagnostic_variables(jn)%name, model%get_horizontal_diagnostic_data(jn))
147      END DO
148
149      IF( l_trdtrc ) THEN      ! Save the trends in the mixed layer
150          DO jn = jp_fabm0, jp_fabm1
151            ztrfabm(:,:,:) = tra(:,:,:,jn)
152            CALL trd_trc( ztrfabm, jn, jptra_sms, kt )   ! save trends
153          END DO
154          CALL wrk_dealloc( jpi, jpj, jpk, ztrfabm )
155      END IF
156      !
157      IF( nn_timing == 1 )  CALL timing_stop('trc_sms_fabm')
158
159   END SUBROUTINE trc_sms_fabm
160   
161   SUBROUTINE asmdiags_fabm( kt )
162      INTEGER, INTENT(IN) :: kt
163      INTEGER :: ji,jj,jk,jkmax
164      REAL(wp), DIMENSION(jpi,jpj) :: zmld
165      REAL(wp), DIMENSION(jpi,jpj,jpk) :: pgrow_3d, ploss_3d
166     
167      IF (kt == nittrc000) THEN
168         MLD_MAX(:,:) = 0.0
169      ENDIF
170      PGROW_AVG(:,:) = 0.0
171      PLOSS_AVG(:,:) = 0.0
172      PHYT_AVG(:,:)  = 0.0
173       
174      pgrow_3d(:,:,:) = model%get_interior_diagnostic_data(jp_fabm_pgrow)
175      ploss_3d(:,:,:) = model%get_interior_diagnostic_data(jp_fabm_ploss)
176     
177      SELECT CASE( mld_choice_bgc )
178      CASE ( 1 )                   ! Turbocline/mixing depth [W points]
179         zmld(:,:) = hmld(:,:)
180      CASE ( 2 )                   ! Density criterion (0.01 kg/m^3 change from 10m) [W points]
181         zmld(:,:) = hmlp(:,:)
182      CASE ( 3 )                   ! Kara MLD [Interpolated]
183#if defined key_karaml
184         IF ( ln_kara ) THEN
185            zmld(:,:) = hmld_kara(:,:)
186         ELSE
187            CALL ctl_stop( ' Kara mixed layer requested for BGC assimilation,', &
188               &           ' but ln_kara=.false.' )
189         ENDIF
190#else
191         CALL ctl_stop( ' Kara mixed layer requested for BGC assimilation,', &
192            &           ' but is not defined' )
193#endif
194      CASE ( 4 )                   ! Temperature criterion (0.2 K change from surface) [T points]
195         zmld(:,:) = hmld_tref(:,:)
196      CASE ( 5 )                   ! Density criterion (0.01 kg/m^3 change from 10m) [T points]
197         zmld(:,:) = hmlpt(:,:)
198      END SELECT
199   
200      DO jj = 2, jpjm1
201         DO ji = 2, jpim1
202            !
203            jkmax = jpk-1
204            DO jk = jpk-1, 1, -1
205               IF ( ( zmld(ji,jj) >  gdepw_n(ji,jj,jk)   ) .AND. &
206                  & ( zmld(ji,jj) <= gdepw_n(ji,jj,jk+1) ) ) THEN
207                  zmld(ji,jj) = gdepw_n(ji,jj,jk+1)
208                  jkmax = jk
209               ENDIF
210            END DO
211            !
212            DO jk = 1, jkmax
213               PHYT_AVG(ji,jj) = PHYT_AVG(ji,jj) + &
214                  &              trn(ji,jj,jk,jp_fabm_m1+jp_fabm_p1n) + &
215                  &              trn(ji,jj,jk,jp_fabm_m1+jp_fabm_p2n) + &
216                  &              trn(ji,jj,jk,jp_fabm_m1+jp_fabm_p3n) + &
217                  &              trn(ji,jj,jk,jp_fabm_m1+jp_fabm_p4n)
218               IF ( pgrow_3d(ji,jj,jk) .GT. 0.0 ) THEN
219                  PGROW_AVG(ji,jj) = PGROW_AVG(ji,jj) + &
220                     &               pgrow_3d(ji,jj,jk)
221               ENDIF
222               IF ( ploss_3d(ji,jj,jk) .GT. 0.0 ) THEN
223                  PLOSS_AVG(ji,jj) = PLOSS_AVG(ji,jj) + &
224                     &               ploss_3d(ji,jj,jk)
225               ENDIF
226            END DO
227           
228            PHYT_AVG(ji,jj)  = PHYT_AVG(ji,jj)  / REAL(jkmax)
229            PGROW_AVG(ji,jj) = PGROW_AVG(ji,jj) / REAL(jkmax)
230            PLOSS_AVG(ji,jj) = PLOSS_AVG(ji,jj) / REAL(jkmax)
231   
232            IF ( zmld(ji,jj) .GT. MLD_MAX(ji,jj) ) THEN
233               MLD_MAX(ji,jj) = zmld(ji,jj)
234            ENDIF
235            !
236         END DO
237      END DO
238     
239      PHYT_AVG(:,:)  = PHYT_AVG(:,:)  * tmask(:,:,1)
240      PGROW_AVG(:,:) = PGROW_AVG(:,:) * tmask(:,:,1)
241      PLOSS_AVG(:,:) = PLOSS_AVG(:,:) * tmask(:,:,1)
242      MLD_MAX(:,:)   = MLD_MAX(:,:)   * tmask(:,:,1)
243     
244   END SUBROUTINE asmdiags_fabm
245
246   SUBROUTINE compute_fabm( kt )
247      INTEGER, INTENT(in) :: kt   ! ocean time-step index
248
249      INTEGER :: ji,jj,jk,jn
250      LOGICAL :: valid, repaired
251      REAL(wp) :: zalfg,zztmpx,zztmpy
252
253      ! Validate current model state (setting argument to .TRUE. enables repair=clipping)
254      CALL check_state(.TRUE., valid, repaired)
255      IF (.NOT. valid) THEN
256         WRITE(numout,*) "Invalid value in FABM encountered in area ",narea,"!!!"
257#if defined key_iomput
258         CALL xios_finalize                ! end mpp communications with xios
259         IF( lk_oasis ) CALL cpl_finalize    ! end coupling and mpp communications with OASIS
260#else
261         IF( lk_oasis ) THEN
262            CALL cpl_finalize              ! end coupling and mpp communications with OASIS
263         ELSE
264            IF( lk_mpp )   CALL mppstop    ! end mpp communications
265         ENDIF
266#endif
267      END IF
268      IF (repaired) THEN
269         WRITE(numout,*) "Total interior repairs up to now on process",narea,":",repair_interior_count
270         WRITE(numout,*) "Total surface repairs up to now on process",narea,":",repair_surface_count
271         WRITE(numout,*) "Total bottom repairs up to now on process",narea,":",repair_bottom_count
272      ENDIF
273
274      daynumber_in_year = fjulday - fjulstartyear + 1
275
276      ! Compute the now hydrostatic pressure
277      ! copied from istate.F90
278      ! ------------------------------------
279
280      IF (ALLOCATED(rho)) rho = rau0 * ( 1._wp + rhd )
281
282      IF (ALLOCATED(prn)) THEN
283         zalfg = 0.5e-4_wp * grav ! FABM wants dbar, convert from Pa (and multiply with 0.5 to average 2 cell thicknesses below)
284         prn(:,:,1) = 10.1325_wp + zalfg * fse3t(:,:,1) * rho(:,:,1)
285         DO jk = 2, jpkm1                                              ! Vertical integration from the surface
286            prn(:,:,jk) = prn(:,:,jk-1) + zalfg * ( &
287                        fse3t(:,:,jk-1) * rho(:,:,jk-1)  &
288                        + fse3t(:,:,jk) * rho(:,:,jk) )
289         END DO
290      END IF
291
292      ! Compute the bottom stress
293      ! copied from diawri.F90
294      ! ------------------------------------
295
296      IF (ALLOCATED(taubot)) THEN
297         taubot(:,:) = 0._wp
298         DO jj = 2, jpjm1
299            DO ji = fs_2, fs_jpim1   ! vector opt.
300                  zztmpx = (  bfrua(ji  ,jj) * un(ji  ,jj,mbku(ji  ,jj))  &
301                        &  + bfrua(ji-1,jj) * un(ji-1,jj,mbku(ji-1,jj))  )
302                  zztmpy = (  bfrva(ji,  jj) * vn(ji,jj  ,mbkv(ji,jj  ))  &
303                        &  + bfrva(ji,jj-1) * vn(ji,jj-1,mbkv(ji,jj-1))  )
304                  taubot(ji,jj) = 0.5_wp * rau0 * SQRT( zztmpx * zztmpx + zztmpy * zztmpy ) * tmask(ji,jj,1)
305                  !
306            END DO
307         END DO
308      END IF
309
310      CALL model%prepare_inputs(real(kt, wp),nyear,nmonth,nday,REAL(nsec_day,wp))
311
312      ! Retrieve 3D shortwave and store in etot3
313      IF (ln_qsr_spec) etot3(:,:,:) = swrad(:,:,:)
314
315      ! Zero rate array of interface-attached state variables
316      fabm_st2Da = 0._wp
317
318      ! Compute interfacial source terms and fluxes
319      DO jj=2,jpjm1
320         ! Process bottom (get_bottom_sources increments rather than sets, so zero flux array first)
321         flux = 0._wp
322         CALL model%get_bottom_sources(fs_2,fs_jpim1,jj,flux,fabm_st2Da(fs_2:fs_jpim1,jj,jp_fabm_surface+1:))
323         DO jn=1,jp_fabm
324            ! Divide bottom fluxes by height of bottom layer and add to source terms.
325            DO ji=fs_2,fs_jpim1
326               tra(ji,jj,mbkt(ji,jj),jp_fabm_m1+jn) = tra(ji,jj,mbkt(ji,jj),jp_fabm_m1+jn) + flux(ji,jn)/fse3t(ji,jj,mbkt(ji,jj))
327            END DO
328         END DO
329
330         ! Process surface (get_surface_sources increments rather than sets, so zero flux array first)
331         flux = 0._wp
332         CALL model%get_surface_sources(fs_2,fs_jpim1,jj,flux,fabm_st2Da(fs_2:fs_jpim1,jj,1:jp_fabm_surface))
333         ! Divide surface fluxes by height of surface layer and add to source terms.
334         DO jn=1,jp_fabm
335            DO ji=fs_2,fs_jpim1
336               tra(ji,jj,1,jp_fabm_m1+jn) = tra(ji,jj,1,jp_fabm_m1+jn) + flux(ji,jn)/fse3t(ji,jj,1)
337            END DO
338         END DO
339      END DO
340
341      ! Compute interior source terms (NB get_interior_sources increments rather than sets)
342      DO jk=1,jpkm1
343         DO jj=2,jpjm1
344            CALL model%get_interior_sources(fs_2,fs_jpim1,jj,jk,tra(fs_2:fs_jpim1,jj,jk,jp_fabm0:jp_fabm1))
345         END DO
346      END DO
347
348      CALL model%finalize_outputs()
349   END SUBROUTINE compute_fabm
350
351   SUBROUTINE check_state(repair, valid, repaired)
352      LOGICAL, INTENT(IN)  :: repair
353      LOGICAL, INTENT(OUT) :: valid, repaired
354
355      INTEGER :: jj, jk
356      LOGICAL :: valid_int, valid_sf, valid_bt
357
358      valid = .TRUE.     ! Whether the model state is valid after this subroutine returns
359      repaired = .FALSE. ! Whether the model state has been repaired by this subroutine
360      DO jk=1,jpkm1
361         DO jj=2,jpjm1
362            CALL model%check_interior_state(fs_2, fs_jpim1, jj, jk, repair, valid_int)
363            IF (repair .AND. .NOT. valid_int) THEN
364               repair_interior_count = repair_interior_count + 1
365               repaired = .TRUE.
366            END IF
367            IF (.NOT. (valid_int .OR. repair)) valid = .FALSE.
368         END DO
369      END DO
370      DO jj=2,jpjm1
371         CALL model%check_surface_state(fs_2, fs_jpim1, jj, repair, valid_sf)
372         IF (repair .AND. .NOT. valid_sf) THEN
373            repair_surface_count = repair_surface_count + 1
374            repaired = .TRUE.
375         END IF
376         IF (.NOT. (valid_sf .AND. valid_bt) .AND. .NOT. repair) valid = .FALSE.
377         CALL model%check_bottom_state(fs_2, fs_jpim1, jj, repair, valid_bt)
378         IF (repair .AND. .NOT. valid_bt) THEN
379            repair_bottom_count = repair_bottom_count + 1
380            repaired = .TRUE.
381         END IF
382         IF (.NOT. (valid_sf .AND. valid_bt) .AND. .NOT. repair) valid = .FALSE.
383      END DO
384   END SUBROUTINE
385
386   SUBROUTINE trc_sms_fabm_check_mass()
387      REAL(wp) :: total(SIZE(model%conserved_quantities))
388      INTEGER :: ji,jk,jj,jn
389
390      total = 0._wp
391
392      IF (.NOT. started) CALL nemo_fabm_start
393
394      DO jk=1,jpkm1
395         DO jj=2,jpjm1
396            CALL model%get_interior_conserved_quantities(fs_2,fs_jpim1,jj,jk,current_total)
397            DO jn=1,SIZE(model%conserved_quantities)
398               DO ji=fs_2,fs_jpim1
399                  total(jn) = total(jn) + cvol(ji,jj,jk) * current_total(ji,jn) * tmask_i(ji,jj)
400               END DO
401            END DO
402         END DO
403      END DO
404
405      DO jj=2,jpjm1
406         CALL model%get_horizontal_conserved_quantities(fs_2,fs_jpim1,jj,current_total)
407         DO jn=1,SIZE(model%conserved_quantities)
408            DO ji=fs_2,fs_jpim1
409               total(jn) = total(jn) + e1e2t(ji,jj) * current_total(ji,jn) * tmask_i(ji,jj)
410            END DO
411         END DO
412      END DO
413
414      IF( lk_mpp ) CALL mpp_sum(total,SIZE(model%conserved_quantities))
415
416      DO jn=1,SIZE(model%conserved_quantities)
417         IF(lwp) WRITE(numout,*) 'FABM '//TRIM(model%conserved_quantities(jn)%name),total(jn),TRIM(model%conserved_quantities(jn)%units)//'*m3'
418      END DO
419
420   END SUBROUTINE trc_sms_fabm_check_mass
421
422   SUBROUTINE st2d_fabm_nxt( kt )
423      !!----------------------------------------------------------------------
424      !!                     ***  st2d_fabm_nxt  ***
425      !!
426      !! ** Purpose :   routine to integrate 2d states in time
427      !!
428      !! ** Method  :   based on integration of 3D passive tracer fields
429      !!                implemented in TOP_SRC/TRP/trcnxt.F90, plus
430      !!                tra_nxt_fix in OPA_SRC/TRA/tranxt.F90. Similar to
431      !!                time integration of sea surface height in
432      !!                OPA_SRC/DYN/sshwzv.F90.
433      !!----------------------------------------------------------------------
434      !
435      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
436      REAL(wp) :: z2dt
437      INTEGER :: jn
438
439!!----------------------------------------------------------------------
440      !
441      IF( neuler == 0 .AND. kt == nittrc000 ) THEN
442          z2dt = rdt                  ! set time step size (Euler)
443      ELSE
444          z2dt = 2._wp * rdt          ! set time step size (Leapfrog)
445      ENDIF
446
447      ! Forward Euler time step to compute "now"
448      DO jn=1,jp_fabm_surface+jp_fabm_bottom
449         fabm_st2Da(:,:,jn) = (fabm_st2db(:,:,jn) + z2dt * fabm_st2da(:,:,jn)) * tmask(:,:,1)
450      ENDDO
451
452      IF( neuler == 0 .AND. kt == nittrc000 ) THEN        ! Euler time-stepping at first time-step
453         !                                                ! (only swap)
454         fabm_st2Dn(:,:,:) = fabm_st2Da(:,:,:)
455         !
456      ELSE
457         ! Update now state + Asselin filter time stepping
458         fabm_st2Db(:,:,:) = (1._wp - 2._wp*atfp) * fabm_st2Dn(:,:,:) + &
459             atfp * ( fabm_st2Db(:,:,:) + fabm_st2Da(:,:,:) )
460         fabm_st2Dn(:,:,:) = fabm_st2Da(:,:,:)
461      ENDIF
462
463   END SUBROUTINE st2d_fabm_nxt
464
465   INTEGER FUNCTION trc_sms_fabm_alloc()
466      INTEGER :: jn
467      !!----------------------------------------------------------------------
468      !!              ***  ROUTINE trc_sms_fabm_alloc  ***
469      !!----------------------------------------------------------------------
470      !
471      ! ALLOCATE here the arrays specific to FABM
472      ALLOCATE( lk_rad_fabm(jp_fabm))
473      IF (model%variable_needs_values(fabm_standard_variables%pressure)) ALLOCATE(prn(jpi, jpj, jpk))
474      IF (ALLOCATED(prn) .or. model%variable_needs_values(fabm_standard_variables%density)) ALLOCATE(rho(jpi, jpj, jpk))
475      IF (model%variable_needs_values(fabm_standard_variables%bottom_stress)) ALLOCATE(taubot(jpi, jpj))
476      ! ALLOCATE( tab(...) , STAT=trc_sms_fabm_alloc )
477
478      ! Allocate arrays to hold state for surface-attached and bottom-attached state variables
479      ALLOCATE(fabm_st2Dn(jpi, jpj, jp_fabm_surface+jp_fabm_bottom))
480      ALLOCATE(fabm_st2Da(jpi, jpj, jp_fabm_surface+jp_fabm_bottom))
481      ALLOCATE(fabm_st2Db(jpi, jpj, jp_fabm_surface+jp_fabm_bottom))
482
483      ! Work array to hold surface and bottom fluxes
484      ALLOCATE(flux(fs_2:fs_jpim1,jp_fabm))
485
486      ! Allocate work arrays for vertical movement
487      ALLOCATE(w_ct(fs_2:fs_jpim1,1:jpkm1,jp_fabm))
488      ALLOCATE(current_total(fs_2:fs_jpim1,SIZE(model%conserved_quantities)))
489#if defined key_trdtrc && defined key_iomput
490      IF( lk_trdtrc ) ALLOCATE(tr_vmv(jpi,jpj,jpk,jp_fabm))
491      IF( lk_trdtrc ) ALLOCATE(tr_inp(jpi,jpj,jpk))
492#endif
493
494      trc_sms_fabm_alloc = 0      ! set to zero if no array to be allocated
495      !
496      IF( trc_sms_fabm_alloc /= 0 ) CALL ctl_warn('trc_sms_fabm_alloc : failed to allocate arrays')
497      !
498
499      ! Provide FABM with domain extents
500      CALL model%set_domain(jpi, jpj, jpk, rdt)
501      CALL model%set_domain_start(fs_2, 2, 1)
502      CALL model%set_domain_stop(fs_jpim1, jpjm1, jpkm1)
503
504      ! Provide FABM with the vertical indices of the bottom, and the land-sea mask.
505      CALL model%set_bottom_index(mbkt)  ! NB mbkt extents should match dimension lengths provided to set_domain
506      CALL model%set_mask(tmask,tmask(:,:,1)) ! NB tmask extents should match dimension lengths provided to set_domain
507
508      ! Initialize state and send pointers to state data to FABM
509      ! We mask land points in states with zeros, as per with NEMO "convention"
510      ! NB we cannot call model%initialize_*_state at this point, because model%start has not been called yet.
511      DO jn=1,jp_fabm
512         trn(:,:,:,jp_fabm_m1+jn) = model%interior_state_variables(jn)%initial_value * tmask
513         CALL model%link_interior_state_data(jn,trn(:,:,:,jp_fabm_m1+jn))
514      END DO
515      DO jn=1,jp_fabm_surface
516         fabm_st2Dn(:,:,jn) = model%surface_state_variables(jn)%initial_value * tmask(:,:,1)
517         CALL model%link_surface_state_data(jn,fabm_st2Dn(:,:,jn))
518      END DO
519      DO jn=1,jp_fabm_bottom
520         fabm_st2Dn(:,:,jp_fabm_surface+jn) = model%bottom_state_variables(jn)%initial_value * tmask(:,:,1)
521         CALL model%link_bottom_state_data(jn,fabm_st2Dn(:,:,jp_fabm_surface+jn))
522      END DO
523
524      ! Send pointers to environmental data to FABM
525      CALL model%link_interior_data(fabm_standard_variables%depth, fsdept(:,:,:))
526      CALL model%link_interior_data(fabm_standard_variables%temperature, tsn(:,:,:,jp_tem))
527      CALL model%link_interior_data(fabm_standard_variables%practical_salinity, tsn(:,:,:,jp_sal))
528      IF (ALLOCATED(rho)) CALL model%link_interior_data(fabm_standard_variables%density, rho(:,:,:))
529      IF (ALLOCATED(prn)) CALL model%link_interior_data(fabm_standard_variables%pressure, prn)
530      IF (ALLOCATED(taubot)) CALL model%link_horizontal_data(fabm_standard_variables%bottom_stress, taubot(:,:))
531      CALL model%link_interior_data(fabm_standard_variables%cell_thickness, fse3t(:,:,:))
532      CALL model%link_horizontal_data(fabm_standard_variables%latitude, gphit)
533      CALL model%link_horizontal_data(fabm_standard_variables%longitude, glamt)
534      CALL model%link_scalar(fabm_standard_variables%number_of_days_since_start_of_the_year, daynumber_in_year)
535      CALL model%link_horizontal_data(fabm_standard_variables%wind_speed, wndm(:,:))
536      CALL model%link_horizontal_data(fabm_standard_variables%surface_downwelling_shortwave_flux, qsr(:,:))
537      CALL model%link_horizontal_data(fabm_standard_variables%bottom_depth_below_geoid, bathy(:,:))
538      CALL model%link_horizontal_data(fabm_standard_variables%ice_area_fraction, fr_i(:,:))
539
540      ! Obtain user-specified input variables (read from NetCDF file)
541      CALL link_inputs
542      CALL update_inputs(nit000, .FALSE.)
543
544      ! Set mask for negativity corrections to the relevant states
545      lk_rad_fabm(:) = .FALSE.
546      DO jn=1,jp_fabm
547        IF (model%interior_state_variables(jn)%minimum >= 0._wp) THEN
548          lk_rad_fabm(jn) = .TRUE.
549          IF(lwp) WRITE(numout,*) 'FABM clipping for '//TRIM(model%interior_state_variables(jn)%name)//' activated.'
550        END IF
551      END DO
552
553
554      ! Copy initial condition for interface-attached state variables to "previous" state field
555      ! NB NEMO does this itself for pelagic state variables (trb) in TOP_SRC/trcini.F90.
556      fabm_st2Db = fabm_st2Dn
557
558      ! Pointer to spectral heating term
559      swrad => model%get_data(model%get_interior_variable_id(standard_variables%net_rate_of_absorption_of_shortwave_energy_in_layer))
560
561
562   END FUNCTION trc_sms_fabm_alloc
563
564   SUBROUTINE nemo_fabm_start()
565      INTEGER :: jn
566
567      ! Make FABM aware of diagnostics that are not needed [not included in output]
568      ! This works only after iom has completely initialised, because it depends on iom_use
569      DO jn=1,size(model%interior_diagnostic_variables)
570         model%interior_diagnostic_variables(jn)%save = iom_use(model%interior_diagnostic_variables(jn)%name) &
571            .or. iom_use(TRIM(model%interior_diagnostic_variables(jn)%name)//'_VINT')
572      END DO
573      DO jn=1,size(model%horizontal_diagnostic_variables)
574         model%horizontal_diagnostic_variables(jn)%save = iom_use(model%horizontal_diagnostic_variables(jn)%name)
575      END DO
576
577      ! Check whether FABM has all required data
578      ! [after this, the save attribute of diagnostic variables can no longer change!]
579      CALL model%start()
580
581      started = .TRUE.
582   END SUBROUTINE
583
584#else
585   !!----------------------------------------------------------------------
586   !!   Dummy module                                        No FABM model
587   !!----------------------------------------------------------------------
588CONTAINS
589   SUBROUTINE trc_sms_fabm( kt )             ! Empty routine
590      INTEGER, INTENT( in ) ::   kt
591      WRITE(*,*) 'trc_sms_fabm: You should not have seen this print! error?', kt
592   END SUBROUTINE trc_sms_fabm
593#endif
594
595   !!======================================================================
596END MODULE trcsms_fabm
Note: See TracBrowser for help on using the repository browser.