source: trunk/LMDZ.TITAN/libf/phytitan/cond_muphy.F90 @ 3529

Last change on this file since 3529 was 3090, checked in by slebonnois, 14 months ago

BdeBatz? : Cleans microphysics and makes few corrections for physics

File size: 2.0 KB
Line 
1subroutine cond_muphy(ngrid,nlayer,pt,dqmuficond,dtlc)
2
3use comcstfi_mod, only: cpp ! Heat capacity of the atmosphere [J.kg-1.K-1]
4use tracer_h
5
6implicit none
7
8!====================================================================
9!
10!   Purpose
11!   -------
12!   The routine calculates the latent heat
13!   from condensation of species from microphysics
14!   
15!   INPUT
16!   -----
17!   ngrid      = Number of grid points in the chunk [-]
18!   nlayer     = Number of vertical layers [-]
19!   pt         = Temperature [K]
20!   dqmuficond = Trend of gas condensed in the microphysics [kg.kg_air-1.s-1]
21!
22!   OUTPUT
23!   ------
24!   dtlc = Condensation heating rate [K.s-1]
25!
26!   Author
27!   ------
28!   B. de Batz de Trenquelléon (07/2023)
29!====================================================================   
30
31
32!------------------------------------
33! 0. DECLARATIONS
34!------------------------------------
35
36! Inputs :
37!---------
38integer, intent(in) :: ngrid
39integer, intent(in) :: nlayer
40real, intent(in)    :: pt(ngrid,nlayer)
41real, intent(in)    :: dqmuficond(ngrid,nlayer,size(ices_indx))
42
43! Outputs :
44!----------
45real, intent(out) :: dtlc(ngrid,nlayer)
46
47! Local variables :
48!------------------
49integer :: iq
50real*8  :: Lc(ngrid,nlayer,size(ices_indx)) ! Condensation latente heat [J.kg-1]
51
52
53!------------------------------------
54! 1. INITIALISATION
55!------------------------------------
56
57! Variables set to 0 :
58! ~~~~~~~~~~~~~~~~~~~~
59dtlc(:,:) = 0.0D0
60Lc(:,:,:) = 0.0D0
61
62
63! Computes Lc for ices [J.kg-1] :
64! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
65call calc_condlh(ngrid,nlayer,pt,Lc)
66
67
68! Sum of condensation latent heat [J.kg_air-1.s-1] :
69! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
70do iq = 1, size(ices_indx)
71    dtlc(:,:) = dtlc(:,:) + (dqmuficond(:,:,iq) * Lc(:,:,iq))
72enddo
73
74
75! Condensation heating rate :
76! ~~~~~~~~~~~~~~~~~~~~~~~~~~~
77! If ice formation  : dqmuficond > 0 --> dtlc > 0
78! Else vaporisation : dqmuficond < 0 --> dtlc < 0
79dtlc(:,:) = dtlc(:,:) / cpp ! [K.s-1]
80
81return
82
83end subroutine cond_muphy
Note: See TracBrowser for help on using the repository browser.