Skip to content
Snippets Groups Projects
ini_budget.f90 204 KiB
Newer Older
!MNH_LIC Copyright 1995-2023 CNRS, Meteo-France and Universite Paul Sabatier
VIE Benoit's avatar
VIE Benoit committed
!MNH_LIC This is part of the Meso-NH software governed by the CeCILL-C licence
!MNH_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
!MNH_LIC for details. version 1.
!-----------------------------------------------------------------
! Modifications:
!  P. Wautelet 17/08/2020: add Budget_preallocate subroutine
!-----------------------------------------------------------------
module mode_ini_budget

  use mode_msg

  implicit none

  private

  public :: Budget_preallocate, Ini_budget

  integer, parameter :: NSOURCESMAX = 60 !Maximum number of sources in a budget

contains

subroutine Budget_preallocate()

use modd_budget, only: nbudgets, tbudgets,                                         &
                       NBUDGET_U, NBUDGET_V, NBUDGET_W, NBUDGET_TH, NBUDGET_TKE,   &
                       NBUDGET_RV, NBUDGET_RC, NBUDGET_RR, NBUDGET_RI, NBUDGET_RS, &
                       NBUDGET_RG, NBUDGET_RH, NBUDGET_SV1
use modd_nsv,    only: nsv, tsvlist
VIE Benoit's avatar
VIE Benoit committed

integer          :: ibudget
integer          :: jsv

call Print_msg( NVERB_DEBUG, 'BUD', 'Budget_preallocate', 'called' )

if ( allocated( tbudgets ) ) then
  call Print_msg( NVERB_WARNING, 'BUD', 'Budget_preallocate', 'tbudgets already allocated' )
  return
end if

nbudgets = NBUDGET_SV1 - 1 + nsv
allocate( tbudgets( nbudgets ) )

tbudgets(NBUDGET_U)%cname    = "UU"
tbudgets(NBUDGET_U)%ccomment = "Budget for U"
tbudgets(NBUDGET_U)%nid      = NBUDGET_U

tbudgets(NBUDGET_V)%cname    = "VV"
tbudgets(NBUDGET_V)%ccomment = "Budget for V"
tbudgets(NBUDGET_V)%nid      = NBUDGET_V

tbudgets(NBUDGET_W)%cname    = "WW"
tbudgets(NBUDGET_W)%ccomment = "Budget for W"
tbudgets(NBUDGET_W)%nid      = NBUDGET_W

tbudgets(NBUDGET_TH)%cname    = "TH"
tbudgets(NBUDGET_TH)%ccomment = "Budget for potential temperature"
tbudgets(NBUDGET_TH)%nid      = NBUDGET_TH

tbudgets(NBUDGET_TKE)%cname    = "TK"
tbudgets(NBUDGET_TKE)%ccomment = "Budget for turbulent kinetic energy"
tbudgets(NBUDGET_TKE)%nid      = NBUDGET_TKE

tbudgets(NBUDGET_RV)%cname    = "RV"
tbudgets(NBUDGET_RV)%ccomment = "Budget for water vapor mixing ratio"
tbudgets(NBUDGET_RV)%nid      = NBUDGET_RV

tbudgets(NBUDGET_RC)%cname    = "RC"
tbudgets(NBUDGET_RC)%ccomment = "Budget for cloud water mixing ratio"
tbudgets(NBUDGET_RC)%nid      = NBUDGET_RC

tbudgets(NBUDGET_RR)%cname    = "RR"
tbudgets(NBUDGET_RR)%ccomment = "Budget for rain water mixing ratio"
tbudgets(NBUDGET_RR)%nid      = NBUDGET_RR

tbudgets(NBUDGET_RI)%cname    = "RI"
tbudgets(NBUDGET_RI)%ccomment = "Budget for cloud ice mixing ratio"
tbudgets(NBUDGET_RI)%nid      = NBUDGET_RI

tbudgets(NBUDGET_RS)%cname    = "RS"
tbudgets(NBUDGET_RS)%ccomment = "Budget for snow/aggregate mixing ratio"
tbudgets(NBUDGET_RS)%nid      = NBUDGET_RS

tbudgets(NBUDGET_RG)%cname    = "RG"
tbudgets(NBUDGET_RG)%ccomment = "Budget for graupel mixing ratio"
tbudgets(NBUDGET_RG)%nid      = NBUDGET_RG

tbudgets(NBUDGET_RH)%cname    = "RH"
tbudgets(NBUDGET_RH)%ccomment = "Budget for hail mixing ratio"
tbudgets(NBUDGET_RH)%nid      = NBUDGET_RH

do jsv = 1, nsv
  ibudget = NBUDGET_SV1 - 1 + jsv
  tbudgets(ibudget)%cname    = Trim( tsvlist(jsv)%cmnhname )
  tbudgets(ibudget)%ccomment = 'Budget for scalar variable ' // Trim( tsvlist(jsv)%cmnhname )
VIE Benoit's avatar
VIE Benoit committed
  tbudgets(ibudget)%nid      = ibudget
end do


end subroutine Budget_preallocate


!     #################################################################
      SUBROUTINE Ini_budget(KLUOUT,PTSTEP,KSV,KRR,                    &
      ONUMDIFU,ONUMDIFTH,ONUMDIFSV,                                   &
      OHORELAX_UVWTH,OHORELAX_RV,OHORELAX_RC,OHORELAX_RR,             &
      OHORELAX_RI,OHORELAX_RS, OHORELAX_RG, OHORELAX_RH,OHORELAX_TKE, &
      OHORELAX_SV, OVE_RELAX, ove_relax_grd, OCHTRANS,                &
      ONUDGING,ODRAGTREE,ODEPOTREE, OAERO_EOL,                        &
      HRAD,HDCONV,HSCONV,HTURB,HTURBDIM,HCLOUD                        )
