source: trunk/LMDZ.PLUTO/libf/phypluto/generic_tracer_index_mod.F90 @ 3506

Last change on this file since 3506 was 3275, checked in by afalco, 8 months ago

Pluto PCM:
Changed _vap to _gas;
Included surfprop.F90;
callcorrk includes methane
AF

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_gas_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_gas,igcm_generic_ice,call_ice_gas_generic)
20
21        USE tracer_h, only: noms, is_condensable
22        integer nq,iq, igcm_generic_gas, igcm_generic_ice, i_shortname, ii
23        logical call_ice_gas_generic
24        character(len=10) :: shortname
25        if((is_condensable(iq)==1) .and. (index(noms(iq),"gas") .ne. 0) &
26                                                    .and. (index(noms(iq),"h2o") .eq. 0) .and. (index(noms(iq),"n2") .eq. 0)) then
27            ! Let's get the index of our tracers (we look for igcm _generic_gas and igcm_generic_ice)
28
29            igcm_generic_gas = 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_gas_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_gas_generic = .true.
70        else
71
72            call_ice_gas_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.