source: trunk/LMDZ.TITAN/libf/phytitan/calmufi.F90 @ 1926

Last change on this file since 1926 was 1926, checked in by jvatant, 7 years ago

1) Microphysics diags / outputs :


+ Add supplementary diagnostics outputs for microphysics ( precip, flux, rc ... ) ( new muphy_diag.F90 module )
+ Correct the outputs of microphys tracers to be in X/m-3 to be comparable to "standard values"

+ Also update the deftank callphys.def with latest revs modifs for microphysics

2) Condensation / chemistry updates :


+ Moved chemistry AFTER microphysics

  • To have mufi condensation before photochem
  • Chemistry called last coherent with the fact that it brings back fields to equilibrium

+ If 2D chemistry, make zonally averaged fields go through mufi and chem condensation

to have non saturated profiles in input of photochemistry
( other 'short' processes neglected as 2D -> no diurnal cycle, just seasonal evolution )

+ Also corrected the positivity check ( took Mars GCM syntax ) after chemistry ( could previously lead to negs )

3) Noticed a weird behaviour ( bug? ) :


+ If generalize the use of arrays *_indx for tracers, to get rid of ugly "iq+nmicro",

it ends up with weird results / crash in optim mode ( ok in debug ) but didn't find out why ...

--JVO

  • Property svn:executable set to *