!     #################################################################
!
!!****  *INI_BUDGET* - routine to initialize the parameters for the budgets
!!
!!    PURPOSE
!!    -------
!       The purpose of this routine is to set or compute the parameters used
!     by the MESONH budgets. Names of files for budget recording are processed 
!     and storage arrays are initialized.               
!
!!**  METHOD
!!    ------
!!      The essential of information is passed by modules. The choice of budgets 
!!    and processes set by the user as integers is converted in "actions" 
!!    readable  by the subroutine BUDGET under the form of string characters. 
!!    For each complete process composed of several elementary processes, names 
!!    of elementary processes are concatenated in order to have an explicit name
!!    in the comment of the recording file for budget. 
!!
!!      
!!    EXTERNAL
!!    --------   
!!      None
!!
!!    IMPLICIT ARGUMENTS
!!    ------------------ 
!!      Modules MODD_*
!!
!!    REFERENCE
!!    ---------
!!      Book2 of documentation (routine INI_BUDGET)
!!      
!!
!!    AUTHOR
!!    ------
!!	P. Hereil      * Meteo France *
!!
!!    MODIFICATIONS
!!    -------------
!!      Original        01/03/95 
!!      J. Stein        25/06/95  put the sources in phase with the code
!!      J. Stein        20/07/95  reset to FALSE of all the switches when
!!                                CBUTYPE /= MASK or CART
!!      J. Stein        26/06/96  add the new sources + add the increment between
!!                                2 active processes
!!      J.-P. Pinty     13/12/96  Allowance of multiple SVs
!!      J.-P. Pinty     11/01/97  Includes deep convection ice and forcing processes
!!      J.-P. Lafore    10/02/98  Allocation of the RHODJs for budget
!!      V. Ducrocq      04/06/99  //  
!!      N. Asencio      18/06/99  // MASK case : delete KIMAX and KJMAX arguments,
!!                                GET_DIM_EXT_ll initializes the dimensions of the
!!                                extended local domain.
!!                                LBU_MASK and NBUSURF are allocated on the extended
!!                                local domain.
!!                                add 3 local variables IBUDIM1,IBUDIM2,IBUDIM3
!!                                to define the dimensions of the budget arrays
!!                                in the different cases CART and MASK
!!      J.-P. Pinty     23/09/00  add budget for C2R2
!!      V. Masson       18/11/02  add budget for 2way nesting
!!      O.Geoffroy      03/2006   Add KHKO scheme
!!      J.-P. Pinty     22/04/97  add the explicit hail processes
!!      C.Lac           10/08/07  Add ADV for PPM without contribution
!!                                of each direction
!!      C. Barthe       19/11/09  Add atmospheric electricity
!!      C.Lac           01/07/11  Add vegetation drag        
!!      P. Peyrille, M. Tomasini : include in the forcing term the 2D forcing
!!                                terms in term 2DFRC search for modif PP . but Not very clean! 
!!      C .Lac          27/05/14    add negativity corrections for chemical species
!!      C.Lac           29/01/15  Correction for NSV_USER
!!      J.Escobar       02/10/2015 modif for JPHEXT(JPVEXT) variable  
!!      C.Lac           04/12/15  Correction for LSUPSAT 
!  C. Lac         04/2016: negative contribution to the budget split between advection, turbulence and microphysics for KHKO/C2R2
!  C. Barthe      01/2016: add budget for LIMA
!  C. Lac         10/2016: add budget for droplet deposition
!  S. Riette      11/2016: new budgets for ICE3/ICE4
!  P. Wautelet 05/2016-04/2018: new data structures and calls for I/O
!  P. Wautelet 10/04/2019: replace ABORT and STOP calls by Print_msg
!  P. Wautelet 15/11/2019: remove unused CBURECORD variable
!  P. Wautelet 24/02/2020: bugfix: corrected condition for budget NCDEPITH
!  P. Wautelet 26/02/2020: bugfix: rename CEVA->REVA for budget for raindrop evaporation in C2R2 (necessary after commit 4ed805fc)
!  P. Wautelet 26/02/2020: bugfix: add missing condition on OCOLD for NSEDIRH budget in LIMA case
!  P. Wautelet 02-03/2020: use the new data structures and subroutines for budgets
!  B. Vie      02/03/2020: LIMA negativity checks after turbulence, advection and microphysics budgets
!  P .Wautelet 09/03/2020: add missing budgets for electricity
!  P. Wautelet 25/03/2020: add missing ove_relax_grd
!  P. Wautelet 23/04/2020: add nid in tbudgetdata datatype
!  P. Wautelet + Benoit Vié 11/06/2020: improve removal of negative scalar variables + adapt the corresponding budgets
!  P. Wautelet 30/06/2020: use NADVSV when possible
!  P. Wautelet 30/06/2020: add NNETURSV, NNEADVSV and NNECONSV variables
!  P. Wautelet 06/07/2020: bugfix: add condition on HTURB for NETUR sources for SV budgets
!  P. Wautelet 08/12/2020: add nbusubwrite and nbutotwrite
!  P. Wautelet 11/01/2021: ignore xbuwri for cartesian boxes (write at every xbulen interval)
!  P. Wautelet 01/02/2021: bugfix: add missing CEDS source terms for SV budgets
!  P. Wautelet 02/02/2021: budgets: add missing source terms for SV budgets in LIMA
!  P. Wautelet 03/02/2021: budgets: add new source if LIMA splitting: CORR2
!  P. Wautelet 10/02/2021: budgets: add missing sources for NSV_C2R2BEG+3 budget
!  P. Wautelet 11/02/2021: budgets: add missing term SCAV for NSV_LIMA_SCAVMASS budget
!  P. Wautelet 02/03/2021: budgets: add terms for blowing snow
!  P. Wautelet 04/03/2021: budgets: add terms for drag due to buildings
!  P. Wautelet 17/03/2021: choose source terms for budgets with character strings instead of multiple integer variables
!  C. Barthe   14/03/2022: budgets: add terms for CIBU and RDSF in LIMA
!  M. Taufour  01/07/2022: budgets: add concentration for snow, graupel, hail
!-------------------------------------------------------------------------------
!
!*       0.    DECLARATIONS
!              ------------ 
!
use modd_2d_frc,        only: l2d_adv_frc, l2d_rel_frc
use modd_blowsnow,      only: lblowsnow
use modd_blowsnow_n,    only: lsnowsubl
use modd_budget
use modd_ch_aerosol,    only: lorilam
use modd_conf,          only: l1d, lcartesian, lforcing, lthinshell, nmodel
use modd_dim_n,         only: nimax_ll, njmax_ll, nkmax
use modd_dragbldg_n,    only: ldragbldg
use modd_dust,          only: ldust
use modd_dyn,           only: lcorio, xseglen
use modd_dyn_n,         only: xtstep, locean
use modd_elec_descr,    only: linductive, lrelax2fw_ion
use modd_field,         only: TYPEREAL
use modd_fire_n,        only: lblaze
use modd_nsv,           only: nsv_aerbeg, nsv_aerend, nsv_aerdepbeg, nsv_aerdepend, nsv_c2r2beg, nsv_c2r2end,      &
VIE Benoit's avatar
VIE Benoit committed
                              nsv_chembeg, nsv_chemend, nsv_chicbeg, nsv_chicend, nsv_csbeg, nsv_csend,            &
                              nsv_dstbeg, nsv_dstend, nsv_dstdepbeg, nsv_dstdepend, nsv_elecbeg, nsv_elecend,      &
