Skip to content
Snippets Groups Projects
ini_budget.f90 159 KiB
Newer Older
!MNH_LIC Copyright 1995-2020 CNRS, Meteo-France and Universite Paul Sabatier
!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.
!-----------------------------------------------------------------
module mode_ini_budget

  implicit none

  private

  public :: Ini_budget

contains

!     #################################################################
      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,                                   &
      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
!!    ------------------ 
!!      Module MODD_PARAMETERS: JPBUMAX,JPBUPROMAX
!!
!!      Module MODD_CONF: CCONF        
!!
!!      Module MODD_DYN:  XSEGLEN 
!!
!!
!!    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 XBUSURF 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 negative 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 
!!                   04/2016 (C.LAC) negative contribution to the budget splitted between advection, turbulence and microphysics for KHKO/C2R2
Gaelle TANGUY's avatar
Gaelle TANGUY committed
!!      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
!!  Philippe 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 02-03/2020: use the new data structures and subroutines for budgets
!-------------------------------------------------------------------------------
!
!*       0.    DECLARATIONS
!              ------------ 
!
use modd_2d_frc,        only: l2d_adv_frc, l2d_rel_frc
use modd_budget
use modd_ch_aerosol,    only: lorilam
use modd_conf,          only: l1d, lcartesian, lforcing, lthinshell, nmodel
use modd_dust,          only: ldust
use modd_dyn,           only: lcorio
use modd_elec_descr,    only: linductive, lrelax2fw_ion
use modd_field,         only: TYPEREAL
use modd_nsv,           only: nsv_aerbeg, nsv_aerend, nsv_c2r2beg, nsv_c2r2end, nsv_chembeg, nsv_chemend,          &
                              nsv_elecbeg, nsv_elecend,                                                            &
                              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_scavmass,                            &
                              nsv_user
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: lacti_lima => lacti, lcold_lima => lcold, ldepoc_lima => ldepoc, lhail_lima => lhail, &
                             lhhoni_lima => lhhoni, lmeyers_lima => lmeyers, lnucl_lima => lnucl, lptsplit,        &
                             lrain_lima => lrain, lscav_lima => lscav, lsedc_lima => lsedc, lsedi_lima => lsedi,   &
                             lsnow_lima => lsnow, lwarm_lima => lwarm,                                             &
                             nmod_ccn, nmod_ifn, nmod_imm
use modd_salt,         only: lsalt
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
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
!                                                     
INTEGER, DIMENSION(JPBUMAX,JPBUPROMAX+1) :: IPROACTV      ! switches set by the
                                                          ! user for process 
                                                          ! activation
INTEGER :: JI, JJ, JK , JJJ                               ! 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
                                                          ! = NBUWRNB in MASK case
INTEGER :: IBUDIM3                                        ! third dimension of the budget arrays
                                                          ! = NBUKMAX in CART case
                                                          ! = NBUMASK in MASK case
LOGICAL :: GERROR                                         ! switch for error in
                                                          ! budget specifcation
CHARACTER(LEN=7), DIMENSION(JPBUMAX) :: YEND_COMMENT      ! last part of comment
                                                          ! for budgets records
CHARACTER(LEN=6), DIMENSION(JPBUMAX,JPBUPROMAX) :: YWORK2 ! used for
                                                          ! concatenattion of  
                                                          ! comments for budgets
CHARACTER(LEN=40)                  :: YSTRING
character(len=3)    :: ybudgetnum
INTEGER             :: ILEN
INTEGER             :: JSV               ! loop indice for the SVs
INTEGER             :: IBUPROCNBR_SV_MAX ! Max number of processes for the SVs
INTEGER             :: IINFO_ll ! return status of the interface routine
integer             :: ibudget
integer             :: isourcesmax          ! Maximum number of source terms in a budget
integer             :: igroup
logical             :: gcond
logical             :: gtmp
type(tbusourcedata) :: tzsource ! Used to prepare metadate of source terms

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

