source: LMDZ6/trunk/libf/phylmd/tracco2i_mod.F90 @ 3523

Last change on this file since 3523 was 3453, checked in by oboucher, 6 years ago

Adding some diagnostics for type_trac=co2i

File size: 8.6 KB
RevLine 
[3361]1MODULE tracco2i_mod
2!
3! This module does the work for the interactive CO2 tracers
4!
5CONTAINS
6
7  SUBROUTINE tracco2i(pdtphys, debutphy, &
8       xlat, xlon, pphis, pphi, &
9       t_seri, pplay, paprs, tr_seri, source)
10
11    USE dimphy
[3435]12    USE infotrac_phy
[3450]13    USE geometry_mod, ONLY: cell_area
[3453]14    USE carbon_cycle_mod, ONLY: id_CO2, nbcf_in, fields_in, cfname_in
15    USE carbon_cycle_mod, ONLY: fco2_ocn_day, fco2_ff, fco2_bb, fco2_land, fco2_ocean
[3450]16    USE carbon_cycle_mod, ONLY: carbon_cycle_tr, carbon_cycle_rad, RCO2_glo, RCO2_tot
[3361]17    USE mod_grid_phy_lmdz
[3450]18    USE mod_phys_lmdz_mpi_data, ONLY: is_mpi_root
[3361]19    USE mod_phys_lmdz_para, ONLY: gather, bcast, scatter
20    USE phys_cal_mod
[3450]21    USE phys_state_var_mod, ONLY: pctsrf
22    USE indice_sol_mod, ONLY: nbsrf, is_ter, is_lic, is_oce, is_sic
[3361]23
24    IMPLICIT NONE
25
26    INCLUDE "clesphys.h"
27    INCLUDE "YOMCST.h"
28
29! Input argument
30!---------------
31    REAL,INTENT(IN)    :: pdtphys    ! Pas d'integration pour la physique (seconde)
32    LOGICAL,INTENT(IN) :: debutphy   ! le flag de l'initialisation de la physique
33
34    REAL,DIMENSION(klon),INTENT(IN)        :: xlat    ! latitudes pour chaque point
35    REAL,DIMENSION(klon),INTENT(IN)        :: xlon    ! longitudes pour chaque point
36    REAL,DIMENSION(klon),INTENT(IN)        :: pphis   ! geopotentiel du sol
37    REAL,DIMENSION(klon,klev),INTENT(IN)   :: pphi    ! geopotentiel de chaque couche
38
39    REAL,DIMENSION(klon,klev),INTENT(IN)   :: t_seri  ! Temperature
40    REAL,DIMENSION(klon,klev),INTENT(IN)   :: pplay   ! pression pour le mileu de chaque couche (en Pa)
41    REAL,DIMENSION(klon,klev+1),INTENT(IN) :: paprs   ! pression pour chaque inter-couche (en Pa)
42    REAL,DIMENSION(klon,nbtr),INTENT(INOUT):: source  ! flux de traceur [U/m2/s]
43
44! Output argument
45!----------------
46    REAL,DIMENSION(klon,klev,nbtr),INTENT(INOUT)  :: tr_seri ! Concentration Traceur [U/kgA] 
47
48! Local variables
49!----------------
50
[3421]51    INTEGER                        :: it, k, i, nb
52    REAL, DIMENSION(klon,klev)     :: m_air     ! mass of air in every grid box [kg]
[3361]53    REAL, DIMENSION(klon_glo,klev) :: co2_glo   ! variable temporaire sur la grille global
54    REAL, DIMENSION(klon_glo,klev) :: m_air_glo ! variable temporaire sur la grille global
55
[3450]56    INTEGER, SAVE :: mth_pre=0, day_pre=0
57!$OMP THREADPRIVATE(mth_pre, day_pre)
[3361]58
59    IF (is_mpi_root) THEN
60      PRINT *,'in tracco2i: date from phys_cal_mod =',year_cur,'-',mth_cur,'-',day_cur,'-',hour
61    ENDIF
62
63!--initialisation of CO2 field if not in restart file
64!--dirty way of doing, do it better later
65!--convert 280 ppm into kg CO2 / kg air
66    IF (debutphy) THEN
67      IF (MAXVAL(tr_seri(:,:,id_CO2)).LT.1.e-15) THEN
[3450]68        !!tr_seri(:,:,id_CO2)=280.e-6/RMD*RMCO2
69        tr_seri(:,:,id_CO2)=400.e-6/RMD*RMCO2 !--initialised to 400 ppm for a test
[3361]70      ENDIF
71    ENDIF
72
[3421]73!--calculate mass of air in every grid box in kg air
[3361]74    DO i=1, klon
75    DO k=1, klev
76      m_air(i,k)=(paprs(i,k)-paprs(i,k+1))/RG*cell_area(i)
77    ENDDO
78    ENDDO
79
80!--call CO2 emission routine
81!--co2bb is zero for now
[3421]82!--unit kg CO2 m-2 s-1
83    CALL co2_emissions(debutphy)
[3361]84
[3421]85!--retrieving land and ocean CO2 flux
[3453]86    fco2_land(:)=0.0
87    fco2_ocean(:)=0.0
[3421]88    DO nb=1, nbcf_in
[3450]89!--fCO2_nep comes in unit of kg C m-2 s-1
90!--converting to kg CO2 m-2 s-1
[3453]91      IF (cfname_in(nb) == "fCO2_nep" )   fco2_land(:)=fields_in(:,nb)*RMCO2/RMC*pctsrf(:,is_ter)
[3450]92!--fCO2_fgco2 comes in unit of mol C02 m-2 s-1
93!--converting to kg CO2 m-2 s-1 + change sign
[3453]94      IF (cfname_in(nb) == "fCO2_fgco2" ) fco2_ocean(:)=-1.*fco2_ocn_day(:)*RMCO2/1.e3*(pctsrf(:,is_oce)+pctsrf(:,is_sic))
[3421]95    ENDDO
96
[3453]97!--preparing the net anthropogenic flux at the surface for mixing layer
[3361]98!--unit kg CO2 / m2 / s
[3453]99    source(:,id_CO2)=fco2_ff(:)+fco2_bb(:)+fco2_land(:)+fco2_ocean(:)
[3361]100
101!--computing global mean CO2 for radiation
[3450]102!--every timestep for now but enough every day
[3361]103!    IF (debutphy.OR.mth_cur.NE.mth_pre) THEN
[3450]104!    IF (debutphy.OR.day_cur.NE.day_pre) THEN
[3361]105      CALL gather(tr_seri(:,:,id_CO2),co2_glo)
106      CALL gather(m_air,m_air_glo)
107!$OMP MASTER
[3450]108!--compute a global mean CO2 value and print its value in ppm
[3361]109       IF (is_mpi_root) THEN
[3450]110         RCO2_tot=SUM(co2_glo*m_air_glo)  !--unit kg CO2
111         RCO2_glo=RCO2_tot/SUM(m_air_glo) !--unit kg CO2 / kg air
112         PRINT *,'tracco2i: global CO2 in ppm =', RCO2_glo*1.e6*RMD/RMCO2
113         PRINT *,'tracco2i: total CO2 in kg =', RCO2_tot
[3361]114       ENDIF
115!$OMP END MASTER
116       CALL bcast(RCO2_glo)
117       mth_pre=mth_cur
[3450]118       day_pre=day_cur
119!--if not carbon_cycle_tr, then we reinitialize the CO2 each day to its global mean value
120       IF (.NOT.carbon_cycle_tr) THEN
121         tr_seri(:,:,id_CO2)=RCO2_glo
122       ENDIF
[3361]123!    ENDIF
124
125  END SUBROUTINE tracco2i
126
[3421]127  SUBROUTINE co2_emissions(debutphy)
[3361]128
129    USE dimphy
[3435]130    USE infotrac_phy
[3361]131    USE geometry_mod, ONLY : cell_area
132    USE mod_grid_phy_lmdz
[3366]133    USE mod_phys_lmdz_mpi_data, ONLY: is_mpi_root
[3361]134    USE mod_phys_lmdz_para, ONLY: gather, scatter
135    USE phys_cal_mod
136
137    USE netcdf95, ONLY: nf95_close, nf95_gw_var, nf95_inq_varid, nf95_open
138    USE netcdf, ONLY: nf90_get_var, nf90_noerr, nf90_nowrite
139
[3421]140    USE carbon_cycle_mod, ONLY : fco2_ff, fco2_bb
141
[3361]142    IMPLICIT NONE
143
[3368]144    INCLUDE "YOMCST.h"
[3361]145    LOGICAL,INTENT(IN) :: debutphy
146
147! For NetCDF:
[3366]148    INTEGER ncid_in  ! IDs for input files
149    INTEGER varid, ncerr
[3361]150
151    INTEGER :: n_glo, n_month
152    REAL, POINTER:: vector(:), time(:)
[3366]153    REAL,ALLOCATABLE       :: flx_co2ff_glo(:,:) !  fossil-fuel CO2
[3421]154    REAL,ALLOCATABLE       :: flx_co2bb_glo(:,:) !  biomass-burning CO2
[3366]155    REAL,ALLOCATABLE, SAVE :: flx_co2ff(:,:)     !  fossil-fuel CO2
156    REAL,ALLOCATABLE, SAVE :: flx_co2bb(:,:)     !  biomass-burning CO2
[3421]157!$OMP THREADPRIVATE(flx_co2ff,flx_co2bb)
[3361]158
[3383]159!! may be controlled via the .def later on
160!! also co2bb for now comes from ORCHIDEE
161    LOGICAL, PARAMETER :: readco2ff=.TRUE., readco2bb=.FALSE.
[3361]162
[3366]163    IF (debutphy) THEN
[3361]164
[3366]165    ALLOCATE(flx_co2ff(klon,12))
166    ALLOCATE(flx_co2bb(klon,12))
167
[3361]168!$OMP MASTER
[3366]169    IF (is_mpi_root) THEN
170   
[3383]171      IF (.NOT.ALLOCATED(flx_co2ff_glo)) ALLOCATE(flx_co2ff_glo(klon_glo,12))
172      IF (.NOT.ALLOCATED(flx_co2bb_glo)) ALLOCATE(flx_co2bb_glo(klon_glo,12))
[3361]173
[3366]174!--reading CO2 fossil fuel emissions
175      IF (readco2ff) THEN
[3361]176
[3366]177        ! ... Open the COZff file
178        CALL nf95_open("sflx_lmdz_co2_ff.nc", nf90_nowrite, ncid_in)
[3361]179
[3366]180        CALL nf95_inq_varid(ncid_in, "vector", varid)
181        CALL nf95_gw_var(ncid_in, varid, vector)
182        n_glo = size(vector)
183        IF (n_glo.NE.klon_glo) THEN
184           PRINT *,'sflx_lmdz_co2_ff: le nombre de points n est pas egal a klon_glo'
185           STOP
186        ENDIF
[3361]187
[3366]188        CALL nf95_inq_varid(ncid_in, "time", varid)
189        CALL nf95_gw_var(ncid_in, varid, time)
190        n_month = size(time)
191        IF (n_month.NE.12) THEN
192           PRINT *,'sflx_lmdz_co2_ff: le nombre de month n est pas egal a 12'
193           STOP
194        ENDIF
[3361]195
[3366]196!--reading flx_co2 for fossil fuel
197        CALL nf95_inq_varid(ncid_in, "flx_co2", varid)
198        ncerr = nf90_get_var(ncid_in, varid, flx_co2ff_glo)
199
200        CALL nf95_close(ncid_in)
201   
202      ELSE  !--co2ff not to be read
203        flx_co2ff_glo(:,:)=0.0
204      ENDIF
205
206!--reading CO2 biomass burning emissions
207      IF (readco2bb) THEN
208
209      ! ... Open the CO2bb file
210      CALL nf95_open("sflx_lmdz_co2_bb.nc", nf90_nowrite, ncid_in)
211
212      CALL nf95_inq_varid(ncid_in, "vector", varid)
213      CALL nf95_gw_var(ncid_in, varid, vector)
214      n_glo = size(vector)
215      IF (n_glo.NE.klon_glo) THEN
216         PRINT *,'sflx_lmdz_co2_bb: le nombre de points n est pas egal a klon_glo'
217         STOP
218      ENDIF
219
220      CALL nf95_inq_varid(ncid_in, "time", varid)
221      CALL nf95_gw_var(ncid_in, varid, time)
222      n_month = size(time)
223      IF (n_month.NE.12) THEN
224         PRINT *,'sflx_lmdz_co2_bb: le nombre de month n est pas egal a 12'
225         STOP
226      ENDIF
227
228!--reading flx_co2 for biomass burning
229      CALL nf95_inq_varid(ncid_in, "flx_co2", varid)
230      ncerr = nf90_get_var(ncid_in, varid, flx_co2bb_glo)
231
232      CALL nf95_close(ncid_in)
233   
234      ELSE  !--co2bb not to be read
235        flx_co2bb_glo(:,:)=0.0
236      ENDIF
237
238    ENDIF
[3361]239!$OMP END MASTER
240
241!--scatter on all proc
[3366]242    CALL scatter(flx_co2ff_glo,flx_co2ff)
243    CALL scatter(flx_co2bb_glo,flx_co2bb)
[3361]244
245!$OMP MASTER
[3366]246    IF (is_mpi_root) THEN
247       DEALLOCATE(flx_co2ff_glo)
248       DEALLOCATE(flx_co2bb_glo)
249    ENDIF
[3361]250!$OMP END MASTER
251
252  ENDIF !--end debuthy
253
254!---select the correct month
255  IF (mth_cur.LT.1.OR.mth_cur.GT.12) THEN
[3366]256    PRINT *,'probleme avec le mois dans co2_ini =', mth_cur
[3361]257  ENDIF
[3421]258  IF (.NOT.ALLOCATED(fco2_ff)) ALLOCATE(fco2_ff(klon))
259  IF (.NOT.ALLOCATED(fco2_bb)) ALLOCATE(fco2_bb(klon))
260  fco2_ff(:) = flx_co2ff(:,mth_cur)
261  fco2_bb(:) = flx_co2bb(:,mth_cur)
[3361]262
263  END SUBROUTINE co2_emissions
264
265END MODULE tracco2i_mod
Note: See TracBrowser for help on using the repository browser.