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

Last change on this file since 3966 was 3961, checked in by jbclement, 4 months ago

PEM:
Clearing all the tags and code related to the Generic PCM.
JBC

File size: 14.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
[3571]17SUBROUTINE read_data_PCM(filename1,filename2,timelen,iim_input,jjm_input,ngrid,nslope,vmr_co2_PCM,ps_timeseries,ps_avg,tsurf_avg_yr1,tsurf_avg, &
18                         tsoil_avg,tsoil_timeseries,min_co2_ice,min_h2o_ice,q_co2,q_h2o,watersurf_density_avg,watersoil_density_timeseries)
[3076]19use comsoil_h,             only: nsoilmx
20use comsoil_h_PEM,         only: soil_pem
21use constants_marspem_mod, only: m_co2, m_noco2
22
23implicit none
24
25!=======================================================================
[2794]26!
[3096]27! Purpose: Read initial confitions file from the PCM
[2794]28!
[3602]29! Authors: JBC
[2794]30!=======================================================================
31
[3076]32include "dimensions.h"
[2794]33
[3076]34!=======================================================================
[3571]35! Inputs:
36character(*), intent(in) :: filename1, filename2                ! File names
37integer,      intent(in) :: timelen                             ! Number of times stored in the file
38integer,      intent(in) :: iim_input, jjm_input, ngrid, nslope ! Number of points in the lat x lon dynamical grid, number of subgrid slopes
[2855]39! Ouputs
[3571]40real, dimension(ngrid,nslope,2),               intent(out) :: min_co2_ice                  ! Minimum of co2 ice per slope of the year [kg/m^2]
41real, dimension(ngrid,nslope,2),               intent(out) :: min_h2o_ice                  ! Minimum of h2o ice per slope of the year [kg/m^2]
42real, dimension(ngrid,timelen),                intent(out) :: vmr_co2_PCM                  ! Grid points x Times co2 volume mixing ratio retrieve from the PCM [m^3/m^3]
43real, dimension(ngrid,timelen),                intent(out) :: q_co2                        ! CO2 mass mixing ratio in the first layer [kg/m^3]
44real, dimension(ngrid,timelen),                intent(out) :: q_h2o                        ! H2O mass mixing ratio in the first layer [kg/m^3]
45real, dimension(ngrid,timelen),                intent(out) :: ps_timeseries                ! Surface pressure timeseries [Pa]
46real, dimension(ngrid),                        intent(out) :: ps_avg                       ! Averaged surface pressure [K]
47real, dimension(ngrid,nslope),                 intent(out) :: tsurf_avg                    ! Averaged surface temperature [K]
48real, dimension(ngrid,nslope),                 intent(out) :: tsurf_avg_yr1                ! Averaged surface temperature for year 1 [K]
49real, dimension(ngrid,nsoilmx,nslope),         intent(out) :: tsoil_avg                    ! Averaged soil temperature for year 2 [K]
50real, dimension(ngrid,nsoilmx,nslope,timelen), intent(out) :: tsoil_timeseries             ! Soil temperature timeseries [K]
51real, dimension(ngrid,nslope),                 intent(out) :: watersurf_density_avg        ! Water density at the surface [kg/m^3]
52real, dimension(ngrid,nsoilmx,nslope,timelen), intent(out) :: watersoil_density_timeseries ! Water density timeseries in the soil layer [kg/m^3]
53! Local variables
[3620]54integer                               :: i, j, t, islope                       ! Loop variables
[3609]55real                                  :: A, B                                  ! Intermediate variables to compute the mean molar mass of the layer
56character(:), allocatable             :: num                                   ! For reading sloped variables
[3602]57real, dimension(:,:),     allocatable :: var2_read                             ! Variables for reading (2 dimensions)
58real, dimension(:,:,:),   allocatable :: var3_read_1, var3_read_2, var3_read_3 ! Variables for reading (3 dimensions)
59real, dimension(:,:,:,:), allocatable :: var4_read                             ! Variables for reading (4 dimensions)
[2794]60!-----------------------------------------------------------------------
[3571]61
[3602]62!------------------------------- Year 1 --------------------------------
[3122]63! Open the NetCDF file
[3571]64write(*,*) "Opening "//filename1//"..."
65call error_msg(NF90_OPEN(filename1,NF90_NOWRITE,fID),"open",filename1)
[2794]66
[3598]67! Download the data from the file
[3602]68allocate(var2_read(iim_input + 1,jjm_input + 1),var3_read_1(iim_input + 1,jjm_input + 1,timelen),var3_read_2(iim_input + 1,jjm_input + 1,timelen))
69
[3609]70if (nslope == 1) then ! No slope
71    allocate(character(0) :: num)
72else ! We use slopes
73    allocate(character(8) :: num)
74endif
75
76do islope = 1,nslope
77    if (nslope /= 1) then
[3611]78        num='  '
[3609]79        write(num,'(i2.2)') islope
80        num = '_slope'//num
81    endif
82
[3602]83    ! CO2 ice
84    !--------
[3609]85    call get_var3("co2ice"//num,var3_read_1)
[3602]86    where (var3_read_1 < 0.) var3_read_1 = 0.
[3609]87    write(*,*) "Data for co2_ice"//num//" downloaded."
88    call get_var3("perennial_co2ice"//num,var3_read_2)
89    write(*,*) "Data for perennial_co2ice"//num//" downloaded."
[3571]90
[3602]91    ! Compute the minimum over the year for each point
92    var2_read = minval(var3_read_1 + var3_read_2,3)
93#ifndef CPP_1D
[3609]94    call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_co2_ice(:,islope,1))
[3602]95#else
[3609]96    min_co2_ice(1,islope,1) = var2_read(1,1)
[3602]97#endif
[3609]98    write(*,*) "Min of co2_ice"//num//" computed."
[3602]99
100    ! H2O ice
101    !--------
[3609]102    call get_var3("h2o_ice_s"//num,var3_read_1)
[3602]103    where (var3_read_1 < 0.) var3_read_1 = 0.
[3609]104    write(*,*) "Data for h2o_ice_s"//num//" downloaded."
105    call get_var3("watercap"//num,var3_read_2)
106    write(*,*) "Data for watercap"//num//" downloaded."
[3571]107
[3602]108    ! Compute the minimum over the year for each point
109    var2_read = minval(var3_read_1 + var3_read_2,3)
110#ifndef CPP_1D
[3609]111    call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_h2o_ice(:,islope,1))
[3602]112#else
[3609]113    min_h2o_ice(1,islope,1) = var2_read(1,1)
[3571]114#endif
[3609]115    write(*,*) "Min of h2o_ice"//num//" computed."
[3602]116
117    ! Tsurf
118    !------
[3609]119    call get_var3("tsurf"//num,var3_read_1)
120    write(*,*) "Data for tsurf"//num//" downloaded."
[3602]121
122    ! Compute average over the year for each point
123    var2_read = sum(var3_read_1,3)/timelen
124#ifndef CPP_1D
[3609]125    call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,tsurf_avg_yr1(:,islope))
[3602]126#else
[3609]127    tsurf_avg_yr1(1,islope) = var2_read(1,1)
[3602]128#endif
[3609]129    write(*,*) "Average of tsurf"//num//" computed."
130enddo
[3602]131
132! Close the NetCDF file
133call error_msg(nf90_close(fID),"close",filename1)
[3603]134write(*,*) '> "'//filename1//'" processed!'
[3571]135
[3602]136!------------------------------- Year 2 --------------------------------
[3571]137! Open the NetCDF file
138write(*,*) "Opening "//filename2//"..."
139call error_msg(NF90_OPEN(filename2,NF90_NOWRITE,fID),"open",filename2)
140
[3363]141! Download the data from the file
[3602]142allocate(var3_read_3(iim_input + 1,jjm_input + 1,nsoilmx),var4_read(iim_input + 1,jjm_input + 1,nsoilmx,timelen))
143
144! Surface pressure
145!-----------------
146call get_var3("ps",var3_read_1)
[3591]147write(*,*) "Data for surface pressure downloaded."
[2794]148
[3602]149! Compute average over the year for each point
150var2_read = sum(var3_read_1,3)/timelen
151#ifndef CPP_1D
152    call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,var3_read_1,ps_timeseries)
153    call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,ps_avg)
154#else
[3609]155    ps_timeseries(1,:) = var3_read_1(1,1,:)
[3602]156    ps_avg(1) = var2_read(1,1)
157#endif
158write(*,*) "Average of surface pressure computed."
[2794]159
[3602]160! CO2 vmr
161!--------
162call get_var3("co2_layer1",var3_read_1)
163where (var3_read_1 < 0.)
164    var3_read_1 = 1.e-10
165else where (var3_read_1 > 1.)
166    var3_read_1 = 1.
167end where
168#ifndef CPP_1D
169    call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,var3_read_1,q_co2)
170#else
171    q_co2(1,:) = var3_read_1(1,1,:)
172#endif
173A = (1./m_co2 - 1./m_noco2)
174B = 1./m_noco2
175vmr_co2_PCM = q_co2/(A*q_co2 + B)/m_co2
176write(*,*) "Data for CO2 vmr downloaded."
[2794]177
[3602]178! H2O vmr
179!--------
180call get_var3("h2o_layer1",var3_read_1)
181where (var3_read_1 < 0.)
182    var3_read_1 = 1.e-10
183else where (var3_read_1 > 1.)
184    var3_read_1 = 1.
185end where
186#ifndef CPP_1D
187    call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,var3_read_1,q_h2o)
188#else
189    q_h2o(1,:) = var3_read_1(1,1,:)
190#endif
191write(*,*) "Data for H2O vmr downloaded."
192
[3609]193do islope = 1,nslope
194    if (nslope /= 1) then
[3611]195        num='  '
[3609]196        write(num,'(i2.2)') islope
197        num = '_slope'//num
198    endif
199
[3602]200    ! CO2 ice
201    !--------
[3609]202    call get_var3("co2ice"//num,var3_read_1)
[3602]203    where (var3_read_1 < 0.) var3_read_1 = 0.
[3609]204    write(*,*) "Data for co2_ice"//num//" downloaded."
205    call get_var3("perennial_co2ice"//num,var3_read_2)
206    write(*,*) "Data for perennial_co2ice"//num//" downloaded."
[2835]207
[3602]208    ! Compute the minimum over the year for each point
209    var2_read = minval(var3_read_1 + var3_read_2,3)
210#ifndef CPP_1D
[3609]211    call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_co2_ice(:,islope,2))
[3602]212#else
[3609]213    min_co2_ice(1,islope,2) = var2_read(1,1)
[3602]214#endif
[3609]215    write(*,*) "Min of co2_ice"//num//" computed."
[3602]216
217    ! H2O ice
218    !--------
[3609]219    call get_var3("h2o_ice_s"//num,var3_read_1)
[3602]220    where (var3_read_1 < 0.) var3_read_1 = 0.
[3609]221    write(*,*) "Data for h2o_ice_s"//num//" downloaded."
222    call get_var3("watercap"//num,var3_read_2)
223    write(*,*) "Data for watercap"//num//" downloaded."
[2835]224
[3602]225    ! Compute the minimum over the year for each point
226    var2_read = minval(var3_read_1 + var3_read_2,3)
227#ifndef CPP_1D
[3609]228    call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_h2o_ice(:,islope,2))
[3602]229#else
[3609]230    min_h2o_ice(1,islope,2) = var2_read(1,1)
[3602]231#endif
[3609]232    write(*,*) "Min of h2o_ice"//num//" computed."
[3602]233
234    ! Tsurf
235    !------
[3609]236    call get_var3("tsurf"//num,var3_read_1)
[3620]237    write(*,*) "Data for tsurf"//num//" downloaded."
[3143]238
[3602]239    ! Compute average over the year for each point
240    var2_read = sum(var3_read_1,3)/timelen
241#ifndef CPP_1D
[3609]242    call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,tsurf_avg(:,islope))
[3602]243#else
[3609]244    tsurf_avg(1,islope) = var2_read(1,1)
[3602]245#endif
[3609]246    write(*,*) "Average of tsurf"//num//" computed."
[3130]247
[3122]248    if (soil_pem) then
[3602]249        ! Tsoil
250        !------
[3609]251        call get_var4("soiltemp"//num,var4_read)
252        write(*,*) "Data for soiltemp"//num//" downloaded."
[2794]253
[3602]254        ! Compute average over the year for each point
255        var3_read_3 = sum(var4_read,4)/timelen
256#ifndef CPP_1D
[3620]257        do t = 1,timelen
258            call gr_dyn_fi(nsoilmx,iim_input + 1,jjm_input + 1,ngrid,var4_read(:,:,:,t),tsoil_timeseries(:,:,islope,t))
[3602]259        enddo
[3609]260        call gr_dyn_fi(nsoilmx,iim_input + 1,jjm_input + 1,ngrid,var3_read_3,tsoil_avg(:,:,islope))
[3602]261#else
[3609]262        tsoil_timeseries(1,:,islope,:) = var4_read(1,1,:,:)
263        tsoil_avg(1,:,islope) = var3_read_3(1,1,:)
[3602]264#endif
[3609]265        write(*,*) "Average of tsoil"//num//" computed."
[3602]266
267        ! Soil water density
268        !-------------------
[3609]269        call get_var4("waterdensity_soil"//num,var4_read)
[3602]270#ifndef CPP_1D
[3620]271        do t = 1,timelen
272            call gr_dyn_fi(nsoilmx,iim_input + 1,jjm_input + 1,ngrid,var4_read(:,:,:,t),watersoil_density_timeseries(:,:,islope,t))
[3602]273        enddo
274#else
[3609]275        watersoil_density_timeseries(1,:,islope,:) = var4_read(1,1,:,:)
[3602]276#endif
[3609]277        write(*,*) "Data for waterdensity_soil"//num//" downloaded."
[2794]278
[3602]279        ! Surface water density
280        !----------------------
[3609]281        call get_var3("waterdensity_surface"//num,var3_read_1)
282        write(*,*) "Data for waterdensity_surface"//num//" downloaded."
[3602]283
284        ! Compute average over the year for each point
285        var2_read = sum(var3_read_1,3)/timelen
286#ifndef CPP_1D
[3609]287        call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,watersurf_density_avg(:,islope))
[3602]288#else
[3609]289        watersurf_density_avg(1,islope) = var2_read(1,1)
[3122]290#endif
[3609]291        write(*,*) "Average of waterdensity_surface"//num//" computed."
[3602]292    endif
[3609]293enddo
294deallocate(var2_read,var3_read_1,var3_read_2,var3_read_3,var4_read,num)
[3602]295
296! Close the NetCDF file
297call error_msg(nf90_close(fID),"close",filename2)
[3616]298write(*,*) '> "'//filename2//'" processed!'
[3602]299
[3096]300END SUBROUTINE read_data_PCM
[2897]301
[3122]302!=======================================================================
303
[3076]304SUBROUTINE check_dim(n1,n2,str1,str2)
[2980]305
[3076]306implicit none
[2980]307
[3076]308integer,            intent(in) :: n1, n2
309character(len = *), intent(in) :: str1, str2
[2980]310
[3076]311character(256) :: s1, s2
[2794]312
[3076]313if (n1 /= n2) then
314    s1 = 'value of '//trim(str1)//' ='
315    s2 = ' read in starting file differs from parametrized '//trim(str2)//' ='
316    write(msg,'(10x,a,i4,2x,a,i4)')trim(s1),n1,trim(s2),n2
317    call abort_gcm(trim(modname),trim(msg),1)
318endif
319
[2794]320END SUBROUTINE check_dim
321
[3076]322!=======================================================================
[2794]323
[3130]324SUBROUTINE get_var1(var,v)
[3076]325
326implicit none
327
[3130]328character(len = *), intent(in)  :: var
[3076]329real, dimension(:), intent(out) :: v
330
[3130]331call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var)
332call error_msg(NF90_GET_VAR(fID,vID,v),"get",var)
[3076]333
[2794]334END SUBROUTINE get_var1
335
[3076]336!=======================================================================
[2794]337
[3130]338SUBROUTINE get_var3(var,v) ! on U grid
[2794]339
[3076]340implicit none
341
[3130]342character(len = *),     intent(in)  :: var
[3076]343real, dimension(:,:,:), intent(out) :: v
344
[3130]345call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var)
346call error_msg(NF90_GET_VAR(fID,vID,v),"get",var)
[3076]347
[2794]348END SUBROUTINE get_var3
349
[3076]350!=======================================================================
351
[3130]352SUBROUTINE get_var4(var,v)
[3076]353
354implicit none
355
[3130]356character(len = *),       intent(in)  :: var
[3076]357real, dimension(:,:,:,:), intent(out) :: v
358
[3130]359call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var)
360call error_msg(NF90_GET_VAR(fID,vID,v),"get",var)
[3076]361
[2794]362END SUBROUTINE get_var4
363
[3076]364!=======================================================================
[2794]365
[3076]366SUBROUTINE error_msg(ierr,typ,nam)
367
368implicit none
369
370integer,            intent(in) :: ierr !--- NetCDF ERROR CODE
371character(len = *), intent(in) :: typ  !--- TYPE OF OPERATION
372character(len = *), intent(in) :: nam  !--- FIELD/FILE NAME
373
374if (ierr == nf90_noerr) return
375select case(typ)
376    case('inq');   msg="Field <"//trim(nam)//"> is missing"
377    case('get');   msg="Reading failed for <"//trim(nam)//">"
[3363]378    case('put');   msg="Writing failed for <"//trim(nam)//">"
[3076]379    case('open');  msg="File opening failed for <"//trim(nam)//">"
380    case('close'); msg="File closing failed for <"//trim(nam)//">"
381    case default
382        write(*,*) 'There is no message for this error.'
383        error stop
384end select
385call abort_gcm(trim(modname),trim(msg),ierr)
386
387END SUBROUTINE error_msg
388
[3096]389END MODULE read_data_PCM_mod
[3602]390
Note: See TracBrowser for help on using the repository browser.