source: trunk/LMDZ.COMMON/libf/evolution/read_data_PCM_mod.F90 @ 3563

Last change on this file since 3563 was 3367, checked in by jbclement, 6 months ago

PEM:
Removal of useless condition and variable + some updates.
JBC

File size: 15.2 KB
RevLine 
[3096]1MODULE read_data_PCM_mod
[2794]2
[3076]3use netcdf, only: nf90_open, NF90_NOWRITE, nf90_noerr, nf90_strerror, &
4                  nf90_get_var, nf90_inq_varid, nf90_inq_dimid,       &
5                  nf90_inquire_dimension, nf90_close
[2794]6
[3076]7implicit none
8
[3122]9character(13), parameter :: modname = 'read_data_PCM'
10character(256)           :: msg      ! for reading
11integer                  :: fID, vID ! for reading
[3076]12
[2794]13!=======================================================================
[3076]14contains
15!=======================================================================
16
[3367]17SUBROUTINE read_data_PCM(filename,timelen,iim_input,jjm_input,ngrid,nslope,vmr_co2_PCM_phys,ps_timeseries, &
18                         min_co2_ice,min_h2o_ice,tsurf_avg,tsoil_avg,tsurf_PCM,tsoil_PCM,q_co2,q_h2o,      &
[3149]19                         watersurf_density_avg,watersoil_density)
[3076]20use comsoil_h,             only: nsoilmx
21use comsoil_h_PEM,         only: soil_pem
22use constants_marspem_mod, only: m_co2, m_noco2
23
24implicit none
25
26!=======================================================================
[2794]27!
[3096]28! Purpose: Read initial confitions file from the PCM
[2794]29!
[2855]30! Authors: RV & LL
[2794]31!=======================================================================
32
[3076]33include "dimensions.h"
[2794]34
[3076]35!=======================================================================
[2794]36! Arguments:
[3149]37character(*), intent(in) :: filename                            ! File name
38integer,      intent(in) :: timelen                             ! number of times stored in the file
39integer,      intent(in) :: iim_input, jjm_input, ngrid, nslope ! number of points in the lat x lon dynamical grid, number of subgrid slopes
[2855]40! Ouputs
[3181]41real, dimension(ngrid,nslope),         intent(out) :: min_co2_ice      ! Minimum of co2 ice per slope of the year [kg/m^2]
[3076]42real, dimension(ngrid,nslope),         intent(out) :: min_h2o_ice      ! Minimum of h2o ice per slope of the year [kg/m^2]
[3181]43real, dimension(ngrid,timelen),        intent(out) :: vmr_co2_PCM_phys ! Physics x Times co2 volume mixing ratio retrieve from the PCM [m^3/m^3]
[3076]44real, dimension(ngrid,timelen),        intent(out) :: ps_timeseries    ! Surface Pressure [Pa]
45real, dimension(ngrid,timelen),        intent(out) :: q_co2            ! CO2 mass mixing ratio in the first layer [kg/m^3]
46real, dimension(ngrid,timelen),        intent(out) :: q_h2o            ! H2O mass mixing ratio in the first layer [kg/m^3]
[2794]47!SOIL
[3149]48real, dimension(ngrid,nslope),                 intent(out) :: tsurf_avg             ! Average surface temperature of the concatenated file [K]
49real, dimension(ngrid,nsoilmx,nslope),         intent(out) :: tsoil_avg             ! Average soil temperature of the concatenated file [K]
50real, dimension(ngrid,nslope,timelen),         intent(out) :: tsurf_PCM             ! Surface temperature of the concatenated file, time series [K]
51real, dimension(ngrid,nsoilmx,nslope,timelen), intent(out) :: tsoil_PCM             ! Soil temperature of the concatenated file, time series [K]
52real, dimension(ngrid,nslope),                 intent(out) :: watersurf_density_avg ! Water density at the surface [kg/m^3]
[3076]53real, dimension(ngrid,nsoilmx,nslope,timelen), intent(out) :: watersoil_density     ! Water density in the soil layer, time series [kg/m^3]
[3122]54!=======================================================================
55! Local Variables
[3210]56integer                         :: i, j, l, t                                                    ! loop variables
[3206]57real                            :: A, B, mmean                                                   ! Molar Mass of co2 and no co2, A;B intermediate variables to compute the mean molar mass of the layer
58integer                         :: islope                                                        ! loop for variables
59character(2)                    :: num                                                           ! for reading sloped variables
[3123]60real, dimension(iim_input + 1,jjm_input + 1,nslope,timelen)         :: h2o_ice_s_dyn             ! h2o ice per slope of the concatenated file [kg/m^2]
[3122]61real, dimension(iim_input + 1,jjm_input + 1,nslope,timelen)         :: watercap
[3130]62real, dimension(iim_input + 1,jjm_input + 1,nslope,timelen)         :: perennial_co2ice
[3149]63real, dimension(iim_input + 1,jjm_input + 1,timelen)                :: vmr_co2_PCM               ! CO2 volume mixing ratio in the first layer [mol/m^3]
64real, dimension(iim_input + 1,jjm_input + 1,timelen)                :: ps_PCM                    ! Surface Pressure [Pa]
[3076]65real, dimension(iim_input + 1,jjm_input + 1,nslope)                 :: min_co2_ice_dyn
66real, dimension(iim_input + 1,jjm_input + 1,nslope)                 :: min_h2o_ice_dyn
[3149]67real, dimension(iim_input + 1,jjm_input + 1,nslope)                 :: tsurf_avg_dyn             ! Average surface temperature of the concatenated file [K]
68real, dimension(iim_input + 1,jjm_input + 1,nsoilmx,nslope)         :: tsoil_avg_dyn             ! Average soil temperature of the concatenated file [K]
69real, dimension(iim_input + 1,jjm_input + 1,nslope,timelen)         :: tsurf_PCM_dyn             ! Surface temperature of the concatenated file, time series [K]
70real, dimension(iim_input + 1,jjm_input + 1,nsoilmx,nslope,timelen) :: tsoil_PCM_dyn             ! Soil temperature of the concatenated file, time series [K]
[3123]71real, dimension(iim_input + 1,jjm_input + 1,timelen)                :: q_co2_dyn                 ! CO2 mass mixing ratio in the first layer [kg/m^3]
72real, dimension(iim_input + 1,jjm_input + 1,timelen)                :: q_h2o_dyn                 ! H2O mass mixing ratio in the first layer [kg/m^3]
[3367]73real, dimension(iim_input + 1,jjm_input + 1,nslope,timelen)         :: co2_ice_slope_dyn         ! co2 ice amount per slope of the year [kg/m^2]
[3123]74real, dimension(iim_input + 1,jjm_input + 1,nslope,timelen)         :: watersurf_density_dyn     ! Water density at the surface, time series [kg/m^3]
[3149]75real, dimension(iim_input + 1,jjm_input + 1,nslope)                 :: watersurf_density_dyn_avg ! Water density at the surface, dynamic grid, yearly averaged [kg/m^3]
[3123]76real, dimension(iim_input + 1,jjm_input + 1,nsoilmx,nslope,timelen) :: watersoil_density_dyn     ! Water density in the soil layer, time series [kg/m^3]
[2794]77!-----------------------------------------------------------------------
[3122]78! Open the NetCDF file
79write(*,*) "Opening "//filename//"..."
80call error_msg(NF90_OPEN(filename,NF90_NOWRITE,fID),"open",filename)
[2794]81
[3363]82! Download the data from the file
[3149]83call get_var3("ps",ps_PCM)
[3122]84write(*,*) "Data for surface pressure downloaded!"
[2794]85
[3122]86call get_var3("co2_layer1",q_co2_dyn)
87write(*,*) "Data for vmr co2 downloaded!"
[2794]88
[3122]89call get_var3("h2o_layer1",q_h2o_dyn)
90write(*,*) "Data for vmr h2o downloaded!"
[2794]91
[3122]92if (nslope == 1) then ! There is no slope
93    call get_var3("co2ice",co2_ice_slope_dyn(:,:,1,:))
94    write(*,*) "Data for co2_ice downloaded!"
[2835]95
[3122]96    call get_var3("h2o_ice_s",h2o_ice_s_dyn(:,:,1,:))
97    write(*,*) "Data for h2o_ice_s downloaded!"
[2835]98
[3149]99    call get_var3("tsurf",tsurf_PCM_dyn(:,:,1,:))
[3143]100    write(*,*) "Data for tsurf downloaded!"
101
[3122]102#ifndef CPP_STD
103    call get_var3("watercap",watercap(:,:,1,:))
104    write(*,*) "Data for watercap downloaded!"
[3130]105
106    call get_var3("perennial_co2ice",perennial_co2ice(:,:,1,:))
107    write(*,*) "Data for perennial_co2ice downloaded!"
[2835]108
[3122]109    if (soil_pem) then
[3149]110        call get_var4("soiltemp",tsoil_PCM_dyn(:,:,:,1,:))
[3122]111        write(*,*) "Data for soiltemp downloaded!"
[2794]112
[3152]113        call get_var4("waterdensity_soil",watersoil_density_dyn(:,:,:,1,:))
[3122]114        write(*,*) "Data for waterdensity_soil downloaded!"
[2794]115
[3152]116        call get_var3("waterdensity_surface",watersurf_density_dyn(:,:,1,:))
[3122]117        write(*,*) "Data for waterdensity_surface downloaded!"
118    endif !soil_pem
119#endif
120else ! We use slopes
[3076]121    do islope = 1,nslope
122        write(num,'(i2.2)') islope
123        call get_var3("co2ice_slope"//num,co2_ice_slope_dyn(:,:,islope,:))
124    enddo
[3122]125    write(*,*) "Data for co2_ice downloaded!"
126
[3076]127    do islope = 1,nslope
128        write(num,'(i2.2)') islope
129        call get_var3("h2o_ice_s_slope"//num,h2o_ice_s_dyn(:,:,islope,:))
130    enddo
[3122]131    write(*,*) "Data for h2o_ice_s downloaded!"
[2980]132
[3143]133    do islope = 1,nslope
134        write(num,'(i2.2)') islope
[3149]135        call get_var3("tsurf_slope"//num,tsurf_PCM_dyn(:,:,islope,:))
[3143]136    enddo
137    write(*,*) "Data for tsurf downloaded!"
138
[2985]139#ifndef CPP_STD
[3076]140    do islope = 1,nslope
141        write(num,'(i2.2)') islope
[3122]142        call get_var3("watercap_slope"//num,watercap(:,:,islope,:))
[3076]143    enddo
[3122]144    write(*,*) "Data for watercap downloaded!"
[3130]145
146    do islope = 1,nslope
147        write(num,'(i2.2)') islope
148        call get_var3("perennial_co2ice_slope"//num,perennial_co2ice(:,:,islope,:))
149    enddo
150    write(*,*) "Data for perennial_co2ice downloaded!"
[2985]151
[3076]152    if (soil_pem) then
153        do islope = 1,nslope
154            write(num,'(i2.2)') islope
[3149]155            call get_var4("soiltemp_slope"//num,tsoil_PCM_dyn(:,:,:,islope,:))
[3076]156        enddo
[3122]157        write(*,*) "Data for soiltemp downloaded!"
[2835]158
[3076]159        do islope = 1,nslope
160            write(num,'(i2.2)') islope
[3122]161            call get_var4("waterdensity_soil_slope"//num,watersoil_density_dyn(:,:,:,islope,:))
[3076]162        enddo
[3122]163        write(*,*) "Data for waterdensity_soil downloaded!"
[2794]164
[3076]165        do islope = 1,nslope
166            write(num,'(i2.2)') islope
[3122]167            call get_var3("waterdensity_surface"//num,watersurf_density_dyn(:,:,islope,:))
[3076]168        enddo
[3122]169        write(*,*) "Data for waterdensity_surface downloaded!"
[3076]170    endif !soil_pem
[2985]171#endif
[3122]172endif
[2897]173
[3363]174! Close the NetCDF file
175write(*,*) "Closing "//filename//"..."
176call error_msg(nf90_close(fID),"close",filename)
177
[2794]178! Compute the minimum over the year for each point
[3122]179write(*,*) "Computing the min of h2o_ice..."
180where (h2o_ice_s_dyn < 0.) h2o_ice_s_dyn = 0.
[3149]181min_h2o_ice_dyn = minval(h2o_ice_s_dyn + watercap,4)
[3122]182write(*,*) "Computing the min of co2_ice..."
183where (co2_ice_slope_dyn < 0.) co2_ice_slope_dyn = 0.
[3149]184min_co2_ice_dyn = minval(co2_ice_slope_dyn + perennial_co2ice,4)
[2794]185
[3122]186! Compute averages over the year for each point
187write(*,*) "Computing the average of tsurf..."
[3149]188tsurf_avg_dyn = sum(tsurf_PCM_dyn,4)/timelen
[2794]189
[3076]190if (soil_pem) then
[3122]191    write(*,*) "Computing average of tsoil..."
[3149]192    tsoil_avg_dyn = sum(tsoil_PCM_dyn,5)/timelen
[3122]193    write(*,*) "Computing average of waterdensity_surface..."
[3149]194    watersurf_density_dyn_avg = sum(watersurf_density_dyn,4)/timelen
[3076]195endif
[2794]196
[3122]197! By definition, we get rid of the negative values
198A = (1./m_co2 - 1./m_noco2)
199B = 1./m_noco2
[3076]200do i = 1,iim + 1
201    do j = 1,jjm_input + 1
202        do t = 1, timelen
203            if (q_co2_dyn(i,j,t) < 0) then
204                q_co2_dyn(i,j,t) = 1.e-10
205            else if (q_co2_dyn(i,j,t) > 1) then
206                q_co2_dyn(i,j,t) = 1.
207            endif
208            if (q_h2o_dyn(i,j,t) < 0) then
[3130]209                q_h2o_dyn(i,j,t) = 1.e-10
[3076]210            else if (q_h2o_dyn(i,j,t) > 1) then
211                q_h2o_dyn(i,j,t) = 1.
212            endif
213            mmean = 1/(A*q_co2_dyn(i,j,t) + B)
[3149]214            vmr_co2_PCM(i,j,t) = q_co2_dyn(i,j,t)*mmean/m_co2
[3076]215        enddo
216    enddo
217enddo
[2794]218
[2980]219#ifndef CPP_1D
[3149]220    call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,vmr_co2_PCM,vmr_co2_PCM_phys)
221    call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,ps_PCM,ps_timeseries)
[3076]222    call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,q_co2_dyn,q_co2)
223    call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,q_h2o_dyn,q_h2o)
[2980]224
[3076]225    do islope = 1,nslope
226        call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,min_co2_ice_dyn(:,:,islope),min_co2_ice(:,islope))
227        call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,min_h2o_ice_dyn(:,:,islope),min_h2o_ice(:,islope))
228        if (soil_pem) then
229            do l = 1,nsoilmx
[3149]230                call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,tsoil_avg_dyn(:,:,l,islope),tsoil_avg(:,l,islope))
[3122]231                do t = 1,timelen
[3149]232                    call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,tsoil_PCM_dyn(:,:,l,islope,t),tsoil_PCM(:,l,islope,t))
[3076]233                    call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,watersoil_density_dyn(:,:,l,islope,t),watersoil_density(:,l,islope,t))
234                enddo
235            enddo
[3149]236            call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,watersurf_density_dyn_avg(:,:,islope),watersurf_density_avg(:,islope))
[3122]237        endif ! soil_pem
[3076]238        do t = 1,timelen
[3149]239            call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,tsurf_PCM_dyn(:,:,islope,t),tsurf_PCM(:,islope,t))
[3076]240        enddo
241    enddo
[2849]242
[3149]243    call gr_dyn_fi(nslope,iim_input + 1,jjm_input + 1,ngrid,tsurf_avg_dyn,tsurf_avg)
[3076]244#else
[3149]245    vmr_co2_PCM_phys(1,:) = vmr_co2_PCM(1,1,:)
246    ps_timeseries(1,:) = ps_PCM(1,1,:)
[3076]247    q_co2(1,:) = q_co2_dyn(1,1,:)
248    q_h2o(1,:) = q_h2o_dyn(1,1,:)
249    min_co2_ice(1,:) = min_co2_ice_dyn(1,1,:)
250    min_h2o_ice(1,:) = min_h2o_ice_dyn(1,1,:)
251    if (soil_pem) then
[3149]252        tsoil_avg(1,:,:) = tsoil_avg_dyn(1,1,:,:)
253        tsoil_PCM(1,:,:,:) = tsoil_PCM_dyn(1,1,:,:,:)
[3076]254        watersoil_density(1,:,:,:) = watersoil_density_dyn(1,1,:,:,:)
[3149]255        watersurf_density_avg(1,:) = watersurf_density_dyn_avg(1,1,:)
[3122]256    endif ! soil_pem
[3149]257    tsurf_PCM(1,:,:) = tsurf_PCM_dyn(1,1,:,:)
258    tsurf_avg(1,:) = tsurf_avg_dyn(1,1,:)
[3076]259#endif
[2897]260
[3096]261END SUBROUTINE read_data_PCM
[2897]262
[3122]263!=======================================================================
264
[3076]265SUBROUTINE check_dim(n1,n2,str1,str2)
[2980]266
[3076]267implicit none
[2980]268
[3076]269integer,            intent(in) :: n1, n2
270character(len = *), intent(in) :: str1, str2
[2980]271
[3076]272character(256) :: s1, s2
[2794]273
[3076]274if (n1 /= n2) then
275    s1 = 'value of '//trim(str1)//' ='
276    s2 = ' read in starting file differs from parametrized '//trim(str2)//' ='
277    write(msg,'(10x,a,i4,2x,a,i4)')trim(s1),n1,trim(s2),n2
278    call abort_gcm(trim(modname),trim(msg),1)
279endif
280
[2794]281END SUBROUTINE check_dim
282
[3076]283!=======================================================================
[2794]284
[3130]285SUBROUTINE get_var1(var,v)
[3076]286
287implicit none
288
[3130]289character(len = *), intent(in)  :: var
[3076]290real, dimension(:), intent(out) :: v
291
[3130]292call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var)
293call error_msg(NF90_GET_VAR(fID,vID,v),"get",var)
[3076]294
[2794]295END SUBROUTINE get_var1
296
[3076]297!=======================================================================
[2794]298
[3130]299SUBROUTINE get_var3(var,v) ! on U grid
[2794]300
[3076]301implicit none
302
[3130]303character(len = *),     intent(in)  :: var
[3076]304real, dimension(:,:,:), intent(out) :: v
305
[3130]306call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var)
307call error_msg(NF90_GET_VAR(fID,vID,v),"get",var)
[3076]308
[2794]309END SUBROUTINE get_var3
310
[3076]311!=======================================================================
312
[3130]313SUBROUTINE get_var4(var,v)
[3076]314
315implicit none
316
[3130]317character(len = *),       intent(in)  :: var
[3076]318real, dimension(:,:,:,:), intent(out) :: v
319
[3130]320call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var)
321call error_msg(NF90_GET_VAR(fID,vID,v),"get",var)
[3076]322
[2794]323END SUBROUTINE get_var4
324
[3076]325!=======================================================================
[2794]326
[3076]327SUBROUTINE error_msg(ierr,typ,nam)
328
329implicit none
330
331integer,            intent(in) :: ierr !--- NetCDF ERROR CODE
332character(len = *), intent(in) :: typ  !--- TYPE OF OPERATION
333character(len = *), intent(in) :: nam  !--- FIELD/FILE NAME
334
335if (ierr == nf90_noerr) return
336select case(typ)
337    case('inq');   msg="Field <"//trim(nam)//"> is missing"
338    case('get');   msg="Reading failed for <"//trim(nam)//">"
[3363]339    case('put');   msg="Writing failed for <"//trim(nam)//">"
[3076]340    case('open');  msg="File opening failed for <"//trim(nam)//">"
341    case('close'); msg="File closing failed for <"//trim(nam)//">"
342    case default
343        write(*,*) 'There is no message for this error.'
344        error stop
345end select
346call abort_gcm(trim(modname),trim(msg),ierr)
347
348END SUBROUTINE error_msg
349
[3096]350END MODULE read_data_PCM_mod
Note: See TracBrowser for help on using the repository browser.