#ifdef MNH_FOREFIRE
                              nsv_ffbeg, nsv_ffend,                                                                &
#endif
                              nsv_lgbeg, nsv_lgend,                                                                &
                              nsv_lima_beg, nsv_lima_end, nsv_lima_ccn_acti, nsv_lima_ccn_free, nsv_lima_hom_haze, &
                              nsv_lima_ifn_free, nsv_lima_ifn_nucl, nsv_lima_imm_nucl,                             &
                              nsv_lima_nc, nsv_lima_nr, nsv_lima_ni, nsv_lima_ns, nsv_lima_ng, nsv_lima_nh,        &
                              nsv_lima_scavmass, nsv_lima_spro,                                                    &
VIE Benoit's avatar
VIE Benoit committed
                              nsv_lnoxbeg, nsv_lnoxend, nsv_ppbeg, nsv_ppend,                                      &
                              nsv_sltbeg, nsv_sltend, nsv_sltdepbeg, nsv_sltdepend, nsv_snwbeg, nsv_snwend,        &
VIE Benoit's avatar
VIE Benoit committed
use modd_parameters,   only: jphext
use modd_param_c2r2,   only: ldepoc_c2r2 => ldepoc, lrain_c2r2 => lrain, lsedc_c2r2 => lsedc, lsupsat_c2r2 => lsupsat
use modd_param_ice,    only: ladj_after, ladj_before, ldeposc_ice => ldeposc, lred, lsedic_ice => lsedic, lwarm_ice => lwarm
use modd_param_n,      only: cactccn, celec
use modd_param_lima,   only: laero_mass_lima => laero_mass, lacti_lima => lacti, ldepoc_lima => ldepoc, &
                             lhhoni_lima => lhhoni, lmeyers_lima => lmeyers, lnucl_lima => lnucl,       &
                             lptsplit,                                                                                       &
                             lscav_lima => lscav, lsedc_lima => lsedc, lsedi_lima => lsedi,             &
                             lspro_lima => lspro,  lcibu, lrdsf,           &
                             nmom_c, nmom_r, nmom_i, nmom_s, nmom_g, nmom_h, nmod_ccn, nmod_ifn, nmod_imm
use modd_ref,          only: lcouples
use modd_salt,         only: lsalt
use modd_turb_n,       only: lsubg_cond
use modd_viscosity,    only: lvisc, lvisc_r, lvisc_sv, lvisc_th, lvisc_uvw

USE MODE_ll

IMPLICIT NONE
!
!*       0.1   declarations of argument
!
!
INTEGER,         INTENT(IN) :: KLUOUT   ! Logical unit number for prints
REAL, INTENT(IN)    :: PTSTEP           ! time step
INTEGER, INTENT(IN) :: KSV              ! number of scalar variables
INTEGER, INTENT(IN) :: KRR              ! number of moist variables
LOGICAL, INTENT(IN) :: ONUMDIFU         ! switch to activate the numerical
                                        ! diffusion for momentum
LOGICAL, INTENT(IN) :: ONUMDIFTH        ! for meteorological scalar variables
LOGICAL, INTENT(IN) :: ONUMDIFSV        ! for tracer scalar variables
LOGICAL, INTENT(IN) :: OHORELAX_UVWTH  ! switch for the
                       ! horizontal relaxation for U,V,W,TH
LOGICAL, INTENT(IN) :: OHORELAX_RV     ! switch for the
                       ! horizontal relaxation for Rv
LOGICAL, INTENT(IN) :: OHORELAX_RC     ! switch for the
                       ! horizontal relaxation for Rc
