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 @ 13451

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

Modifications to accommodate interaction between NEMO-FABM coupler only saving diagnostics requested in iodef.xml, and the assimilation and 25h/tmb diagnostic code.

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