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

Last change on this file since 3609 was 3609, checked in by jbclement, 6 days ago

PEM:

  • Simplification of "read_data_PCM_mod.F90" to treat the cases with/without slope.
  • Few bug corrections related to 1D.

JBC

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