nbudgets = NBUDGET_SV1 - 1 + ksv
allocate( tbudgets( nbudgets ) )
!
!*       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') THEN              ! cartesian case only
!
  NBUWRNB = NINT (XBUWRI / XBULEN)  ! only after NBUWRNB budget periods, we write the
                                    ! result on the FM_FILE   
  IF (LBU_ICP) THEN 
    NBUIMAX_ll = 1
  ELSE
    NBUIMAX_ll = NBUIH - NBUIL +1
  END IF
  IF (LBU_JCP) THEN 
    NBUJMAX_ll = 1
  ELSE
    NBUJMAX_ll = NBUJH - NBUJL +1
  END IF
!
  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.
  NBUWRNB = NINT (XBUWRI / XBULEN)  ! only after NBUWRNB budget periods, we write the
                                    ! 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( XBUSURF( IIU, IJU, NBUMASK, NBUWRNB) )
  XBUSURF(:,:,:,:) = 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=NBUWRNB
  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
!              ------------------------------------------------
!
ALLOCATE( NBUPROCNBR(JPBUMAX) )
ALLOCATE( NBUPROCCTR(JPBUMAX) )
ALLOCATE( CBUACTION(JPBUMAX, JPBUPROMAX) )
ALLOCATE( CBUCOMMENT(JPBUMAX, JPBUPROMAX) )
NBUPROCCTR(:) = 0
NBUCTR_ACTV(:) = 0
NBUPROCNBR(:) = 0
CBUACTION(:,:) = 'OF' 
CBUCOMMENT(:,:) = ' '
LBU_BEG =.TRUE. 
!
!-------------------------------------------------------------------------------
!
!*       3.    INITALIZE VARIABLES
!              -------------------
!
IPROACTV(:,:) = 3
IPROACTV(:,4) = 1
IPROACTV(:,JPBUPROMAX+1) = 0
GERROR=.FALSE.
YWORK2(:,:) = ' '
YEND_COMMENT(:) = ' '
!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


tbudgets(NBUDGET_U)%cname    = "BU_RU"
tbudgets(NBUDGET_U)%ccomment = "Budget for U"
  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)
  isourcesmax = 18
  tbudgets(NBUDGET_U)%nsourcesmax = isourcesmax
  allocate( tbudgets(NBUDGET_U)%tsources(isourcesmax) )
  allocate( tbudgets(NBUDGET_U)%xtmpstore(ibudim1, ibudim2, ibudim3) )
  tzsource%ccomment = 'Budget of momentum along X axis'
  tzsource%ngrid    = 2
  tzsource%cmnhname  = 'INIF'
  tzsource%clongname = 'initial state'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, 1, odonotinit = .true., ooverwrite = .true. )
  tzsource%cmnhname  = 'ENDF'
  tzsource%clongname = 'final state'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, 1, odonotinit = .true., ooverwrite = .true. )
  tzsource%cmnhname  = 'AVEF'
  tzsource%clongname = 'averaged state'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, 1, odonotinit = .true., ooverwrite = .false. )
  tzsource%cmnhname  = 'ASSE'
  tzsource%clongname = 'time filter (Asselin)'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, nasseu )
  tzsource%cmnhname  = 'NEST'
  tzsource%clongname = 'nesting'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, nnestu )
  tzsource%cmnhname  = 'FRC'
  tzsource%clongname = 'forcing'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, nfrcu )
  tzsource%cmnhname  = 'NUD'
  tzsource%clongname = 'nudging'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, nnudu )
  gcond = .not.l1d .and. .not.lcartesian
  tzsource%cmnhname  = 'CURV'
  tzsource%clongname = 'curvature'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, ncurvu )
  tzsource%cmnhname  = 'COR'
  tzsource%clongname = 'Coriolis'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, ncoru )
  tzsource%cmnhname  = 'DIF'
  tzsource%clongname = 'numerical diffusion'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, ndifu )
  gcond = ohorelax_uvwth .or. ove_relax .or. ove_relax_grd
  tzsource%cmnhname  = 'REL'
  tzsource%clongname = 'relaxation'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, nrelu )
  tzsource%cmnhname  = 'DRAG'
  tzsource%clongname = 'drag force'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, ndragu )
  tzsource%cmnhname  = 'VTURB'
  tzsource%clongname = 'vertical turbulent diffusion'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, nvturbu )
  gcond = hturb == 'TKEL' .and. HTURBDIM == '3DIM'
  tzsource%cmnhname  = 'HTURB'
  tzsource%clongname = 'horizontal turbulent diffusion'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, nhturbu )
  tzsource%cmnhname  = 'MAFL'
  tzsource%clongname = 'mass flux'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, nmaflu )
  tzsource%cmnhname  = 'VISC'
  tzsource%clongname = 'viscosity'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, nviscu )
  tzsource%cmnhname  = 'ADV'
  tzsource%clongname = 'advection'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, nadvu )
  tzsource%cmnhname  = 'PRES'
  tzsource%clongname = 'pressure'
  call Budget_source_add( tbudgets(NBUDGET_U), tzsource, gcond, npresu )