LOGICAL, INTENT(IN) :: OHORELAX_RR     ! switch for the
                       ! horizontal relaxation for Rr
LOGICAL, INTENT(IN) :: OHORELAX_RI     ! switch for the
                       ! horizontal relaxation for Ri
LOGICAL, INTENT(IN) :: OHORELAX_RS     ! switch for the
                       ! horizontal relaxation for Rs
LOGICAL, INTENT(IN) :: OHORELAX_RG     ! switch for the
                       ! horizontal relaxation for Rg
LOGICAL, INTENT(IN) :: OHORELAX_RH     ! switch for the
                       ! horizontal relaxation for Rh
LOGICAL, INTENT(IN) :: OHORELAX_TKE    ! switch for the
                       ! horizontal relaxation for tke
LOGICAL,DIMENSION(:),INTENT(IN):: OHORELAX_SV     ! switch for the
                       ! horizontal relaxation for scalar variables
LOGICAL, INTENT(IN) :: OVE_RELAX        ! switch to activate the vertical
                                        ! relaxation
logical, intent(in) :: ove_relax_grd    ! switch to activate the vertical
                                        ! relaxation to the lowest verticals
LOGICAL, INTENT(IN) :: OCHTRANS         ! switch to activate convective
                                        !transport for SV
LOGICAL, INTENT(IN) :: ONUDGING         ! switch to activate nudging
LOGICAL, INTENT(IN) :: ODRAGTREE        ! switch to activate vegetation drag
LOGICAL, INTENT(IN) :: ODEPOTREE        ! switch to activate droplet deposition on tree
LOGICAL, INTENT(IN) :: OAERO_EOL        ! switch to activate wind turbine wake
CHARACTER (LEN=*), INTENT(IN) :: HRAD   ! type of the radiation scheme
CHARACTER (LEN=*), INTENT(IN) :: HDCONV ! type of the deep convection scheme
CHARACTER (LEN=*), INTENT(IN) :: HSCONV ! type of the shallow convection scheme
CHARACTER (LEN=*), INTENT(IN) :: HTURB  ! type of the turbulence scheme
CHARACTER (LEN=*), INTENT(IN) :: HTURBDIM! dimensionnality of the turbulence 
                                        ! scheme
CHARACTER (LEN=*), INTENT(IN) :: HCLOUD ! type of microphysical scheme
!
!*       0.2   declarations of local variables
!
real, parameter :: ITOL = 1e-6

INTEGER :: JI, JJ                                         ! loop indices
INTEGER :: IIMAX_ll, IJMAX_ll ! size of the physical global domain
INTEGER :: IIU, IJU                                       ! size along x and y directions
                                                          ! of the extended subdomain
INTEGER :: IBUDIM1                                        ! first dimension of the budget arrays
                                                          ! = NBUIMAX in CART case
                                                          ! = NBUKMAX in MASK case
INTEGER :: IBUDIM2                                        ! second dimension of the budget arrays
                                                          ! = NBUJMAX in CART case
                                                          ! = nbusubwrite in MASK case
INTEGER :: IBUDIM3                                        ! third dimension of the budget arrays
                                                          ! = NBUKMAX in CART case
                                                          ! = NBUMASK in MASK case
INTEGER             :: JSV               ! loop indice for the SVs
INTEGER             :: IINFO_ll ! return status of the interface routine
integer             :: ibudget
logical             :: gtmp
type(tbusourcedata) :: tzsource ! Used to prepare metadate of source terms

call Print_msg( NVERB_DEBUG, 'BUD', 'Ini_budget', 'called' )
!
!*       1.    COMPUTE BUDGET VARIABLES
!              ------------------------
!
NBUSTEP = NINT (XBULEN / PTSTEP)
NBUTSHIFT=0
!
!  common dimension for all CBUTYPE values
!
IF (LBU_KCP) THEN
  NBUKMAX = 1
ELSE
  NBUKMAX = NBUKH - NBUKL +1
END IF
!
if ( cbutype == 'CART' .or. cbutype == 'MASK' ) then
  !Check if xbulen is a multiple of xtstep (within tolerance)
  if ( Abs( Nint( xbulen / xtstep ) * xtstep - xbulen ) > ( ITOL * xtstep ) ) &
    call Print_msg( NVERB_WARNING, 'BUD', 'Ini_budget', 'xbulen is not a multiple of xtstep' )

  if ( cbutype == 'CART' ) then
    !Check if xseglen is a multiple of xbulen (within tolerance)
    if ( Abs( Nint( xseglen / xbulen ) * xbulen - xseglen ) > ( ITOL * xseglen ) ) &
      call Print_msg( NVERB_WARNING, 'BUD', 'Ini_budget', 'xseglen is not a multiple of xbulen' )

    !Write cartesian budgets every xbulen time period (do not take xbuwri into account)
    xbuwri = xbulen

    nbusubwrite = 1                                      !Number of budget time average periods for each write
    nbutotwrite = nbusubwrite * Nint( xseglen / xbulen ) !Total number of budget time average periods
  else if ( cbutype == 'MASK' ) then
    !Check if xbuwri is a multiple of xtstep (within tolerance)
    if ( Abs( Nint( xbuwri / xtstep ) * xtstep - xbuwri ) > ( ITOL * xtstep ) ) &
      call Print_msg( NVERB_WARNING, 'BUD', 'Ini_budget', 'xbuwri is not a multiple of xtstep' )

    !Check if xbuwri is a multiple of xbulen (within tolerance)
    if ( Abs( Nint( xbuwri / xbulen ) * xbulen - xbuwri ) > ( ITOL * xbulen ) ) &
      call Print_msg( NVERB_WARNING, 'BUD', 'Ini_budget', 'xbuwri is not a multiple of xbulen' )

    !Check if xseglen is a multiple of xbuwri (within tolerance)
    if ( Abs( Nint( xseglen / xbuwri ) * xbuwri - xseglen ) > ( ITOL * xseglen ) ) &
      call Print_msg( NVERB_WARNING, 'BUD', 'Ini_budget', 'xseglen is not a multiple of xbuwri' )

    nbusubwrite = Nint ( xbuwri / xbulen )               !Number of budget time average periods for each write
    nbutotwrite = nbusubwrite * Nint( xseglen / xbuwri ) !Total number of budget time average periods
  end if
