MODULE read_data_PCM_mod use netcdf, only: nf90_open, NF90_NOWRITE, nf90_noerr, nf90_strerror, & nf90_get_var, nf90_inq_varid, nf90_inq_dimid, & nf90_inquire_dimension, nf90_close implicit none character(13), parameter :: modname = 'read_data_PCM' character(256) :: msg ! for reading integer :: fID, vID ! for reading !======================================================================= contains !======================================================================= SUBROUTINE read_data_PCM(filename,timelen,iim_input,jjm_input,ngrid,nslope,vmr_co2_PCM_phys,ps_timeseries, & min_co2_ice,min_h2o_ice,tsurf_avg,tsoil_avg,tsurf_PCM,tsoil_PCM,q_co2,q_h2o, & watersurf_density_avg,watersoil_density) use comsoil_h, only: nsoilmx use comsoil_h_PEM, only: soil_pem use constants_marspem_mod, only: m_co2, m_noco2 implicit none !======================================================================= ! ! Purpose: Read initial confitions file from the PCM ! ! Authors: RV & LL !======================================================================= include "dimensions.h" !======================================================================= ! Arguments: character(*), intent(in) :: filename ! File name integer, intent(in) :: timelen ! number of times stored in the file integer, intent(in) :: iim_input, jjm_input, ngrid, nslope ! number of points in the lat x lon dynamical grid, number of subgrid slopes ! Ouputs real, dimension(ngrid,nslope), intent(out) :: min_co2_ice ! Minimum of co2 ice per slope of the year [kg/m^2] real, dimension(ngrid,nslope), intent(out) :: min_h2o_ice ! Minimum of h2o ice per slope of the year [kg/m^2] real, dimension(ngrid,timelen), intent(out) :: vmr_co2_PCM_phys ! Physics x Times co2 volume mixing ratio retrieve from the PCM [m^3/m^3] real, dimension(ngrid,timelen), intent(out) :: ps_timeseries ! Surface Pressure [Pa] real, dimension(ngrid,timelen), intent(out) :: q_co2 ! CO2 mass mixing ratio in the first layer [kg/m^3] real, dimension(ngrid,timelen), intent(out) :: q_h2o ! H2O mass mixing ratio in the first layer [kg/m^3] !SOIL real, dimension(ngrid,nslope), intent(out) :: tsurf_avg ! Average surface temperature of the concatenated file [K] real, dimension(ngrid,nsoilmx,nslope), intent(out) :: tsoil_avg ! Average soil temperature of the concatenated file [K] real, dimension(ngrid,nslope,timelen), intent(out) :: tsurf_PCM ! Surface temperature of the concatenated file, time series [K] real, dimension(ngrid,nsoilmx,nslope,timelen), intent(out) :: tsoil_PCM ! Soil temperature of the concatenated file, time series [K] real, dimension(ngrid,nslope), intent(out) :: watersurf_density_avg ! Water density at the surface [kg/m^3] real, dimension(ngrid,nsoilmx,nslope,timelen), intent(out) :: watersoil_density ! Water density in the soil layer, time series [kg/m^3] !======================================================================= ! Local Variables integer :: i, j, l, t ! loop variables real :: A, B, mmean ! Molar Mass of co2 and no co2, A;B intermediate variables to compute the mean molar mass of the layer integer :: islope ! loop for variables character(2) :: num ! for reading sloped variables real, dimension(iim_input + 1,jjm_input + 1,nslope,timelen) :: h2o_ice_s_dyn ! h2o ice per slope of the concatenated file [kg/m^2] real, dimension(iim_input + 1,jjm_input + 1,nslope,timelen) :: watercap real, dimension(iim_input + 1,jjm_input + 1,nslope,timelen) :: perennial_co2ice real, dimension(iim_input + 1,jjm_input + 1,timelen) :: vmr_co2_PCM ! CO2 volume mixing ratio in the first layer [mol/m^3] real, dimension(iim_input + 1,jjm_input + 1,timelen) :: ps_PCM ! Surface Pressure [Pa] real, dimension(iim_input + 1,jjm_input + 1,nslope) :: min_co2_ice_dyn real, dimension(iim_input + 1,jjm_input + 1,nslope) :: min_h2o_ice_dyn real, dimension(iim_input + 1,jjm_input + 1,nslope) :: tsurf_avg_dyn ! Average surface temperature of the concatenated file [K] real, dimension(iim_input + 1,jjm_input + 1,nsoilmx,nslope) :: tsoil_avg_dyn ! Average soil temperature of the concatenated file [K] real, dimension(iim_input + 1,jjm_input + 1,nslope,timelen) :: tsurf_PCM_dyn ! Surface temperature of the concatenated file, time series [K] real, dimension(iim_input + 1,jjm_input + 1,nsoilmx,nslope,timelen) :: tsoil_PCM_dyn ! Soil temperature of the concatenated file, time series [K] real, dimension(iim_input + 1,jjm_input + 1,timelen) :: q_co2_dyn ! CO2 mass mixing ratio in the first layer [kg/m^3] real, dimension(iim_input + 1,jjm_input + 1,timelen) :: q_h2o_dyn ! H2O mass mixing ratio in the first layer [kg/m^3] real, dimension(iim_input + 1,jjm_input + 1,nslope,timelen) :: co2_ice_slope_dyn ! co2 ice amount per slope of the year [kg/m^2] real, dimension(iim_input + 1,jjm_input + 1,nslope,timelen) :: watersurf_density_dyn ! Water density at the surface, time series [kg/m^3] real, dimension(iim_input + 1,jjm_input + 1,nslope) :: watersurf_density_dyn_avg ! Water density at the surface, dynamic grid, yearly averaged [kg/m^3] real, dimension(iim_input + 1,jjm_input + 1,nsoilmx,nslope,timelen) :: watersoil_density_dyn ! Water density in the soil layer, time series [kg/m^3] !----------------------------------------------------------------------- ! Open the NetCDF file write(*,*) "Opening "//filename//"..." call error_msg(NF90_OPEN(filename,NF90_NOWRITE,fID),"open",filename) ! Download the data from the file call get_var3("ps",ps_PCM) write(*,*) "Data for surface pressure downloaded!" call get_var3("co2_layer1",q_co2_dyn) write(*,*) "Data for vmr co2 downloaded!" call get_var3("h2o_layer1",q_h2o_dyn) write(*,*) "Data for vmr h2o downloaded!" if (nslope == 1) then ! There is no slope call get_var3("co2ice",co2_ice_slope_dyn(:,:,1,:)) write(*,*) "Data for co2_ice downloaded!" call get_var3("h2o_ice_s",h2o_ice_s_dyn(:,:,1,:)) write(*,*) "Data for h2o_ice_s downloaded!" call get_var3("tsurf",tsurf_PCM_dyn(:,:,1,:)) write(*,*) "Data for tsurf downloaded!" #ifndef CPP_STD call get_var3("watercap",watercap(:,:,1,:)) write(*,*) "Data for watercap downloaded!" call get_var3("perennial_co2ice",perennial_co2ice(:,:,1,:)) write(*,*) "Data for perennial_co2ice downloaded!" if (soil_pem) then call get_var4("soiltemp",tsoil_PCM_dyn(:,:,:,1,:)) write(*,*) "Data for soiltemp downloaded!" call get_var4("waterdensity_soil",watersoil_density_dyn(:,:,:,1,:)) write(*,*) "Data for waterdensity_soil downloaded!" call get_var3("waterdensity_surface",watersurf_density_dyn(:,:,1,:)) write(*,*) "Data for waterdensity_surface downloaded!" endif !soil_pem #endif else ! We use slopes do islope = 1,nslope write(num,'(i2.2)') islope call get_var3("co2ice_slope"//num,co2_ice_slope_dyn(:,:,islope,:)) enddo write(*,*) "Data for co2_ice downloaded!" do islope = 1,nslope write(num,'(i2.2)') islope call get_var3("h2o_ice_s_slope"//num,h2o_ice_s_dyn(:,:,islope,:)) enddo write(*,*) "Data for h2o_ice_s downloaded!" do islope = 1,nslope write(num,'(i2.2)') islope call get_var3("tsurf_slope"//num,tsurf_PCM_dyn(:,:,islope,:)) enddo write(*,*) "Data for tsurf downloaded!" #ifndef CPP_STD do islope = 1,nslope write(num,'(i2.2)') islope call get_var3("watercap_slope"//num,watercap(:,:,islope,:)) enddo write(*,*) "Data for watercap downloaded!" do islope = 1,nslope write(num,'(i2.2)') islope call get_var3("perennial_co2ice_slope"//num,perennial_co2ice(:,:,islope,:)) enddo write(*,*) "Data for perennial_co2ice downloaded!" if (soil_pem) then do islope = 1,nslope write(num,'(i2.2)') islope call get_var4("soiltemp_slope"//num,tsoil_PCM_dyn(:,:,:,islope,:)) enddo write(*,*) "Data for soiltemp downloaded!" do islope = 1,nslope write(num,'(i2.2)') islope call get_var4("waterdensity_soil_slope"//num,watersoil_density_dyn(:,:,:,islope,:)) enddo write(*,*) "Data for waterdensity_soil downloaded!" do islope = 1,nslope write(num,'(i2.2)') islope call get_var3("waterdensity_surface"//num,watersurf_density_dyn(:,:,islope,:)) enddo write(*,*) "Data for waterdensity_surface downloaded!" endif !soil_pem #endif endif ! Close the NetCDF file write(*,*) "Closing "//filename//"..." call error_msg(nf90_close(fID),"close",filename) ! Compute the minimum over the year for each point write(*,*) "Computing the min of h2o_ice..." where (h2o_ice_s_dyn < 0.) h2o_ice_s_dyn = 0. min_h2o_ice_dyn = minval(h2o_ice_s_dyn + watercap,4) write(*,*) "Computing the min of co2_ice..." where (co2_ice_slope_dyn < 0.) co2_ice_slope_dyn = 0. min_co2_ice_dyn = minval(co2_ice_slope_dyn + perennial_co2ice,4) ! Compute averages over the year for each point write(*,*) "Computing the average of tsurf..." tsurf_avg_dyn = sum(tsurf_PCM_dyn,4)/timelen if (soil_pem) then write(*,*) "Computing average of tsoil..." tsoil_avg_dyn = sum(tsoil_PCM_dyn,5)/timelen write(*,*) "Computing average of waterdensity_surface..." watersurf_density_dyn_avg = sum(watersurf_density_dyn,4)/timelen endif ! By definition, we get rid of the negative values A = (1./m_co2 - 1./m_noco2) B = 1./m_noco2 do i = 1,iim + 1 do j = 1,jjm_input + 1 do t = 1, timelen if (q_co2_dyn(i,j,t) < 0) then q_co2_dyn(i,j,t) = 1.e-10 else if (q_co2_dyn(i,j,t) > 1) then q_co2_dyn(i,j,t) = 1. endif if (q_h2o_dyn(i,j,t) < 0) then q_h2o_dyn(i,j,t) = 1.e-10 else if (q_h2o_dyn(i,j,t) > 1) then q_h2o_dyn(i,j,t) = 1. endif mmean = 1/(A*q_co2_dyn(i,j,t) + B) vmr_co2_PCM(i,j,t) = q_co2_dyn(i,j,t)*mmean/m_co2 enddo enddo enddo #ifndef CPP_1D call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,vmr_co2_PCM,vmr_co2_PCM_phys) call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,ps_PCM,ps_timeseries) call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,q_co2_dyn,q_co2) call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,q_h2o_dyn,q_h2o) do islope = 1,nslope call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,min_co2_ice_dyn(:,:,islope),min_co2_ice(:,islope)) call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,min_h2o_ice_dyn(:,:,islope),min_h2o_ice(:,islope)) if (soil_pem) then do l = 1,nsoilmx call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,tsoil_avg_dyn(:,:,l,islope),tsoil_avg(:,l,islope)) do t = 1,timelen call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,tsoil_PCM_dyn(:,:,l,islope,t),tsoil_PCM(:,l,islope,t)) call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,watersoil_density_dyn(:,:,l,islope,t),watersoil_density(:,l,islope,t)) enddo enddo call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,watersurf_density_dyn_avg(:,:,islope),watersurf_density_avg(:,islope)) endif ! soil_pem do t = 1,timelen call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,tsurf_PCM_dyn(:,:,islope,t),tsurf_PCM(:,islope,t)) enddo enddo call gr_dyn_fi(nslope,iim_input + 1,jjm_input + 1,ngrid,tsurf_avg_dyn,tsurf_avg) #else vmr_co2_PCM_phys(1,:) = vmr_co2_PCM(1,1,:) ps_timeseries(1,:) = ps_PCM(1,1,:) q_co2(1,:) = q_co2_dyn(1,1,:) q_h2o(1,:) = q_h2o_dyn(1,1,:) min_co2_ice(1,:) = min_co2_ice_dyn(1,1,:) min_h2o_ice(1,:) = min_h2o_ice_dyn(1,1,:) if (soil_pem) then tsoil_avg(1,:,:) = tsoil_avg_dyn(1,1,:,:) tsoil_PCM(1,:,:,:) = tsoil_PCM_dyn(1,1,:,:,:) watersoil_density(1,:,:,:) = watersoil_density_dyn(1,1,:,:,:) watersurf_density_avg(1,:) = watersurf_density_dyn_avg(1,1,:) endif ! soil_pem tsurf_PCM(1,:,:) = tsurf_PCM_dyn(1,1,:,:) tsurf_avg(1,:) = tsurf_avg_dyn(1,1,:) #endif END SUBROUTINE read_data_PCM !======================================================================= SUBROUTINE check_dim(n1,n2,str1,str2) implicit none integer, intent(in) :: n1, n2 character(len = *), intent(in) :: str1, str2 character(256) :: s1, s2 if (n1 /= n2) then s1 = 'value of '//trim(str1)//' =' s2 = ' read in starting file differs from parametrized '//trim(str2)//' =' write(msg,'(10x,a,i4,2x,a,i4)')trim(s1),n1,trim(s2),n2 call abort_gcm(trim(modname),trim(msg),1) endif END SUBROUTINE check_dim !======================================================================= SUBROUTINE get_var1(var,v) implicit none character(len = *), intent(in) :: var real, dimension(:), intent(out) :: v call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var) call error_msg(NF90_GET_VAR(fID,vID,v),"get",var) END SUBROUTINE get_var1 !======================================================================= SUBROUTINE get_var3(var,v) ! on U grid implicit none character(len = *), intent(in) :: var real, dimension(:,:,:), intent(out) :: v call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var) call error_msg(NF90_GET_VAR(fID,vID,v),"get",var) END SUBROUTINE get_var3 !======================================================================= SUBROUTINE get_var4(var,v) implicit none character(len = *), intent(in) :: var real, dimension(:,:,:,:), intent(out) :: v call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var) call error_msg(NF90_GET_VAR(fID,vID,v),"get",var) END SUBROUTINE get_var4 !======================================================================= SUBROUTINE error_msg(ierr,typ,nam) implicit none integer, intent(in) :: ierr !--- NetCDF ERROR CODE character(len = *), intent(in) :: typ !--- TYPE OF OPERATION character(len = *), intent(in) :: nam !--- FIELD/FILE NAME if (ierr == nf90_noerr) return select case(typ) case('inq'); msg="Field <"//trim(nam)//"> is missing" case('get'); msg="Reading failed for <"//trim(nam)//">" case('put'); msg="Writing failed for <"//trim(nam)//">" case('open'); msg="File opening failed for <"//trim(nam)//">" case('close'); msg="File closing failed for <"//trim(nam)//">" case default write(*,*) 'There is no message for this error.' error stop end select call abort_gcm(trim(modname),trim(msg),ierr) END SUBROUTINE error_msg END MODULE read_data_PCM_mod