! Budget of RV
tbudgets(NBUDGET_V)%cname    = "BU_RV"
tbudgets(NBUDGET_V)%ccomment = "Budget for V"
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)
  isourcesmax = 18
  tbudgets(NBUDGET_V)%nsourcesmax = isourcesmax
  allocate( tbudgets(NBUDGET_V)%tsources(isourcesmax) )
  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
  gcond = .true.
  tzsource%cmnhname  = 'INIF'
  tzsource%clongname = 'initial state'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, 1, odonotinit = .true., ooverwrite = .true. )
  gcond = .true.
  tzsource%cmnhname  = 'ENDF'
  tzsource%clongname = 'final state'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, 1, odonotinit = .true., ooverwrite = .true. )
  gcond = .true.
  tzsource%cmnhname  = 'AVEF'
  tzsource%clongname = 'averaged state'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, 1, odonotinit = .true., ooverwrite = .false. )
  gcond = .true.
  tzsource%cmnhname  = 'ASSE'
  tzsource%clongname = 'time filter (Asselin)'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, nassev )
  gcond = nmodel > 1
  tzsource%cmnhname  = 'NEST'
  tzsource%clongname = 'nesting'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, nnestv )
  gcond = lforcing
  tzsource%cmnhname  = 'FRC'
  tzsource%clongname = 'forcing'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, nfrcv )
  gcond = onudging
  tzsource%cmnhname  = 'NUD'
  tzsource%clongname = 'nudging'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, nnudv )
  gcond = .not.l1d .and. .not.lcartesian
  tzsource%cmnhname  = 'CURV'
  tzsource%clongname = 'curvature'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, ncurvv )
  gcond = lcorio
  tzsource%cmnhname  = 'COR'
  tzsource%clongname = 'Coriolis'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, ncorv )
  gcond = onumdifu
  tzsource%cmnhname  = 'DIF'
  tzsource%clongname = 'numerical diffusion'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, ndifv )
  gcond = ohorelax_uvwth .or. ove_relax .or. ove_relax_grd
  tzsource%cmnhname  = 'REL'
  tzsource%clongname = 'relaxation'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, nrelv )
  gcond = odragtree
  tzsource%cmnhname  = 'DRAG'
  tzsource%clongname = 'drag force'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, ndragv )
  gcond = hturb == 'TKEL'
  tzsource%cmnhname  = 'VTURB'
  tzsource%clongname = 'vertical turbulent diffusion'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, nvturbv )
  gcond = hturb == 'TKEL' .and. HTURBDIM == '3DIM'
  tzsource%cmnhname  = 'HTURB'
  tzsource%clongname = 'horizontal turbulent diffusion'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, nhturbv )
  gcond = hsconv == 'EDKF'
  tzsource%cmnhname  = 'MAFL'
  tzsource%clongname = 'mass flux'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, nmaflv )
  gcond = lvisc .and. lvisc_uvw
  tzsource%cmnhname  = 'VISC'
  tzsource%clongname = 'viscosity'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, nviscv )
  gcond = .true.
  tzsource%cmnhname  = 'ADV'
  tzsource%clongname = 'advection'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, nadvv )
  gcond = .true.
  tzsource%cmnhname  = 'PRES'
  tzsource%clongname = 'pressure'
  call Budget_source_add( tbudgets(NBUDGET_V), tzsource, gcond, npresv )