end if

IF (CBUTYPE=='CART') THEN              ! cartesian case only
!
  IF ( NBUIL < 1 )        CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUIL too small (<1)' )
  IF ( NBUIL > NIMAX_ll ) CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUIL too large (>NIMAX)' )
  IF ( NBUIH < 1 )        CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUIH too small (<1)' )
  IF ( NBUIH > NIMAX_ll ) CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUIH too large (>NIMAX)' )
  IF ( NBUIH < NBUIL )    CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUIH < NBUIL' )
  IF (LBU_ICP) THEN
    NBUIMAX_ll = 1
  ELSE
    NBUIMAX_ll = NBUIH - NBUIL +1
  END IF

  IF ( NBUJL < 1 )        CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUJL too small (<1)' )
  IF ( NBUJL > NJMAX_ll ) CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUJL too large (>NJMAX)' )
  IF ( NBUJH < 1 )        CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUJH too small (<1)' )
  IF ( NBUJH > NJMAX_ll ) CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUJH too large (>NJMAX)' )
  IF ( NBUJH < NBUJL )    CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUJH < NBUJL' )
  IF (LBU_JCP) THEN
    NBUJMAX_ll = 1
  ELSE
    NBUJMAX_ll = NBUJH - NBUJL +1
  END IF

  IF ( NBUKL < 1 )     CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUKL too small (<1)' )
  IF ( NBUKL > NKMAX ) CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUKL too large (>NKMAX)' )
  IF ( NBUKH < 1 )     CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUKH too small (<1)' )
  IF ( NBUKH > NKMAX ) CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUKH too large (>NKMAX)' )
  IF ( NBUKH < NBUKL ) CALL Print_msg( NVERB_ERROR, 'BUD', 'Ini_budget', 'NBUKH < NBUKL' )

  CALL GET_INTERSECTION_ll(NBUIL+JPHEXT,NBUJL+JPHEXT,NBUIH+JPHEXT,NBUJH+JPHEXT, &
      NBUSIL,NBUSJL,NBUSIH,NBUSJH,"PHYS",IINFO_ll)
  IF ( IINFO_ll /= 1 ) THEN ! 
    IF (LBU_ICP) THEN 
      NBUIMAX = 1
    ELSE
      NBUIMAX = NBUSIH - NBUSIL +1
    END IF
    IF (LBU_JCP) THEN 
      NBUJMAX = 1
    ELSE
      NBUJMAX =  NBUSJH - NBUSJL +1
    END IF
  ELSE ! the intersection is void 
    CBUTYPE='SKIP'  ! no budget on this processor       
    NBUIMAX = 0     ! in order to allocate void arrays
    NBUJMAX = 0
  ENDIF
! three first dimensions of budget arrays in cart and skip cases
   IBUDIM1=NBUIMAX
   IBUDIM2=NBUJMAX
   IBUDIM3=NBUKMAX
! these variables are not be used 
   NBUMASK=-1
!
ELSEIF (CBUTYPE=='MASK') THEN          ! mask case only 
!
  LBU_ENABLE=.TRUE.
                                    ! result on the FM_FILE
  NBUTIME = 1

  CALL GET_DIM_EXT_ll ('B', IIU,IJU)
  ALLOCATE( LBU_MASK( IIU ,IJU, NBUMASK) )
  LBU_MASK(:,:,:)=.FALSE.
  ALLOCATE( NBUSURF( IIU, IJU, NBUMASK, nbusubwrite) )
  NBUSURF(:,:,:,:) = 0
!
! three first dimensions of budget arrays in mask case
!  the order of the dimensions are the order expected in WRITE_DIACHRO routine:
!  x,y,z,time,mask,processus  and in this case x and y are missing
!  first dimension of the arrays : dimension along K
!  second dimension of the arrays : number of the budget time period
!  third dimension of the arrays : number of the budget masks zones
  IBUDIM1=NBUKMAX
  IBUDIM2=nbusubwrite
  IBUDIM3=NBUMASK
! these variables are not used in this case
  NBUIMAX=-1
  NBUJMAX=-1
! the beginning and the end along x and y direction : global extended domain
 ! get dimensions of the physical global domain
   CALL GET_GLOBALDIMS_ll (IIMAX_ll,IJMAX_ll)
   NBUIL=1
   NBUIH=IIMAX_ll + 2 * JPHEXT
   NBUJL=1 
   NBUJH=IJMAX_ll + 2 * JPHEXT
!
ELSE                      ! default case
!
  LBU_ENABLE=.FALSE.
  NBUIMAX = -1
  NBUJMAX = -1
  LBU_RU = .FALSE.
  LBU_RV = .FALSE.
  LBU_RW = .FALSE.
  LBU_RTH= .FALSE.
  LBU_RTKE= .FALSE.
  LBU_RRV= .FALSE.
  LBU_RRC= .FALSE.
  LBU_RRR= .FALSE.
  LBU_RRI= .FALSE.
  LBU_RRS= .FALSE.
  LBU_RRG= .FALSE.
  LBU_RRH= .FALSE.
  LBU_RSV= .FALSE.
