source: trunk/LMDZ.MARS/libf/phymars/dyn1d/read_profile_mod.F90 @ 2997

Last change on this file since 2997 was 2991, checked in by jbclement, 18 months ago

Minor fix of initialization of tracers indices + rework of atm_wat_profile. JBC

File size: 8.5 KB
RevLine 
[2355]1module read_profile_mod
2
3implicit none
4
5contains
6!=====================================================================================================================!
7! SUBROUTINE: read_profile ===========================================================================================!
8!=====================================================================================================================!
9! Author: Christophe Mathe
10! Date: 09/06/2020
11!---------------------------------------------------------------------------------------------------------------------!
12! Subject:
13!---------
14!   Read input profile of tracers listed in traceur.def
15!---------------------------------------------------------------------------------------------------------------------!
16! Comments:
17!----------
18!     - If no input profile found, q(traceur) and qsurf(traceur) set to 0, except:
19!          - q(co2) = 0.95
20!          - q(hdo) = q(h2o)*2*155.76e-6*5
21!
22!     - igcm is not available at this part of the code
23!
24!     - To ensure that major isotopologue input profile is read before compute minor isotopologue profile, the
25!       calculation is performed at the end of the subroutine
26!---------------------------------------------------------------------------------------------------------------------!
27! Algorithm:
28!-----------
29!   1. Initialization of q and qsurf to 0
30!   2. Get indices of some tracers
31!   3. Main
32!     3.1. Try to open the input profile
33!       3.1.1. Succeed
34!       3.1.2. Fail
35!         3.1.2.a. Some cases require a special initialization
36!   4. Traitment for minor isotopologue if the major isotopologue input profile does not exist
37!=====================================================================================================================!
38  subroutine read_profile(nb_tracer, nb_layer, qsurf, q)
39
40  use infotrac, only: tname
[2991]41  use tracer_mod, only: igcm_co2, igcm_h2o_vap, igcm_h2o_ice, &
42                        igcm_dust_number, igcm_dust_mass,     &
43                        igcm_ccn_number, igcm_ccn_mass
[2355]44  implicit none
45!---------------------------------------------------------------------------------------------------------------------!
46! VARIABLE DECLARATION
47!---------------------------------------------------------------------------------------------------------------------!
48!  Input arguments:
49!------------------
50  integer, intent(in) :: &
51     nb_layer, & ! number of layer
52     nb_tracer  ! number of traceur read from traceur.def
53!---------------------------------------------------------------------------------------------------------------------!
54!  Output arguments:
55!-------------------
56  real, intent(out) :: &
57     qsurf(nb_tracer), &    ! kg/m2
58     q(nb_layer, nb_tracer) ! kg/kg of atmosphere
59!---------------------------------------------------------------------------------------------------------------------!
60!  Local:
61!--------
62  integer :: &
63     iq,             & ! loop over nb_tracer
64     ilayer,         & ! loop over nb_layer
65     ierr,           & ! open file iostat
66     indice_h2o_vap, & ! indice of h2o_vap tracer
67     indice_h2o_ice, & ! indice of h2o_ice tracer
68     indice_hdo_vap, & ! indice of hdo_vap tracer
69     indice_hdo_ice    ! indice of hdo_ice tracer
70
71  character(len=80), dimension(nb_tracer) :: &
72     name_tracer ! array of all tracers already read in traceur.def
73
74  logical :: &
75     hdo_vap = .false., & ! used to compute hdo_vap profile if its input profile is missing (= .true)
76     hdo_ice = .false.    ! used to compute hdo_ice profile if its input profile is missing (= .true)
77!=====================================================================================================================!
78!=== BEGIN                                                                                                            !
79!=====================================================================================================================!
80! 1. Initialization of q and qsurf to 0
81!---------------------------------------------------------------------------------------------------------------------!
82  q(1:nb_layer,1:nb_tracer) = 0.
83  qsurf(1:nb_tracer) = 0.
84  name_tracer(1:nb_tracer) = ""
85!---------------------------------------------------------------------------------------------------------------------!
86! 2. Get indices of some tracers: igcm is not available at this part of the code
87!---------------------------------------------------------------------------------------------------------------------!
88  do iq = 1, nb_tracer
89    write(name_tracer(iq),"(a)") tname(iq)
90    if (trim(name_tracer(iq)) == 'h2o_vap') then
91      indice_h2o_vap = iq
92    else if (trim(name_tracer(iq)) == 'h2o_ice') then
93      indice_h2o_ice = iq
94    else if (trim(name_tracer(iq)) == 'hdo_vap') then
95      indice_hdo_vap = iq
96    else if (trim(name_tracer(iq)) == 'hdo_ice') then
97      indice_hdo_ice = iq
98    end if
99  end do
100!---------------------------------------------------------------------------------------------------------------------!
101! 3. Main
102!---------------------------------------------------------------------------------------------------------------------!
103  do iq = 1, nb_tracer
104    write(*,*)"  tracer:",trim(name_tracer(iq))
105!---------------------------------------------------------------------------------------------------------------------!
106! 3.1. Try to open the input profile
107!---------------------------------------------------------------------------------------------------------------------!
108    open(91,file='profile_'//trim(name_tracer(iq)), status='old', form='formatted', iostat=ierr)
109!---------------------------------------------------------------------------------------------------------------------!
110! 3.1.1. Succeed: the input file exists
111!---------------------------------------------------------------------------------------------------------------------!
112    if (ierr.eq.0) then
113      read(91,*)qsurf(iq)
[2991]114      if (trim(name_tracer(iq)) == 'co2') igcm_co2 = iq
115      if (trim(name_tracer(iq)) == 'h2o_vap') igcm_h2o_vap = iq
116      if (trim(name_tracer(iq)) == 'h2o_ice') igcm_h2o_ice = iq
117      if (trim(name_tracer(iq)) == 'dust_number') igcm_dust_number = iq
118      if (trim(name_tracer(iq)) == 'dust_mass') igcm_dust_mass = iq
119      if (trim(name_tracer(iq)) == 'ccn_number') igcm_ccn_number = iq
120      if (trim(name_tracer(iq)) == 'ccn_mass') igcm_ccn_mass = iq
[2355]121      do ilayer = 1, nb_layer
122        read(91,*)q(ilayer,iq)
123      end do
124!---------------------------------------------------------------------------------------------------------------------!
125! 3.1.2. Fail: the input file does not exist
126!---------------------------------------------------------------------------------------------------------------------!
127    else
128      write(*,*)"File profile_"//trim(name_tracer(iq))//" not found!"
129!---------------------------------------------------------------------------------------------------------------------!
130! 3.1.2.a. Some cases require a special initialization
131!---------------------------------------------------------------------------------------------------------------------!
132      select case(trim(name_tracer(iq)))
133        case("co2")
134          q(1:nb_layer, iq) = 0.95
135
136        case("hdo_vap")
137          hdo_vap = .true.
138
139        case("hdo_ice")
140          hdo_ice = .true.
141
142       end select
143    end if
144    close(91)
145  end do ! of do iq = 1, nq
146!---------------------------------------------------------------------------------------------------------------------!
147! 4. Traitment for minor isotopologue if the major isotopologue input profile does not exist. At this part, ensure that
148!    the main isotopologue profile is already read before minors isotopologues
149!---------------------------------------------------------------------------------------------------------------------!
150  if (hdo_vap.eqv..true. .and. indice_h2o_vap.ne.0) then
151    do ilayer = 1, nb_layer
152      q(ilayer, indice_hdo_vap) = q(ilayer, indice_h2o_vap)*2*155.76e-6*5
153    end do
154  end if
155
156  if (hdo_ice.eqv..true. .and. indice_h2o_ice.ne.0) then
157    qsurf(indice_hdo_ice) = qsurf(indice_h2o_ice) * 2*155.76e-6*5
158    do ilayer = 1, nb_layer
159      q(ilayer, indice_hdo_ice) = q(ilayer, indice_h2o_ice) * 2*155.76e-6*5
160    end do
161  end if
162!=====================================================================================================================!
163!=== END                                                                                                              !
164!=====================================================================================================================!
165  end subroutine read_profile
166end module read_profile_mod
167
Note: See TracBrowser for help on using the repository browser.