end if
! Budget of RW
tbudgets(NBUDGET_W)%cname    = "BU_RW"
tbudgets(NBUDGET_W)%ccomment = "Budget for W"
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 Y 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)
  isourcesmax = 17
  tbudgets(NBUDGET_W)%nsourcesmax = isourcesmax
  allocate( tbudgets(NBUDGET_W)%tsources(isourcesmax) )
  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
  gcond = .true.
  tzsource%cmnhname  = 'INIF'
  tzsource%clongname = 'initial state'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, 1, odonotinit = .true., ooverwrite = .true. )
  gcond = .true.
  tzsource%cmnhname  = 'ENDF'
  tzsource%clongname = 'final state'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, 1, odonotinit = .true., ooverwrite = .true. )
  gcond = .true.
  tzsource%cmnhname  = 'AVEF'
  tzsource%clongname = 'averaged state'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, 1, odonotinit = .true., ooverwrite = .false. )
  gcond = .true.
  tzsource%cmnhname  = 'ASSE'
  tzsource%clongname = 'time filter (Asselin)'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, nassew )
  gcond = nmodel > 1
  tzsource%cmnhname  = 'NEST'
  tzsource%clongname = 'nesting'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, nnestw )
  gcond = lforcing
  tzsource%cmnhname  = 'FRC'
  tzsource%clongname = 'forcing'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, nfrcw )
  gcond = onudging
  tzsource%cmnhname  = 'NUD'
  tzsource%clongname = 'nudging'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, nnudw )
  gcond = .not.l1d .and. .not.lcartesian .and. .not.lthinshell
  tzsource%cmnhname  = 'CURV'
  tzsource%clongname = 'curvature'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, ncurvw )
  gcond = lcorio .and. .not.l1d .and. .not.lthinshell
  tzsource%cmnhname  = 'COR'
  tzsource%clongname = 'Coriolis'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, ncorw )
  gcond = onumdifu
  tzsource%cmnhname  = 'DIF'
  tzsource%clongname = 'numerical diffusion'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, ndifw )
  gcond = ohorelax_uvwth .or. ove_relax .or. ove_relax_grd
  tzsource%cmnhname  = 'REL'
  tzsource%clongname = 'relaxation'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, nrelw )
Gaelle TANGUY's avatar
Gaelle TANGUY committed

  gcond = hturb == 'TKEL'
  tzsource%cmnhname  = 'VTURB'
  tzsource%clongname = 'vertical turbulent diffusion'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, nvturbw )
  gcond = hturb == 'TKEL' .and. HTURBDIM == '3DIM'
  tzsource%cmnhname  = 'HTURB'
  tzsource%clongname = 'horizontal turbulent diffusion'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, nhturbw )
  gcond = lvisc .and. lvisc_uvw
  tzsource%cmnhname  = 'VISC'
  tzsource%clongname = 'viscosity'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, nviscw )
  gcond = .true.
  tzsource%cmnhname  = 'GRAV'
  tzsource%clongname = 'gravity'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, ngravw )
  gcond = .true.
  tzsource%cmnhname  = 'ADV'
  tzsource%clongname = 'advection'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, nadvw )
  gcond = .true.
  tzsource%cmnhname  = 'PRES'
  tzsource%clongname = 'pressure'
  call Budget_source_add( tbudgets(NBUDGET_W), tzsource, gcond, npresw )