!
! three first dimensions of budget arrays in default case
  IBUDIM1=0
  IBUDIM2=0
  IBUDIM3=0
!
END IF  
!
!
!-------------------------------------------------------------------------------
!
!*       2.    ALLOCATE MEMORY FOR BUDGET ARRAYS AND INITIALIZE
!              ------------------------------------------------
!
LBU_BEG =.TRUE. 
!
!-------------------------------------------------------------------------------
!
!*       3.    INITALIZE VARIABLES
!              -------------------
!
!Create intermediate variable to store rhodj for scalar variables
if ( lbu_rth .or. lbu_rtke .or. lbu_rrv .or. lbu_rrc .or. lbu_rrr .or. &
     lbu_rri .or. lbu_rrs  .or. lbu_rrg .or. lbu_rrh .or. lbu_rsv      ) then
  allocate( tburhodj )

  tburhodj%cmnhname  = 'RhodJS'
  tburhodj%cstdname  = ''
  tburhodj%clongname = 'RhodJS'
  tburhodj%cunits    = 'kg'
  tburhodj%ccomment  = 'RhodJ for Scalars variables'
  tburhodj%ngrid     = 1
  tburhodj%ntype     = TYPEREAL
  tburhodj%ndims     = 3

  allocate( tburhodj%xdata(ibudim1, ibudim2, ibudim3) )
  tburhodj%xdata(:, :, :) = 0.
end if

VIE Benoit's avatar
VIE Benoit committed
531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000
tzsource%ntype    = TYPEREAL
tzsource%ndims    = 3

! Budget of RU
tbudgets(NBUDGET_U)%lenabled = lbu_ru

if ( lbu_ru ) then
  allocate( tbudgets(NBUDGET_U)%trhodj )

  tbudgets(NBUDGET_U)%trhodj%cmnhname  = 'RhodJX'
  tbudgets(NBUDGET_U)%trhodj%cstdname  = ''
  tbudgets(NBUDGET_U)%trhodj%clongname = 'RhodJX'
  tbudgets(NBUDGET_U)%trhodj%cunits    = 'kg'
  tbudgets(NBUDGET_U)%trhodj%ccomment  = 'RhodJ for momentum along X axis'
  tbudgets(NBUDGET_U)%trhodj%ngrid     = 2
  tbudgets(NBUDGET_U)%trhodj%ntype     = TYPEREAL
  tbudgets(NBUDGET_U)%trhodj%ndims     = 3

  allocate( tbudgets(NBUDGET_U)%trhodj%xdata(ibudim1, ibudim2, ibudim3) )
  tbudgets(NBUDGET_U)%trhodj%xdata(:, :, :) = 0.

  !Allocate all basic source terms (used or not)
  !The size should be large enough (bigger than necessary is OK)
  tbudgets(NBUDGET_U)%nsourcesmax = NSOURCESMAX
  allocate( tbudgets(NBUDGET_U)%tsources(NSOURCESMAX) )

  allocate( tbudgets(NBUDGET_U)%xtmpstore(ibudim1, ibudim2, ibudim3) )

  tbudgets(NBUDGET_U)%tsources(:)%ngroup = 0

  tzsource%ccomment = 'Budget of momentum along X axis'
  tzsource%ngrid    = 2

  tzsource%cunits   = 'm s-1'

  tzsource%cmnhname   = 'INIF'
  tzsource%clongname  = 'initial state'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, odonotinit = .true., ooverwrite = .true. )

  tzsource%cmnhname   = 'ENDF'
  tzsource%clongname  = 'final state'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, odonotinit = .true., ooverwrite = .true. )

  tzsource%cmnhname   = 'AVEF'
  tzsource%clongname  = 'averaged state'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, odonotinit = .true., ooverwrite = .false. )

  tzsource%cunits   = 'm s-2'

  tzsource%cmnhname   = 'ASSE'
  tzsource%clongname  = 'time filter (Asselin)'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'NEST'
  tzsource%clongname  = 'nesting'
  tzsource%lavailable = nmodel > 1
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'FRC'
  tzsource%clongname  = 'forcing'
  tzsource%lavailable = lforcing
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'NUD'
  tzsource%clongname  = 'nudging'
  tzsource%lavailable = onudging
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'CURV'
  tzsource%clongname  = 'curvature'
  tzsource%lavailable = .not.l1d .and. .not.lcartesian
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'COR'
  tzsource%clongname  = 'Coriolis'
  tzsource%lavailable = lcorio
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'DIF'
  tzsource%clongname  = 'numerical diffusion'
  tzsource%lavailable = onumdifu
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'REL'
  tzsource%clongname  = 'relaxation'
  tzsource%lavailable = ohorelax_uvwth .or. ove_relax .or. ove_relax_grd
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'DRAG'
  tzsource%clongname  = 'drag force due to trees'
  tzsource%lavailable = odragtree
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'DRAGEOL'
  tzsource%clongname  = 'drag force due to wind turbine'
  tzsource%lavailable = OAERO_EOL
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'DRAGB'
  tzsource%clongname  = 'drag force due to buildings'
  tzsource%lavailable = ldragbldg
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'VTURB'
  tzsource%clongname  = 'vertical turbulent diffusion'
  tzsource%lavailable = hturb == 'TKEL'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'HTURB'
  tzsource%clongname  = 'horizontal turbulent diffusion'
  tzsource%lavailable = hturb == 'TKEL' .and. HTURBDIM == '3DIM'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'MAFL'
  tzsource%clongname  = 'mass flux'
  tzsource%lavailable = hsconv == 'EDKF'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'VISC'
  tzsource%clongname  = 'viscosity'
  tzsource%lavailable = lvisc .and. lvisc_uvw
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'ADV'
  tzsource%clongname  = 'advection'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )

  tzsource%cmnhname   = 'PRES'
  tzsource%clongname  = 'pressure'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource )


  call Sourcelist_sort_compact( tbudgets(NBUDGET_U) )

  call Sourcelist_scan( tbudgets(NBUDGET_U), cbulist_ru )
