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/CO6_KD490_amm7_oper_fabm_hem08/NEMOGCM/NEMO/TOP_SRC/FABM – NEMO

source: branches/UKMO/CO6_KD490_amm7_oper_fabm_hem08/NEMOGCM/NEMO/TOP_SRC/FABM/trcsms_fabm.F90 @ 9330

Last change on this file since 9330 was 9330, checked in by dford, 6 years ago

Add required diagnostics.

File size: 20.2 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   !!----------------------------------------------------------------------
8#if defined key_fabm
9   !!----------------------------------------------------------------------
10   !!   'key_fabm'                                               FABM tracers
11   !!----------------------------------------------------------------------
12   !! trc_sms_fabm       : FABM model main routine
13   !! trc_sms_fabm_alloc : allocate arrays specific to FABM sms
14   !!----------------------------------------------------------------------
15   USE par_trc         ! TOP parameters
16   USE oce_trc         ! Ocean variables
17   USE trc             ! TOP variables
18   USE trcbc
19   USE trd_oce
20   USE trdtrc
21#if defined key_trdtrc && defined key_iomput
22   USE trdtrc_oce
23#endif
24
25   USE oce, only: tsn  ! Needed?
26   USE sbc_oce, only: lk_oasis
27   USE dom_oce
28   USE zdf_oce
29   !USE iom
30   USE xios
31   USE cpl_oasis3
32   USE st2D_fabm
33   USE inputs_fabm
34   USE vertical_movement_fabm
35   USE zdfmxl
36
37   !USE fldread         !  time interpolation
38
39   USE fabm
40
41   IMPLICIT NONE
42
43#  include "domzgr_substitute.h90"
44#  include "vectopt_loop_substitute.h90"
45
46   PRIVATE
47
48   PUBLIC   trc_sms_fabm       ! called by trcsms.F90 module
49   PUBLIC   trc_sms_fabm_alloc ! called by trcini_fabm.F90 module
50   PUBLIC   trc_sms_fabm_check_mass
51   PUBLIC   st2d_fabm_nxt ! 2D state intergration
52   PUBLIC   compute_fabm ! Compute FABM sources, sinks and diagnostics
53
54   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: flux    ! Cross-interface flux of pelagic variables (# m-2 s-1)
55
56   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:)   :: ext     ! Light extinction coefficient (m-1)
57
58   ! Work array for mass aggregation
59   REAL(wp), ALLOCATABLE, SAVE, TARGET, DIMENSION(:,:)   :: current_total
60
61
62   ! Arrays for environmental variables
63   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, TARGET, DIMENSION(:,:,:) :: prn,rho
64   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, TARGET, DIMENSION(:,:) :: taubot
65
66   ! repair counters
67   INTEGER :: repair_interior_count,repair_surface_count,repair_bottom_count
68
69   ! state check type
70   TYPE type_state
71      LOGICAL             :: valid
72      LOGICAL             :: repaired
73   END TYPE
74
75   REAL(wp), PUBLIC :: daynumber_in_year
76
77   TYPE (type_bulk_variable_id),SAVE :: swr_id
78
79   !!----------------------------------------------------------------------
80   !! NEMO/TOP 3.3 , NEMO Consortium (2010)
81   !! $Id$
82   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
83   !!----------------------------------------------------------------------
84CONTAINS
85
86   SUBROUTINE trc_sms_fabm( kt )
87      !!----------------------------------------------------------------------
88      !!                     ***  trc_sms_fabm  ***
89      !!
90      !! ** Purpose :   main routine of FABM model
91      !!
92      !! ** Method  : -
93      !!----------------------------------------------------------------------
94      !
95      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
96      INTEGER :: jn
97      REAL(wp), POINTER, DIMENSION(:,:,:) :: ztrfabm
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
109      CALL update_inputs( kt )
110
111      CALL compute_fabm
112
113      CALL compute_vertical_movement( kt )
114
115      CALL st2d_fabm_nxt( kt )
116     
117      CALL asmdiags_fabm
118
119      IF( l_trdtrc )  CALL wrk_alloc( jpi, jpj, jpk, ztrfabm )
120
121      CALL trc_bc_read  ( kt )       ! tracers: surface and lateral Boundary Conditions
122      CALL trc_rnf_fabm ( kt ) ! River forcings
123
124      IF( l_trdtrc ) THEN      ! Save the trends in the mixed layer
125          DO jn = jp_fabm0, jp_fabm1
126            ztrfabm(:,:,:) = tra(:,:,:,jn)
127            CALL trd_trc( ztrfabm, jn, jptra_sms, kt )   ! save trends
128          END DO
129          CALL wrk_dealloc( jpi, jpj, jpk, ztrfabm )
130      END IF
131      !
132      IF( nn_timing == 1 )  CALL timing_stop('trc_sms_fabm')
133
134   END SUBROUTINE trc_sms_fabm
135   
136   SUBROUTINE asmdiags_fabm
137      INTEGER :: ji,jj,jk,jkmax
138   
139      PGROW_AVG(:,:) = 1.0
140      PLOSS_AVG(:,:) = 1.0
141      PHYT_AVG(:,:)  = 0.0
142
143      DO jj = 1, jpj
144         DO ji = 1, jpi
145            !
146            jkmax = jpk-1
147            DO jk = jpk-1, 1, -1
148               IF ( ( hmld_tref(ji,jj) >  gdepw_n(ji,jj,jk)   ) .AND. &
149                  & ( hmld_tref(ji,jj) <= gdepw_n(ji,jj,jk+1) ) ) THEN
150                  hmld_tref(ji,jj) = gdepw_n(ji,jj,jk+1)
151                  jkmax = jk
152               ENDIF
153            END DO
154            !
155            DO jk = 1, jkmax
156               PHYT_AVG(ji,jj) = PHYT_AVG(ji,jj) + &
157                  &              trn(ji,jj,jk,jp_fabm_p1n) + &
158                  &              trn(ji,jj,jk,jp_fabm_p2n) + &
159                  &              trn(ji,jj,jk,jp_fabm_p3n) + &
160                  &              trn(ji,jj,jk,jp_fabm_p4n)
161            END DO
162           
163            PHYT_AVG(ji,jj) = PHYT_AVG(ji,jj) / REAL(jkmax)
164
165            IF ( hmld_tref(ji,jj) .GT. MLD_MAX(ji,jj) ) THEN
166               MLD_MAX(ji,jj) = hmld_tref(ji,jj)
167            ENDIF
168            !
169         END DO
170      END DO
171   
172   END SUBROUTINE asmdiags_fabm
173
174   SUBROUTINE compute_fabm()
175      INTEGER :: ji,jj,jk,jn
176      TYPE(type_state) :: valid_state
177      REAL(wp) :: zalfg,zztmpx,zztmpy
178
179      ! Validate current model state (setting argument to .TRUE. enables repair=clipping)
180      valid_state = check_state(.TRUE.)
181      IF (.NOT.valid_state%valid) THEN
182         WRITE(numout,*) "Invalid value in FABM encountered in area ",narea,"!!!"
183#if defined key_iomput
184         CALL xios_finalize                ! end mpp communications with xios
185         IF( lk_oasis ) CALL cpl_finalize    ! end coupling and mpp communications with OASIS
186#else
187         IF( lk_oasis ) THEN
188            CALL cpl_finalize              ! end coupling and mpp communications with OASIS
189         ELSE
190            IF( lk_mpp )   CALL mppstop    ! end mpp communications
191         ENDIF
192#endif
193      END IF
194      IF (valid_state%repaired) THEN
195         WRITE(numout,*) "Total interior repairs up to now on process",narea,":",repair_interior_count
196         WRITE(numout,*) "Total surface repairs up to now on process",narea,":",repair_surface_count
197         WRITE(numout,*) "Total bottom repairs up to now on process",narea,":",repair_bottom_count
198      ENDIF
199
200      ! Compute the now hydrostatic pressure
201      ! copied from istate.F90
202      ! ------------------------------------
203
204      zalfg = 0.5e-4 * grav ! FABM wants dbar, convert from Pa
205
206      rho = rau0 * ( 1. + rhd )
207
208      prn(:,:,1) = 10.1325 + zalfg * fse3t(:,:,1) * rho(:,:,1)
209
210      daynumber_in_year=(fjulday-fjulstartyear+1)*1._wp
211
212      DO jk = 2, jpk                                              ! Vertical integration from the surface
213         prn(:,:,jk) = prn(:,:,jk-1) + zalfg * ( &
214                     fse3t(:,:,jk-1) * rho(:,:,jk-1)  &
215                     + fse3t(:,:,jk) * rho(:,:,jk) )
216      END DO
217
218      ! Bottom stress
219      taubot(:,:) = 0._wp
220      DO jj = 2, jpjm1
221         DO ji = fs_2, fs_jpim1   ! vector opt.
222               zztmpx = (  bfrua(ji  ,jj) * un(ji  ,jj,mbku(ji  ,jj))  &
223                      &  + bfrua(ji-1,jj) * un(ji-1,jj,mbku(ji-1,jj))  )
224               zztmpy = (  bfrva(ji,  jj) * vn(ji,jj  ,mbkv(ji,jj  ))  &
225                      &  + bfrva(ji,jj-1) * vn(ji,jj-1,mbkv(ji,jj-1))  )
226               taubot(ji,jj) = 0.5_wp * rau0 * SQRT( zztmpx * zztmpx + zztmpy * zztmpy ) * tmask(ji,jj,1)
227               !
228         ENDDO
229      ENDDO
230      ! Compute light extinction
231      DO jk=1,jpk
232          DO jj=1,jpj
233            call fabm_get_light_extinction(model,1,jpi,jj,jk,ext)
234         END DO
235      END DO
236
237      ! Compute light field (stored among FABM's internal diagnostics)
238      DO jj=1,jpj
239          DO ji=1,jpi
240            call fabm_get_light(model,1,jpk,ji,jj)
241         END DO
242      END DO
243
244      ! TODO: retrieve 3D shortwave and store in etot3
245
246      ! Zero rate array of interface-attached state variables
247      fabm_st2Da = 0._wp
248
249      ! Compute interfacial source terms and fluxes
250      DO jj=1,jpj
251         ! Process bottom (fabm_do_bottom increments rather than sets, so zero flux array first)
252         flux = 0._wp
253         CALL fabm_do_bottom(model,1,jpi,jj,flux,fabm_st2Da(:,jj,jp_fabm_surface+1:))
254         DO jn=1,jp_fabm
255             DO ji=1,jpi
256                 ! Divide bottom fluxes by height of bottom layer and add to source terms.
257                 ! TODO: is there perhaps an existing variable for fse3t(ji,jj,mbkt(ji,jj))??
258                 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))
259             END DO
260         END DO
261
262         ! Process surface (fabm_do_surface increments rather than sets, so zero flux array first)
263         flux = 0._wp
264         CALL fabm_do_surface(model,1,jpi,jj,flux,fabm_st2Da(:,jj,1:jp_fabm_surface))
265         DO jn=1,jp_fabm
266             ! Divide surface fluxes by height of surface layer and add to source terms.
267             tra(:,jj,1,jp_fabm_m1+jn) = tra(:,jj,1,jp_fabm_m1+jn) + flux(:,jn)/fse3t(:,jj,1)
268         END DO
269      END DO
270
271      ! Compute interior source terms (NB fabm_do increments rather than sets)
272      DO jk=1,jpk
273          DO jj=1,jpj
274              CALL fabm_do(model,1,jpi,jj,jk,tra(:,jj,jk,jp_fabm0:jp_fabm1))
275          END DO
276      END DO
277   END SUBROUTINE compute_fabm
278
279   FUNCTION check_state(repair) RESULT(exit_state)
280      LOGICAL, INTENT(IN) :: repair
281      TYPE(type_state) :: exit_state
282
283      INTEGER             :: jj,jk
284      LOGICAL             :: valid_int,valid_sf,valid_bt
285
286      exit_state%valid = .TRUE.
287      exit_state%repaired =.FALSE.
288      DO jk=1,jpk
289         DO jj=1,jpj
290            CALL fabm_check_state(model,1,jpi,jj,jk,repair,valid_int)
291            IF (repair.AND..NOT.valid_int) THEN
292               repair_interior_count = repair_interior_count + 1
293               exit_state%repaired = .TRUE.
294            END IF
295            IF (.NOT.(valid_int.OR.repair)) exit_state%valid = .FALSE.
296         END DO
297      END DO
298      DO jj=1,jpj
299         CALL fabm_check_surface_state(model,1,jpi,jj,repair,valid_sf)
300         IF (repair.AND..NOT.valid_sf) THEN
301            repair_surface_count = repair_surface_count + 1
302            exit_state%repaired = .TRUE.
303         END IF
304         IF (.NOT.(valid_sf.AND.valid_bt).AND..NOT.repair) exit_state%valid = .FALSE.
305         CALL fabm_check_bottom_state(model,1,jpi,jj,repair,valid_bt)
306         IF (repair.AND..NOT.valid_bt) THEN
307            repair_bottom_count = repair_bottom_count + 1
308            exit_state%repaired = .TRUE.
309         END IF
310         IF (.NOT.(valid_sf.AND.valid_bt).AND..NOT.repair) exit_state%valid = .FALSE.
311      END DO
312   END FUNCTION
313
314   SUBROUTINE trc_sms_fabm_check_mass()
315      REAL(wp) :: total(SIZE(model%conserved_quantities))
316      INTEGER :: jk,jj,jn
317
318      total = 0._wp
319
320      DO jk=1,jpk
321         DO jj=1,jpj
322            CALL fabm_get_conserved_quantities(model,1,jpi,jj,jk,current_total)
323            DO jn=1,SIZE(model%conserved_quantities)
324               total(jn) = total(jn) + SUM(cvol(:,jj,jk)*current_total(:,jn)*tmask_i(:,jj))
325            END DO
326         END DO
327      END DO
328
329      DO jj=1,jpj
330         CALL fabm_get_horizontal_conserved_quantities(model,1,jpi,jj,current_total)
331         DO jn=1,SIZE(model%conserved_quantities)
332            total(jn) = total(jn) + SUM(e1e2t(:,jj)*current_total(:,jn)*tmask_i(:,jj))
333         END DO
334      END DO
335
336      IF( lk_mpp ) CALL mpp_sum(total,SIZE(model%conserved_quantities))
337
338      DO jn=1,SIZE(model%conserved_quantities)
339         IF(lwp) WRITE(numout,*) 'FABM '//TRIM(model%conserved_quantities(jn)%name),total(jn),TRIM(model%conserved_quantities(jn)%units)//'*m3'
340      END DO
341
342   END SUBROUTINE trc_sms_fabm_check_mass
343
344   SUBROUTINE st2d_fabm_nxt( kt )
345      !!----------------------------------------------------------------------
346      !!                     ***  st2d_fabm_nxt  ***
347      !!
348      !! ** Purpose :   routine to integrate 2d states in time
349      !!
350      !! ** Method  :   based on integration of 3D passive tracer fields
351      !!                implemented in TOP_SRC/TRP/trcnxt.F90, plus
352      !!                tra_nxt_fix in OPA_SRC/TRA/tranxt.F90. Similar to
353      !!                time integration of sea surface height in
354      !!                OPA_SRC/DYN/sshwzv.F90.
355      !!----------------------------------------------------------------------
356      !
357      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
358      REAL(wp) :: z2dt
359      INTEGER :: jn
360
361!!----------------------------------------------------------------------
362      !
363      IF( neuler == 0 .AND. kt == nittrc000 ) THEN
364          z2dt = rdt                  ! set time step size (Euler)
365      ELSE
366          z2dt = 2._wp * rdt          ! set time step size (Leapfrog)
367      ENDIF
368
369      ! Forward Euler time step to compute "now"
370      DO jn=1,jp_fabm_surface+jp_fabm_bottom
371         fabm_st2Da(:,:,jn) = (fabm_st2db(:,:,jn) + z2dt * fabm_st2da(:,:,jn)) * tmask(:,:,1)
372      ENDDO
373
374      IF( neuler == 0 .AND. kt == nittrc000 ) THEN        ! Euler time-stepping at first time-step
375         !                                                ! (only swap)
376         fabm_st2Dn(:,:,:) = fabm_st2Da(:,:,:)
377         !
378      ELSE
379         ! Update now state + Asselin filter time stepping
380         fabm_st2Db(:,:,:) = (1._wp - 2._wp*atfp) * fabm_st2Dn(:,:,:) + &
381             atfp * ( fabm_st2Db(:,:,:) + fabm_st2Da(:,:,:) )
382         fabm_st2Dn(:,:,:) = fabm_st2Da(:,:,:)
383      ENDIF
384
385   END SUBROUTINE st2d_fabm_nxt
386
387   INTEGER FUNCTION trc_sms_fabm_alloc()
388      INTEGER :: jj,jk,jn
389      !!----------------------------------------------------------------------
390      !!              ***  ROUTINE trc_sms_fabm_alloc  ***
391      !!----------------------------------------------------------------------
392      !
393      ! ALLOCATE here the arrays specific to FABM
394      ALLOCATE( lk_rad_fabm(jp_fabm))
395      ALLOCATE( prn(jpi, jpj, jpk))
396      ALLOCATE( rho(jpi, jpj, jpk))
397      ALLOCATE( taubot(jpi, jpj))
398      ! ALLOCATE( tab(...) , STAT=trc_sms_fabm_alloc )
399
400      ! Allocate arrays to hold state for surface-attached and bottom-attached state variables
401      ALLOCATE(fabm_st2Dn(jpi, jpj, jp_fabm_surface+jp_fabm_bottom))
402      ALLOCATE(fabm_st2Da(jpi, jpj, jp_fabm_surface+jp_fabm_bottom))
403      ALLOCATE(fabm_st2Db(jpi, jpj, jp_fabm_surface+jp_fabm_bottom))
404
405      ! Work array to hold surface and bottom fluxes
406      ALLOCATE(flux(jpi,jp_fabm))
407
408      ! Work array to hold extinction coefficients
409      ALLOCATE(ext(jpi))
410      ext=0._wp
411
412      ! Allocate work arrays for vertical movement
413      ALLOCATE(w_ct(jpi,jpk,jp_fabm))
414      ALLOCATE(w_if(jpk,jp_fabm))
415      ALLOCATE(zwgt_if(jpk,jp_fabm))
416      ALLOCATE(flux_if(jpk,jp_fabm))
417      ALLOCATE(current_total(jpi,SIZE(model%conserved_quantities)))
418#if defined key_trdtrc && defined key_iomput
419      IF( lk_trdtrc ) ALLOCATE(tr_vmv(jpi,jpj,jpk,jp_fabm))
420      IF( lk_trdtrc ) ALLOCATE(tr_inp(jpi,jpj,jpk))
421#endif
422
423      trc_sms_fabm_alloc = 0      ! set to zero if no array to be allocated
424      !
425      IF( trc_sms_fabm_alloc /= 0 ) CALL ctl_warn('trc_sms_fabm_alloc : failed to allocate arrays')
426      !
427
428      ! Make FABM aware of diagnostics that are not needed [not included in output]
429      DO jn=1,size(model%diagnostic_variables)
430          !model%diagnostic_variables(jn)%save = iom_use(model%diagnostic_variables(jn)%name)
431      END DO
432      DO jn=1,size(model%horizontal_diagnostic_variables)
433          !model%horizontal_diagnostic_variables(jn)%save = iom_use(model%horizontal_diagnostic_variables(jn)%name)
434      END DO
435
436      ! Provide FABM with domain extents [after this, the save attribute of diagnostic variables can no longe change!]
437      call fabm_set_domain(model,jpi, jpj, jpk)
438
439      ! Provide FABM with the vertical indices of the surface and bottom, and the land-sea mask.
440      call model%set_bottom_index(mbkt)  ! NB mbkt extents should match dimension lengths provided to fabm_set_domain
441      call model%set_surface_index(1)
442      call fabm_set_mask(model,tmask,tmask(:,:,1)) ! NB tmask extents should match dimension lengths provided to fabm_set_domain
443
444      ! Send pointers to state data to FABM
445      do jn=1,jp_fabm
446         call fabm_link_bulk_state_data(model,jn,trn(:,:,:,jp_fabm_m1+jn))
447      end do
448      DO jn=1,jp_fabm_surface
449         CALL fabm_link_surface_state_data(model,jn,fabm_st2Dn(:,:,jn))
450      END DO
451      DO jn=1,jp_fabm_bottom
452         CALL fabm_link_bottom_state_data(model,jn,fabm_st2Dn(:,:,jp_fabm_surface+jn))
453      END DO
454
455      ! Send pointers to environmental data to FABM
456      call fabm_link_bulk_data(model,standard_variables%temperature,tsn(:,:,:,jp_tem))
457      call fabm_link_bulk_data(model,standard_variables%practical_salinity,tsn(:,:,:,jp_sal))
458      call fabm_link_bulk_data(model,standard_variables%density,rho(:,:,:))
459      call fabm_link_bulk_data(model,standard_variables%pressure,prn)
460      call fabm_link_horizontal_data(model,standard_variables%bottom_stress,taubot(:,:))
461      ! correct target for cell thickness depends on NEMO configuration:
462#ifdef key_vvl
463      call fabm_link_bulk_data(model,standard_variables%cell_thickness,e3t_n)
464#else
465      call fabm_link_bulk_data(model,standard_variables%cell_thickness,e3t_0)
466#endif
467      call fabm_link_horizontal_data(model,standard_variables%latitude,gphit)
468      call fabm_link_horizontal_data(model,standard_variables%longitude,glamt)
469      call fabm_link_scalar_data(model,standard_variables%number_of_days_since_start_of_the_year,daynumber_in_year)
470      call fabm_link_horizontal_data(model,standard_variables%wind_speed,wndm(:,:))
471      call fabm_link_horizontal_data(model,standard_variables%surface_downwelling_shortwave_flux,qsr(:,:))
472      call fabm_link_horizontal_data(model,standard_variables%bottom_depth_below_geoid,bathy(:,:))
473
474      swr_id = model%get_bulk_variable_id(standard_variables%downwelling_shortwave_flux)
475
476      ! Obtain user-specified input variables (read from NetCDF file)
477      call link_inputs
478      call update_inputs( nit000, .false. )
479
480      ! Check whether FABM has all required data
481      call fabm_check_ready(model)
482
483      ! Initialize state
484      DO jj=1,jpj
485         CALL fabm_initialize_surface_state(model,1,jpi,jj)
486         CALL fabm_initialize_bottom_state(model,1,jpi,jj)
487      END DO
488      DO jk=1,jpk
489         DO jj=1,jpj
490            CALL fabm_initialize_state(model,1,jpi,jj,jk)
491         END DO
492      END DO
493
494      ! Set mask for negativity corrections to the relevant states
495      DO jn=1,jp_fabm
496        IF (model%state_variables(jn)%minimum.ge.0) THEN
497          lk_rad_fabm(jn)=.TRUE.
498          IF(lwp) WRITE(numout,*) 'FABM clipping for '//TRIM(model%state_variables(jn)%name)//' activated.'
499        END IF
500      END DO
501
502      ! Mask land points in states with zeros, not nice, but coherent
503      ! with NEMO "convention":
504      DO jn=jp_fabm0,jp_fabm1
505        WHERE (tmask==0._wp)
506          trn(:,:,:,jn)=0._wp
507        END WHERE
508      END DO
509      DO jn=1,jp_fabm_surface+jp_fabm_bottom
510        WHERE (tmask(:,:,1)==0._wp)
511          fabm_st2Dn(:,:,jn)=0._wp
512        END WHERE
513      END DO
514
515      ! Copy initial condition for interface-attached state variables to "previous" state field
516      ! NB NEMO does this itself for pelagic state variables (trb) in TOP_SRC/trcini.F90.
517      fabm_st2Db = fabm_st2Dn
518
519      ! Initialise repair counters
520      repair_interior_count = 0
521      repair_surface_count = 0
522      repair_bottom_count = 0
523
524   END FUNCTION trc_sms_fabm_alloc
525
526#else
527   !!----------------------------------------------------------------------
528   !!   Dummy module                                        No FABM model
529   !!----------------------------------------------------------------------
530CONTAINS
531   SUBROUTINE trc_sms_fabm( kt )             ! Empty routine
532      INTEGER, INTENT( in ) ::   kt
533      WRITE(*,*) 'trc_sms_fabm: You should not have seen this print! error?', kt
534   END SUBROUTINE trc_sms_fabm
535#endif
536
537   !!======================================================================
538END MODULE trcsms_fabm
Note: See TracBrowser for help on using the repository browser.