end if
! Budget of RTH
tbudgets(NBUDGET_TH)%cname    = "BU_RTH"
tbudgets(NBUDGET_TH)%ccomment = "Budget for potential temperature"
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)
  isourcesmax = 52
  tbudgets(NBUDGET_TH)%nsourcesmax = isourcesmax
  allocate( tbudgets(NBUDGET_TH)%tsources(isourcesmax) )
  allocate( tbudgets(NBUDGET_TH)%xtmpstore(ibudim1, ibudim2, ibudim3) )
  tbudgets(NBUDGET_TH)%tsources(:)%ngroup = 0
  tzsource%ccomment = 'Budget of potential temperature'
  tzsource%ngrid    = 1
  gcond = .true.
  tzsource%cmnhname  = 'INIF'
  tzsource%clongname = 'initial state'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, 1, odonotinit = .true., ooverwrite = .true. )
  gcond = .true.
  tzsource%cmnhname  = 'ENDF'
  tzsource%clongname = 'final state'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, 1, odonotinit = .true., ooverwrite = .true. )
  gcond = .true.
  tzsource%cmnhname  = 'AVEF'
  tzsource%clongname = 'averaged state'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, 1, odonotinit = .true., ooverwrite = .false. )
  gcond = .true.
  tzsource%cmnhname  = 'ASSE'
  tzsource%clongname = 'time filter (Asselin)'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nasseth )
  gcond = nmodel > 1
  tzsource%cmnhname  = 'NEST'
  tzsource%clongname = 'nesting'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nnestth )
  gcond = lforcing
  tzsource%cmnhname  = 'FRC'
  tzsource%clongname = 'forcing'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nfrcth )
  gcond = l2d_adv_frc
  tzsource%cmnhname  = '2DADV'
  tzsource%clongname = 'advective forcing'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, n2dadvth )
  gcond = l2d_rel_frc
  tzsource%cmnhname  = '2DREL'
  tzsource%clongname = 'relaxation forcing'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, n2drelth )
  gcond = onudging
  tzsource%cmnhname  = 'NUD'
  tzsource%clongname = 'nudging'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nnudth )
  gcond = krr > 0 .and. .not.l1d
  tzsource%cmnhname  = 'PREF'
  tzsource%clongname = 'reference pressure'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nprefth )
  gcond = onumdifth
  tzsource%cmnhname  = 'DIF'
  tzsource%clongname = 'numerical diffusion'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, ndifth )
  gcond = ohorelax_uvwth .or. ove_relax .or. ove_relax_grd
  tzsource%cmnhname  = 'REL'
  tzsource%clongname = 'relaxation'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nrelth )
  gcond = hrad /= 'NONE'
  tzsource%cmnhname  = 'RAD'
  tzsource%clongname = 'radiation'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nradth )
  gcond = hdconv == 'KAFR' .OR. hsconv == 'KAFR'
  tzsource%cmnhname  = 'DCONV'
  tzsource%clongname = 'KAFR convection'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, ndconvth )
  gcond = hturb == 'TKEL'
  tzsource%cmnhname  = 'VTURB'
  tzsource%clongname = 'vertical turbulent diffusion'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nvturbth )
  gcond = hturb == 'TKEL' .and. HTURBDIM == '3DIM'
  tzsource%cmnhname  = 'HTURB'
  tzsource%clongname = 'horizontal turbulent diffusion'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nhturbth )
  gcond = hturb == 'TKEL'
  tzsource%cmnhname  = 'DISSH'
  tzsource%clongname = 'dissipation'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, ndisshth )
  gcond = hturb == 'TKEL' .and. ( hcloud == 'KHKO' .or.  hcloud == 'C2R2' )
  tzsource%cmnhname  = 'NETUR'
  tzsource%clongname = 'negative correction induced by turbulence'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nneturth )
  gcond = hsconv == 'EDKF'
  tzsource%cmnhname  = 'MAFL'
  tzsource%clongname = 'mass flux'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nmaflth )
  gcond = lvisc .and. lvisc_th
  tzsource%cmnhname  = 'VISC'
  tzsource%clongname = 'viscosity'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nviscth )
  gcond = .true.
  tzsource%cmnhname  = 'ADV'
  tzsource%clongname = 'total advection'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nadvth )

  gcond = hcloud == 'KHKO' .or. hcloud == 'C2R2'
  tzsource%cmnhname  = 'NEADV'
  tzsource%clongname = 'negative correction induced by advection'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nneadvth )

  gcond = hcloud /= 'NONE' .and. hcloud /= 'KHKO' .and. hcloud /= 'C2R2'
  tzsource%cmnhname  = 'NEGA'
  tzsource%clongname = 'negative correction'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nnegath )

  gcond = hcloud == 'LIMA' .and. lptsplit
  tzsource%cmnhname  = 'SEDI'
  tzsource%clongname = 'heat transport by hydrometeors sedimentation'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nsedith )

  gtmp = cactccn == 'ABRK' .and. (lorilam .or. ldust .or. lsalt )
  gcond =      ( hcloud      == 'LIMA' .and. lwarm_lima .and. lacti_lima .and. nmod_ccn >= 1 )     &
          .or.   hcloud(1:3) == 'ICE'                                                         &
          .or. ( hcloud      == 'C2R2' .and. ( gtmp .or. ( .not.gtmp .and. .not.lsupsat_c2r2 ) ) ) &
          .or. ( hcloud      == 'KHKO' .and. ( gtmp .or. ( .not.gtmp .and. .not.lsupsat_c2r2 ) ) )
  tzsource%cmnhname  = 'HENU'
  tzsource%clongname = 'heterogeneous nucleation'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nhenuth )

  gcond =      ( hcloud      == 'LIMA' .and. ( ( .not. lptsplit .and. lwarm_lima .and. lrain_lima ) .or. lptsplit ) ) &
          .or. ( hcloud(1:3) == 'ICE'  .and. lwarm_ice )                                                              &
          .or. ( hcloud      == 'C2R2' .and. lrain_c2r2 )                                                             &
          .or. ( hcloud      == 'KHKO' .and. lrain_c2r2 )                                                             &
          .or.   hcloud      == 'KESS'
  tzsource%cmnhname  = 'REVA'
  tzsource%clongname = 'rain evaporation'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nrevath )

  gcond = hcloud == 'LIMA' .and. lcold_lima .and. lnucl_lima
  tzsource%cmnhname  = 'HIND'
  tzsource%clongname = 'heterogeneous nucleation by deposition'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nhindth )

  gcond = hcloud == 'LIMA' .and. lcold_lima .and. lnucl_lima
  tzsource%cmnhname  = 'HINC'
  tzsource%clongname = 'heterogeneous nucleation by contact'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nhincth )

  gcond = hcloud(1:3) == 'ICE'
  tzsource%cmnhname  = 'HON'
  tzsource%clongname = 'homogeneous nucleation'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nhonth )

  gcond = hcloud == 'LIMA' .and. lcold_lima .and. lnucl_lima .and. lhhoni_lima .and. nmod_ccn >= 1
  tzsource%cmnhname  = 'HONH'
  tzsource%clongname = 'haze homogeneous nucleation'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nhonhth )

  gcond = hcloud == 'LIMA' .and. ( lptsplit .or. ( lcold_lima .and. lwarm_lima .and. lnucl_lima ) )
  tzsource%cmnhname  = 'HONC'
  tzsource%clongname = 'droplet homogeneous freezing'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nhoncth )

  gcond = hcloud == 'LIMA' .and. ( lptsplit .or. ( lcold_lima .and. lwarm_lima .and. lnucl_lima .and. lrain_lima ) )
  tzsource%cmnhname  = 'HONR'
  tzsource%clongname = 'raindrop homogeneous freezing'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nhonrth )

  gcond = hcloud(1:3) == 'ICE'
  tzsource%cmnhname  = 'SFR'
  tzsource%clongname = 'spontaneous freezing'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nsfrth )

  gcond = ( hcloud == 'LIMA' .and. ( lptsplit .or. ( lcold_lima .and. lsnow_lima ) ) ) .or. hcloud(1:3) == 'ICE'
  tzsource%cmnhname  = 'DEPS'
  tzsource%clongname = 'deposition on snow'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, ndepsth )

  gcond = ( hcloud == 'LIMA' .and. ( lptsplit .or. ( lcold_lima .and. lwarm_lima .and. lsnow_lima ) ) ) .or. hcloud(1:3) == 'ICE'
  tzsource%cmnhname  = 'DEPG'
  tzsource%clongname = 'deposition on graupel'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, ndepgth )

  gcond = ( hcloud == 'LIMA' .and. ( lptsplit .or. ( lcold_lima .and. lwarm_lima ) ) ) .or. hcloud(1:3) == 'ICE'
  tzsource%cmnhname  = 'IMLT'
  tzsource%clongname = 'ice melting'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nimltth )

  gcond = ( hcloud == 'LIMA' .and. ( lptsplit .or. (lcold_lima .and. lwarm_lima) ) ) .or. hcloud(1:3) == 'ICE'
  tzsource%cmnhname  = 'BERFI'
  tzsource%clongname = 'Bergeron-Findeisen'
  call Budget_source_add( tbudgets(NBUDGET_TH), tzsource, gcond, nberfith )