end if

! Budget of RV
tbudgets(NBUDGET_V)%lenabled = lbu_rv

if ( lbu_rv ) then
  allocate( tbudgets(NBUDGET_V)%trhodj )

  tbudgets(NBUDGET_V)%trhodj%cmnhname  = 'RhodJY'
  tbudgets(NBUDGET_V)%trhodj%cstdname  = ''
  tbudgets(NBUDGET_V)%trhodj%clongname = 'RhodJY'
  tbudgets(NBUDGET_V)%trhodj%cunits    = 'kg'
  tbudgets(NBUDGET_V)%trhodj%ccomment  = 'RhodJ for momentum along Y axis'
  tbudgets(NBUDGET_V)%trhodj%ngrid     = 3
  tbudgets(NBUDGET_V)%trhodj%ntype     = TYPEREAL
  tbudgets(NBUDGET_V)%trhodj%ndims     = 3

  allocate( tbudgets(NBUDGET_V)%trhodj%xdata(ibudim1, ibudim2, ibudim3) )
  tbudgets(NBUDGET_V)%trhodj%xdata(:, :, :) = 0.

  !Allocate all basic source terms (used or not)
  !The size should be large enough (bigger than necessary is OK)
  tbudgets(NBUDGET_V)%nsourcesmax = NSOURCESMAX
  allocate( tbudgets(NBUDGET_V)%tsources(NSOURCESMAX) )

  allocate( tbudgets(NBUDGET_V)%xtmpstore(ibudim1, ibudim2, ibudim3) )

  tbudgets(NBUDGET_V)%tsources(:)%ngroup = 0

  tzsource%ccomment = 'Budget of momentum along Y axis'
  tzsource%ngrid    = 3

  tzsource%cunits   = 'm s-1'

  tzsource%cmnhname   = 'INIF'
  tzsource%clongname  = 'initial state'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, odonotinit = .true., ooverwrite = .true. )

  tzsource%cmnhname   = 'ENDF'
  tzsource%clongname  = 'final state'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, odonotinit = .true., ooverwrite = .true. )

  tzsource%cmnhname   = 'AVEF'
  tzsource%clongname  = 'averaged state'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, odonotinit = .true., ooverwrite = .false. )

  tzsource%cunits   = 'm s-2'

  tzsource%cmnhname   = 'ASSE'
  tzsource%clongname  = 'time filter (Asselin)'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'NEST'
  tzsource%clongname  = 'nesting'
  tzsource%lavailable = nmodel > 1
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'FRC'
  tzsource%clongname  = 'forcing'
  tzsource%lavailable = lforcing
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'NUD'
  tzsource%clongname  = 'nudging'
  tzsource%lavailable = onudging
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'CURV'
  tzsource%clongname  = 'curvature'
  tzsource%lavailable = .not.l1d .and. .not.lcartesian
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'COR'
  tzsource%clongname  = 'Coriolis'
  tzsource%lavailable = lcorio
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'DIF'
  tzsource%clongname  = 'numerical diffusion'
  tzsource%lavailable = onumdifu
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'REL'
  tzsource%clongname  = 'relaxation'
  tzsource%lavailable = ohorelax_uvwth .or. ove_relax .or. ove_relax_grd
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'DRAG'
  tzsource%clongname  = 'drag force due to trees'
  tzsource%lavailable = odragtree
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'DRAGEOL'
  tzsource%clongname  = 'drag force due to wind turbine'
  tzsource%lavailable = OAERO_EOL
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'DRAGB'
  tzsource%clongname  = 'drag force due to buildings'
  tzsource%lavailable = ldragbldg
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'VTURB'
  tzsource%clongname  = 'vertical turbulent diffusion'
  tzsource%lavailable = hturb == 'TKEL'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'HTURB'
  tzsource%clongname  = 'horizontal turbulent diffusion'
  tzsource%lavailable = hturb == 'TKEL' .and. HTURBDIM == '3DIM'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'MAFL'
  tzsource%clongname  = 'mass flux'
  tzsource%lavailable = hsconv == 'EDKF'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'VISC'
  tzsource%clongname  = 'viscosity'
  tzsource%lavailable = lvisc .and. lvisc_uvw
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'ADV'
  tzsource%clongname  = 'advection'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )

  tzsource%cmnhname   = 'PRES'
  tzsource%clongname  = 'pressure'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource )


  call Sourcelist_sort_compact( tbudgets(NBUDGET_V) )

  call Sourcelist_scan( tbudgets(NBUDGET_V), cbulist_rv )
end if

! Budget of RW
tbudgets(NBUDGET_W)%lenabled = lbu_rw

