source: trunk/LMDZ.GENERIC/libf/phystd/generic_tracer_index_mod.F90

Last change on this file was 2804, checked in by aslmd, 2 years ago

This is the RGCS scheme

Now we can use this scheme to take into account radiative effect of a Generic Condensable Specie. One needs to add the keyword is_rgcs=1 to the ModernTraceur?.def to activate the scheme + specify the number of aerogeneric specie in callphys.def. Also, one needs to go into suaer_corrk.F90 to add the Mie scattering file into the code for the specie you want. Currently, MnS, Fe, Cr, Mg2SiO4 are available. This has been tested on a Hot Jupiter case, might have issue in different configuration. Also, the numerical stability of the scheme is not guaranted when using strong scatterers such as silicates. Physical oscillations in the evaporation/condensation can occur, and need fixing. This will be for another commit.

File size: 3.3 KB
Line 
1module generic_tracer_index_mod
2
3!======================================================================C
4!
5!     This module is to get back the index of a generic condensable tracer
6!     for both the vap & ice parts of the tracer
7!
8!     call_ice_vap_generic is a boolean which tells if the tracer is condensable
9!     and allows, if it is condensable, to call the vap/ice pair
10!
11!     Noé Clément & Lucas Teinturier (2022)
12!
13!======================================================================C
14
15    implicit none
16
17    contains
18
19    subroutine generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic)
20
21        USE tracer_h, only: noms, is_condensable
22        integer nq,iq, igcm_generic_vap, igcm_generic_ice, i_shortname, ii
23        logical call_ice_vap_generic
24        character(len=10) :: shortname
25        if((is_condensable(iq)==1) .and. (index(noms(iq),"vap") .ne. 0) &
26                                                    .and. (index(noms(iq),"h2o") .eq. 0) .and. (index(noms(iq),"co2") .eq. 0)) then
27            ! Let's get the index of our tracers (we look for igcm _generic_vap and igcm_generic_ice)
28
29            igcm_generic_vap = iq
30            igcm_generic_ice = -1
31
32            ! ! ! We look for the corresponding ice traceur (before or after in the list of traceurs, maybe could be generalised to the whole list)
33            ! ! if (iq .ne. nq) then
34            ! !     if ((noms(iq)(1:len(trim(noms(iq)))-4) .eq. noms(iq+1)(1:len(trim(noms(iq+1)))-4)) .and. (index(noms(iq+1),"ice") .ne. 0)) then
35            ! !         igcm_generic_ice = iq+1
36            ! !     end if
37            ! ! end if
38            ! ! if ((iq .gt. 1)) then
39            ! !     if ((noms(iq)(1:len(trim(noms(iq)))-4) .eq. noms(iq-1)(1:len(trim(noms(iq-1)))-4)) .and. (index(noms(iq-1),"ice") .ne. 0)) then
40            ! !             igcm_generic_ice = iq-1
41            ! !     end if
42            ! ! end if
43            ! ! if (igcm_generic_ice .eq. -1) then
44            ! !     write(*,*) "ERROR : You set a vap traceur but you forgot to set the corresponding ice traceur, &
45            ! !     or the pair vap/ice is not written one after another in traceur.def"
46            ! ! endif
47
48            ! ! call_ice_vap_generic = .true.
49            i_shortname = index(noms(iq),'_')-1
50            shortname=noms(iq)(1:i_shortname)
51            ! We look for the corresponding ice traceur, first if it's before the vap traceur in pq
52            if (iq>1 ) then
53                do ii=1,iq-1
54                    if (index(noms(ii),adjustl(trim(shortname))) .ne. 0) then
55                        igcm_generic_ice=ii
56                        exit
57                    endif
58                ENDDO
59            endif
60            ! if we didn't find it before, we look after the vap tracer in pq
61            if (igcm_generic_ice .eq. -1 .and. (iq<nq)) then
62                do ii=iq+1,nq
63                    if (index(noms(ii),adjustl(trim(shortname))) .ne. 0) then
64                        igcm_generic_ice=ii
65                        exit
66                    endif
67                enddo
68            endif
69            call_ice_vap_generic = .true.
70        else
71
72            call_ice_vap_generic = .false.
73
74        end if ! if((is_condensable(iq)==1) .and. (index(noms(iq),"vap") .ne. 0))
75
76    end subroutine generic_tracer_index
77
78end module generic_tracer_index_mod
Note: See TracBrowser for help on using the repository browser.