File size: 7.9 KB
RevLine 
[1897]1
2
[1908]3SUBROUTINE calmufi(dt, plev, zlev, play, zlay, temp, pq, zdqfi, zdq)
[1795]4  !! Interface subroutine to YAMMS model for Titan LMDZ GCM.
5  !!
[1908]6  !! The subroutine computes the microphysics processes for a single vertical column.
[1795]7  !!
8  !! - All input vectors are assumed to be defined from GROUND to TOP of the atmosphere.
9  !! - All output vectors are defined from GROUND to TOP of the atmosphere.
10  !! - Only tendencies are returned.
11  !!
12  !! @important
13  !! The method assumes global initialization of YAMMS model (and extras) has been already
14  !! done elsewhere.
15  !!
16  !! Authors : J.Burgalat, J.Vatant d'Ollone - 2017
17  !!
18  USE MMP_GCM
19  USE tracer_h
20  USE comcstfi_mod, only : g
21  USE callkeys_mod, only : callclouds
[1926]22  USE muphy_diag
[1795]23  IMPLICIT NONE
24
[1902]25  REAL(kind=8), INTENT(IN) :: dt  !! Physics timestep (s).
26 
[1795]27  REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: plev  !! Pressure levels (Pa).
28  REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: zlev  !! Altitude levels (m).
29  REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: play  !! Pressure layers (Pa).
30  REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: zlay  !! Altitude at the center of each layer (m).
31  REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: temp  !! Temperature at the center of each layer (K).
32
[1908]33  REAL(kind=8), DIMENSION(:,:,:), INTENT(IN)  :: pq    !! Tracers (\(kg.kg^{-1}}\)).
34  REAL(kind=8), DIMENSION(:,:,:), INTENT(IN)  :: zdqfi !! Tendency from former processes for tracers (\(kg.kg^{-1}}\)).
35  REAL(kind=8), DIMENSION(:,:,:), INTENT(OUT) :: zdq   !! Microphysical tendency for tracers (\(kg.kg^{-1}}\)).
[1795]36 
[1908]37  REAL(kind=8), DIMENSION(:,:,:), ALLOCATABLE :: zq !! Local tracers updated from former processes (\(kg.kg^{-1}}\)).
38 
[1795]39  REAL(kind=8), DIMENSION(:), ALLOCATABLE :: m0as !! 0th order moment of the spherical mode (\(m^{-2}\)).
40  REAL(kind=8), DIMENSION(:), ALLOCATABLE :: m3as !! 3rd order moment of the spherical mode (\(m^{3}.m^{-2}\)).
41  REAL(kind=8), DIMENSION(:), ALLOCATABLE :: m0af !! 0th order moment of the fractal mode (\(m^{-2}\)).
42  REAL(kind=8), DIMENSION(:), ALLOCATABLE :: m3af !! 3rd order moment of the fractal mode (\(m^{3}.m^{-2}\)).
43
44  REAL(kind=8), DIMENSION(:), ALLOCATABLE   :: m0n  !! 0th order moment of the CCN distribution (\(m^{-2}\)).
45  REAL(kind=8), DIMENSION(:), ALLOCATABLE   :: m3n  !! 3rd order moment of the CCN distribution (\(m^{3}.m^{-2}\)).
46  REAL(kind=8), DIMENSION(:,:), ALLOCATABLE :: m3i  !! 3rd order moments of the ice components (\(m^{3}.m^{-2}\)).
47  REAL(kind=8), DIMENSION(:,:), ALLOCATABLE :: gazs !! Condensible species gazs molar fraction (\(mol.mol^{-1}\)).
48
49  REAL(kind=8), DIMENSION(:), ALLOCATABLE   :: dm0as !! Tendency of the 0th order moment of the spherical mode distribution (\(m^{-2}\)).
50  REAL(kind=8), DIMENSION(:), ALLOCATABLE   :: dm3as !! Tendency of the 3rd order moment of the spherical mode distribution (\(m^{3}.m^{-2}\)).
51  REAL(kind=8), DIMENSION(:), ALLOCATABLE   :: dm0af !! Tendency of the 0th order moment of the fractal mode distribution (\(m^{-2}\)).
52  REAL(kind=8), DIMENSION(:), ALLOCATABLE   :: dm3af !! Tendency of the 3rd order moment of the fractal mode distribution (\(m^{3}.m^{-2}\)).
53  REAL(kind=8), DIMENSION(:), ALLOCATABLE   :: dm0n  !! Tendency of the 0th order moment of the _CCN_ distribution (\(m^{-2}\)).
54  REAL(kind=8), DIMENSION(:), ALLOCATABLE   :: dm3n  !! Tendency of the 3rd order moment of the _CCN_ distribution (\(m^{3}.m^{-2}\)).
55  REAL(kind=8), DIMENSION(:,:), ALLOCATABLE :: dm3i  !! Tendencies of the 3rd order moments of each ice components (\(m^{3}.m^{-2}\)).
56  REAL(kind=8), DIMENSION(:,:), ALLOCATABLE :: dgazs !! Tendencies of each condensible gaz species !(\(mol.mol^{-1}\)).
57
58  REAL(kind=8), DIMENSION(:), ALLOCATABLE ::  int2ext
59  TYPE(error) :: err
60
61  INTEGER :: ilon, i,nices
[1908]62  INTEGER :: nq,nlon,nlay
[1795]63
64  ! Read size of arrays
[1908]65  nq    = size(pq,DIM=3)
[1795]66  nlon  = size(play,DIM=1)
67  nlay  = size(play,DIM=2)
68  nices = size(ices_indx)
69  ! Conversion intensive to extensive
70  ALLOCATE( int2ext(nlay) ) 
71
[1908]72  ! Allocate arrays
73  ALLOCATE( zq(nlon,nlay,nq) )
[1795]74
75  ALLOCATE( m0as(nlay) )
76  ALLOCATE( m3as(nlay) )
77  ALLOCATE( m0af(nlay) )
78  ALLOCATE( m3af(nlay) )
79  ALLOCATE( m0n(nlay) )
80  ALLOCATE( m3n(nlay) )
81  ALLOCATE( m3i(nlay,nices) )
82  ALLOCATE( gazs(nlay,nices) )
83
84  ALLOCATE( dm0as(nlay) )
85  ALLOCATE( dm3as(nlay) )
86  ALLOCATE( dm0af(nlay) )
87  ALLOCATE( dm3af(nlay) )
88  ALLOCATE( dm0n(nlay) )
89  ALLOCATE( dm3n(nlay) )
90  ALLOCATE( dm3i(nlay,nices) )
91  ALLOCATE( dgazs(nlay,nices) )
92
93  ! Initialization of zdq here since intent=out and no action performed on every tracers
94  zdq(:,:,:) = 0.0
95
[1908]96  ! Initialize tracers updated with former processes from physics
97  zq(:,:,:) = pq(:,:,:) + zdqfi(:,:,:)*dt
98
99  ! Loop on horizontal grid points
[1795]100  DO ilon = 1, nlon
101    ! Convert tracers to extensive ( except for gazs where we work with molar mass ratio )
102    ! We suppose a given order of tracers !
[1897]103    int2ext(:) = ( plev(ilon,1:nlay)-plev(ilon,2:nlay+1) ) / g
[1908]104    m0as(:) = zq(ilon,:,1) * int2ext(:)
105    m3as(:) = zq(ilon,:,2) * int2ext(:)
106    m0af(:) = zq(ilon,:,3) * int2ext(:)
107    m3af(:) = zq(ilon,:,4) * int2ext(:)
[1795]108   
109    if (callclouds) then ! if call clouds
[1908]110      m0n(:) = zq(ilon,:,5) * int2ext(:)
111      m3n(:) = zq(ilon,:,6) * int2ext(:)
[1795]112      do i=1,nices
[1908]113        m3i(:,nices) = zq(ilon,:,6+i) * int2ext(:)
114        gazs(:,i)    = zq(ilon,:,ices_indx(i)) / rat_mmol(ices_indx(i)) ! For gazs we work on the full tracer array !!
[1795]115        ! We use the molar mass ratio from GCM in case there is discrepancy with the mm one
116      enddo
117    endif
118
119    ! Initialize YAMMS atmospheric column
120    err = mm_column_init(plev(ilon,:),zlev(ilon,:),play(ilon,:),zlay(ilon,:),temp(ilon,:)) ; IF (err /= 0) call abort_program(err)
121    ! Initialize YAMMS aerosols moments column
122    err = mm_aerosols_init(m0as,m3as,m0af,m3af) ; IF (err /= 0) call abort_program(err)
123    IF (callclouds) THEN ! call clouds
124      err = mm_clouds_init(m0n,m3n,m3i,gazs) ; IF (err /= 0) call abort_program(err)
125    ENDIF
126
127    ! Check on size (???)
128
129    ! initializes tendencies:
130    !dm0as(:) = 0._mm_wp ; dm3as(:) = 0._mm_wp ; dm0af(:) = 0._mm_wp ; dm3af(:) = 0._mm_wp
131    !dm0n(:) = 0._mm_wp ; dm3n(:) = 0._mm_wp ; dm3i(:,:) = 0._mm_wp ; dgazs(:,:) = 0._mm_wp
132
133    dm0as(:) = 0.0 ; dm3as(:) = 0.0 ; dm0af(:) = 0.0 ; dm3af(:) = 0.0
134    dm0n(:) = 0.0 ; dm3n(:) = 0.0 ; dm3i(:,:) = 0.0 ; dgazs(:,:) = 0.0
135    ! call microphysics
136
137    IF (callclouds) THEN ! call clouds
138      IF(.NOT.mm_muphys(dm0as,dm3as,dm0af,dm3af,dm0n,dm3n,dm3i,dgazs)) &
139        call abort_program(error("mm_muphys aborted -> initialization not done !",-1))
140    ELSE
141      IF (.NOT.mm_muphys(dm0as,dm3as,dm0af,dm3af)) &
142        call abort_program(error("mm_muphys aborted -> initialization not done !",-1))
143    ENDIF
[1926]144    ! save diags (if no clouds, relevant arrays will be set to 0 !)
145    call mm_diagnostics(mmd_aer_prec(ilon),mmd_aer_s_flux(ilon,:),mmd_aer_f_flux(ilon,:),  &
146                        mmd_ccn_prec(ilon),mmd_ccn_flux(ilon,:), mmd_ice_prec(ilon,:),   &
147                        mmd_ice_fluxes(ilon,:,:),mmd_gazs_sat(ilon,:,:))
148    call mm_get_radii(mmd_rc_sph(ilon,:),mmd_rc_fra(ilon,:),mmd_rc_cld(ilon,:))
[1795]149
150    ! Convert tracers back to intensives ( except for gazs where we work with molar mass ratio )
151    ! We suppose a given order of tracers !
[1897]152 
[1795]153    zdq(ilon,:,1) = dm0as(:) / int2ext(:)
154    zdq(ilon,:,2) = dm3as(:) / int2ext(:)
155    zdq(ilon,:,3) = dm0af(:) / int2ext(:)
156    zdq(ilon,:,4) = dm3af(:) / int2ext(:)
157   
158    if (callclouds) then ! if call clouds
159      zdq(ilon,:,5) = dm0n(:) / int2ext(:)
160      zdq(ilon,:,6) = dm3n(:) / int2ext(:)
161      do i=1,nices
162        zdq(ilon,:,6+i) = dm3i(:,nices) / int2ext(:)
[1904]163        zdq(ilon,:,ices_indx(i)) = dgazs(:,i) * rat_mmol(ices_indx(i)) ! For gazs we work on the full tracer array !!
[1795]164        ! We use the molar mass ratio from GCM in case there is discrepancy with the mm one
165      enddo
166    endif
167  END DO ! loop on ilon
168
[1902]169  ! YAMMS gives a tendency which is integrated for all the timestep but in the GCM
170  ! we want to have routines spitting tendencies in s-1 -> let's divide !
171  zdq(:,:,:) = zdq(:,:,:) / dt
172
[1795]173END SUBROUTINE calmufi
Note: See TracBrowser for help on using the repository browser.