if ( lbu_rw ) then
  allocate( tbudgets(NBUDGET_W)%trhodj )

  tbudgets(NBUDGET_W)%trhodj%cmnhname  = 'RhodJZ'
  tbudgets(NBUDGET_W)%trhodj%cstdname  = ''
  tbudgets(NBUDGET_W)%trhodj%clongname = 'RhodJZ'
  tbudgets(NBUDGET_W)%trhodj%cunits    = 'kg'
  tbudgets(NBUDGET_W)%trhodj%ccomment  = 'RhodJ for momentum along Z axis'
  tbudgets(NBUDGET_W)%trhodj%ngrid     = 4
  tbudgets(NBUDGET_W)%trhodj%ntype     = TYPEREAL
  tbudgets(NBUDGET_W)%trhodj%ndims     = 3

  allocate( tbudgets(NBUDGET_W)%trhodj%xdata(ibudim1, ibudim2, ibudim3) )
  tbudgets(NBUDGET_W)%trhodj%xdata(:, :, :) = 0.

  !Allocate all basic source terms (used or not)
  !The size should be large enough (bigger than necessary is OK)
  tbudgets(NBUDGET_W)%nsourcesmax = NSOURCESMAX
  allocate( tbudgets(NBUDGET_W)%tsources(NSOURCESMAX) )

  allocate( tbudgets(NBUDGET_W)%xtmpstore(ibudim1, ibudim2, ibudim3) )

  tbudgets(NBUDGET_W)%tsources(:)%ngroup = 0

  tzsource%ccomment = 'Budget of momentum along Z axis'
  tzsource%ngrid    = 4

  tzsource%cunits   = 'm s-1'

  tzsource%cmnhname   = 'INIF'
  tzsource%clongname  = 'initial state'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, odonotinit = .true., ooverwrite = .true. )

  tzsource%cmnhname   = 'ENDF'
  tzsource%clongname  = 'final state'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, odonotinit = .true., ooverwrite = .true. )

  tzsource%cmnhname   = 'AVEF'
  tzsource%clongname  = 'averaged state'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, odonotinit = .true., ooverwrite = .false. )

  tzsource%cunits   = 'm s-2'

  tzsource%cmnhname   = 'ASSE'
  tzsource%clongname  = 'time filter (Asselin)'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'NEST'
  tzsource%clongname  = 'nesting'
  tzsource%lavailable = nmodel > 1
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'FRC'
  tzsource%clongname  = 'forcing'
  tzsource%lavailable = lforcing
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'NUD'
  tzsource%clongname  = 'nudging'
  tzsource%lavailable = onudging
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'CURV'
  tzsource%clongname  = 'curvature'
  tzsource%lavailable = .not.l1d .and. .not.lcartesian .and. .not.lthinshell
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'COR'
  tzsource%clongname  = 'Coriolis'
  tzsource%lavailable = lcorio .and. .not.l1d .and. .not.lthinshell
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'DIF'
  tzsource%clongname  = 'numerical diffusion'
  tzsource%lavailable = onumdifu
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'REL'
  tzsource%clongname  = 'relaxation'
  tzsource%lavailable = ohorelax_uvwth .or. ove_relax .or. ove_relax_grd
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'VTURB'
  tzsource%clongname  = 'vertical turbulent diffusion'
  tzsource%lavailable = hturb == 'TKEL'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'HTURB'
  tzsource%clongname  = 'horizontal turbulent diffusion'
  tzsource%lavailable = hturb == 'TKEL' .and. HTURBDIM == '3DIM'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'VISC'
  tzsource%clongname  = 'viscosity'
  tzsource%lavailable = lvisc .and. lvisc_uvw
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'GRAV'
  tzsource%clongname  = 'gravity'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'ADV'
  tzsource%clongname  = 'advection'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'PRES'
  tzsource%clongname  = 'pressure'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  tzsource%cmnhname   = 'DRAGEOL'
  tzsource%clongname  = 'drag force due to wind turbine'
  tzsource%lavailable = OAERO_EOL
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource )

  call Sourcelist_sort_compact( tbudgets(NBUDGET_W) )

  call Sourcelist_scan( tbudgets(NBUDGET_W), cbulist_rw )
end if

! Budget of RTH
tbudgets(NBUDGET_TH)%lenabled = lbu_rth

if ( lbu_rth ) then
  tbudgets(NBUDGET_TH)%trhodj => tburhodj

  !Allocate all basic source terms (used or not)
  !The size should be large enough (bigger than necessary is OK)
  tbudgets(NBUDGET_TH)%nsourcesmax = NSOURCESMAX
  allocate( tbudgets(NBUDGET_TH)%tsources(NSOURCESMAX) )

  allocate( tbudgets(NBUDGET_TH)%xtmpstore(ibudim1, ibudim2, ibudim3) )

  tbudgets(NBUDGET_TH)%tsources(:)%ngroup = 0

  tzsource%ccomment = 'Budget of potential temperature'
  tzsource%ngrid    = 1

  tzsource%cunits   = 'K'

  tzsource%cmnhname   = 'INIF'
  tzsource%clongname  = 'initial state'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, odonotinit = .true., ooverwrite = .true. )

  tzsource%cmnhname   = 'ENDF'
  tzsource%clongname  = 'final state'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, odonotinit = .true., ooverwrite = .true. )

  tzsource%cmnhname   = 'AVEF'
  tzsource%clongname  = 'averaged state'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, odonotinit = .true., ooverwrite = .false. )

  tzsource%cunits   = 'K s-1'

  tzsource%cmnhname   = 'ASSE'
  tzsource%clongname  = 'time filter (Asselin)'
  tzsource%lavailable = .true.
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource )

  tzsource%cmnhname   = 'NEST'
  tzsource%clongname  = 'nesting'
  tzsource%lavailable = nmodel > 1
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource )

  tzsource%cmnhname   = 'FRC'
  tzsource%clongname  = 'forcing'
  tzsource%lavailable = lforcing
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource )

  tzsource%cmnhname   = '2DADV'
  tzsource%clongname  = 'advective forcing'
  tzsource%lavailable = l2d_adv_frc
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource )

  tzsource%cmnhname   = '2DREL'