MODULE regr_horiz_time_climoz_m USE interpolation, ONLY: locate USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat, grid_type, unstructured USE nrtype, ONLY: pi USE netcdf, ONLY: NF90_CLOBBER, NF90_FLOAT, NF90_GET_VAR, NF90_OPEN, & NF90_NOWRITE, NF90_NOERR, NF90_GET_ATT, NF90_GLOBAL USE netcdf95, ONLY: NF95_DEF_DIM, NF95_INQ_DIMID, NF95_INQUIRE_DIMENSION, & NF95_DEF_VAR, NF95_INQ_VARID, NF95_INQUIRE_VARIABLE, & NF95_OPEN, NF95_CREATE, NF95_GET_ATT, NF95_GW_VAR, HANDLE_ERR, & NF95_CLOSE, NF95_ENDDEF, NF95_PUT_ATT, NF95_PUT_VAR, NF95_COPY_ATT USE print_control_mod, ONLY: lunout USE dimphy IMPLICIT NONE PRIVATE PUBLIC :: regr_horiz_time_climoz REAL, PARAMETER :: deg2rad=pi/180. CHARACTER(LEN=13), PARAMETER :: vars_in(2)=['tro3 ','tro3_daylight'] INTEGER :: nlat_ou, nlon_ou REAL, ALLOCATABLE :: latitude_glo(:) !$OMP THREADPRIVATE(latitude_glo) INTEGER, ALLOCATABLE :: ind_cell_glo_glo(:) !$OMP THREADPRIVATE(ind_cell_glo_glo) CONTAINS !------------------------------------------------------------------------------- ! SUBROUTINE regr_horiz_time_climoz(read_climoz,interpt) ! !------------------------------------------------------------------------------- ! Purpose: Regrid horizontally and in time zonal or 3D ozone climatologies. ! * Read ozone climatology from netcdf file ! * Regrid it horizontaly to LMDZ grid (quasi-conservative method) ! * If interpt=T, interpolate linearly in time (one record each day) ! If interpt=F, keep original time sampling (14 months). ! * Save it to a new netcdf file. !------------------------------------------------------------------------------- ! Remarks: ! * Up to 2 variables treated: "tro3" and "tro3_daylight" (if read_climoz=2) ! * Input fields coordinates: (longitudes, latitudes, pressure_levels, time) ! * Output grid cells centers coordinates given by [rlonv,] rlatu. ! * Output grid cells edges coordinates given by [rlonu,] rlatv. ! * Input file [longitudes and] latitudes given in degrees. ! * Input file pressure levels are given in Pa or hPa. ! * All coordinates variables are stricly monotonic. ! * Monthly fields are interpolated linearly in time to get daily values. ! * Fields are known at the middle of the months, so interpolation requires an ! additional record both for 1st half of january and 2nd half of december: ! - For a 14-records "climoz.nc": records 1 and 14. ! - For 12-records files: ! record 12 of "climoz_m.nc" if available, or record 1 of "climoz.nc". ! record 1 of "climoz_p.nc" if available, or record 12 of "climoz.nc". ! * Calendar is taken into account to get one record each day (not 360 always). ! * Missing values are filled in from sky to ground by copying lowest valid one. ! Attribute "missing_value" or "_FillValue" must be present in input file. !------------------------------------------------------------------------------- USE assert_m, ONLY: assert USE cal_tools_m, ONLY: year_len, mid_month !! USE control_mod, ONLY: anneeref USE time_phylmdz_mod, ONLY: annee_ref USE ioipsl, ONLY: ioget_year_len, ioget_calendar USE regr_conserv_m, ONLY: regr_conserv USE regr_lint_m, ONLY: regr_lint USE regular_lonlat_mod, ONLY: boundslon_reg, boundslat_reg, south, west, east USE slopes_m, ONLY: slopes USE xios USE mod_phys_lmdz_para, ONLY: is_mpi_root, is_master, is_omp_master, gather, gather_mpi, bcast_mpi, klon_mpi USE geometry_mod, ONLY : latitude_deg, ind_cell_glo USE mod_grid_phy_lmdz, ONLY: klon_glo !------------------------------------------------------------------------------- ! Arguments: INTEGER, INTENT(IN) :: read_climoz ! read ozone climatology, 1 or 2 ! 1: read a single ozone climatology used day and night ! 2: same + read also a daylight climatology LOGICAL, INTENT(IN) :: interpt ! TRUE => daily interpolation ! FALSE => no interpolation (14 months) !------------------------------------------------------------------------------- ! Local variables: !--- Input files variables INTEGER :: nlon_in ! Number of longitudes INTEGER :: nlat_in ! Number of latitudes INTEGER :: nlev_in ! Number of pressure levels INTEGER :: nmth_in ! Number of months REAL, POINTER :: lon_in(:) ! Longitudes (ascending order, rad) REAL, POINTER :: lat_in(:) ! Latitudes (ascending order, rad) REAL, POINTER :: lev_in(:) ! Pressure levels (ascen. order, hPa) REAL, ALLOCATABLE :: lon_in_edge(:) ! Longitude intervals edges ! (ascending order, / ) REAL, ALLOCATABLE :: sinlat_in_edge(:) ! Sinus of latitude intervals edges ! (ascending order, / ) LOGICAL :: ldec_lon, ldec_lat, ldec_lev ! Decreasing order in input file CHARACTER(LEN=20) :: cal_in ! Calendar REAL, ALLOCATABLE :: o3_in3(:,:,:,:,:) ! Ozone climatologies REAL, ALLOCATABLE :: o3_in3bis(:,:,:,:,:) ! Ozone climatologies REAL, ALLOCATABLE :: o3_in2 (:,:,:,:) ! Ozone climatologies REAL, ALLOCATABLE :: o3_in2bis(:,:,:,:,:) ! Ozone climatologies ! last index: 1 for the day-night average, 2 for the daylight field. REAL :: NaN !--- Partially or totally regridded variables (:,:,nlev_in,:,read_climoz) REAL, ALLOCATABLE :: o3_regr_lon (:,:,:,:,:) ! (nlon_ou,nlat_in,:,0:13 ,:) REAL, ALLOCATABLE :: o3_regr_lonlat(:,:,:,:,:) ! (nlon_ou,nlat_ou,:,0:13 ,:) REAL, ALLOCATABLE :: o3_out3 (:,:,:,:,:) ! (nlon_ou,nlat_ou,:,ntim_ou,:) REAL, ALLOCATABLE :: o3_out3_glo (:,:,:,:) ! (nbp_lat,:,ntim_ou,:) REAL, ALLOCATABLE :: o3_regr_lat (:,:,:,:) ! (nlat_in,:,0:13 ,:) REAL, ALLOCATABLE :: o3_out2 (:,:,:,:) ! (nlat_ou,:,ntim_ou,:) REAL, ALLOCATABLE :: o3_out2_glo (:,:,:,:) ! (nbp_lat,:,ntim_ou,:) REAL, ALLOCATABLE :: o3_out (:,:,:,:) ! (nbp_lat,:,ntim_ou,:) ! Dimension number | Interval | Contains | For variables: ! 1 (longitude) | [rlonu(i-1), rlonu(i)] | rlonv(i) | all ! 2 (latitude) | [rlatv(j), rlatv(j-1)] | rlatu(j) | all but o3_regr_lon ! 3 (press level) | | lev(k) | all ! Note that rlatv(0)=pi/2 and rlatv(nlat_ou)=-pi/2. ! Dimension 4 is: month number (all vars but o3_out) ! days elapsed since Jan. 1st 0h at mid-day (o3_out only) REAL, ALLOCATABLE :: v1(:) !--- For NetCDF: INTEGER :: fID_in_m, fID_in, levID_ou, dimid, vID_in(read_climoz), ntim_ou INTEGER :: fID_in_p, fID_ou, timID_ou, varid, vID_ou(read_climoz), ndims, ncerr INTEGER, POINTER :: dIDs(:) CHARACTER(LEN=20) :: cal_ou !--- Calendar; no time inter => same as input CHARACTER(LEN=80) :: press_unit !--- Pressure unit REAL :: tmidmonth(0:13) !--- Elapsed days since Jan-1 0h at mid-months ! Additional records 0, 13 for interpolation REAL, ALLOCATABLE :: tmidday(:) !--- Output times (mid-days since Jan 1st 0h) LOGICAL :: lprev, lnext !--- Flags: previous/next files are present LOGICAL :: l3D, l2D !--- Flag: input fields are 3D or zonal INTEGER :: ii, i, j, k, l, m, dln, ib, ie, iv, dx1, dx2 INTEGER, ALLOCATABLE :: sta(:), cnt(:) CHARACTER(LEN=80) :: sub, dim_nam, msg REAL :: null_array(0) LOGICAL,SAVE :: first=.TRUE. !$OMP THREADPRIVATE(first) REAL, ALLOCATABLE :: test_o3_in(:,:) REAL, ALLOCATABLE :: test_o3_out(:) IF (grid_type==unstructured) THEN IF (first) THEN IF (is_master) THEN ALLOCATE(latitude_glo(klon_glo)) ALLOCATE(ind_cell_glo_glo(klon_glo)) ELSE ALLOCATE(latitude_glo(0)) ALLOCATE(ind_cell_glo_glo(0)) ENDIF CALL gather(latitude_deg, latitude_glo) CALL gather(ind_cell_glo, ind_cell_glo_glo) ENDIF ENDIF IF (is_omp_master) THEN nlat_ou=nbp_lat nlon_ou=nbp_lon !------------------------------------------------------------------------------- IF (is_mpi_root) THEN sub="regr_horiz_time_climoz" WRITE(lunout,*)"Call sequence information: "//TRIM(sub) CALL assert(read_climoz == 1 .OR. read_climoz == 2, "regr_lat_time_climoz") CALL NF95_OPEN("climoz.nc" , NF90_NOWRITE, fID_in) lprev=NF90_OPEN("climoz_m.nc", NF90_NOWRITE, fID_in_m)==NF90_NOERR lnext=NF90_OPEN("climoz_p.nc", NF90_NOWRITE, fID_in_p)==NF90_NOERR !--- Get coordinates from the input file. Converts lon/lat in radians. ! Few inversions because "regr_conserv" and gcm need ascending vectors. CALL NF95_INQ_VARID(fID_in, vars_in(1), varid) CALL NF95_INQUIRE_VARIABLE(fID_in, varid, dimids=dIDs, ndims=ndims) l3D=ndims==4; l2D=ndims==3 IF(l3D) WRITE(lunout,*)"Input files contain full 3D ozone fields." IF(l2D) WRITE(lunout,*)"Input files contain zonal 2D ozone fields." DO i=1,ndims CALL NF95_INQUIRE_DIMENSION(fID_in, dIDs(i), name=dim_nam, nclen=dln) CALL NF95_INQ_VARID(fID_in, dim_nam, varid) ii=i; IF(l2D) ii=i+1 !--- ndims==3:NO LONGITUDE SELECT CASE(ii) CASE(1) !--- LONGITUDE CALL NF95_GW_VAR(fID_in, varid, lon_in) ldec_lon=lon_in(1)>lon_in(dln); IF(ldec_lon) lon_in=lon_in(dln:1:-1) nlon_in=dln; lon_in=lon_in*deg2rad CASE(2) !--- LATITUDE CALL NF95_GW_VAR(fID_in, varid, lat_in) ldec_lat=lat_in(1)>lat_in(dln); IF(ldec_lat) lat_in=lat_in(dln:1:-1) nlat_in=dln; lat_in=lat_in*deg2rad CASE(3) !--- PRESSURE LEVELS CALL NF95_GW_VAR(fID_in, varid, lev_in) ldec_lev=lev_in(1)>lev_in(dln); IF(ldec_lev) lev_in=lev_in(dln:1:-1) nlev_in=dln CALL NF95_GET_ATT(fID_in, varid, "units", press_unit) k=LEN_TRIM(press_unit) DO WHILE(ICHAR(press_unit(k:k))==0) press_unit(k:k)=' '; k=LEN_TRIM(press_unit) !--- REMOVE NULL END CHAR END DO IF(press_unit == "Pa") THEN lev_in = lev_in/100. !--- CONVERT TO hPa ELSE IF(press_unit /= "hPa") THEN CALL abort_physic(sub, "the only recognized units are Pa and hPa.",1) END IF CASE(4) !--- TIME CALL NF95_INQUIRE_DIMENSION(fID_in, dIDs(i), nclen=nmth_in) cal_in='gregorian' IF(NF90_GET_ATT(fID_in, varid, 'calendar', cal_in)/=NF90_NOERR) & WRITE(lunout,*)'WARNING: missing "calendar" attribute for "'// & TRIM(dim_nam)//'" in "climoz.nc". Choosing default: "gregorian".' k=LEN_TRIM(cal_in) DO WHILE(ICHAR(cal_in(k:k))==0) cal_in(k:k)=' '; k=LEN_TRIM(cal_in) !--- REMOVE NULL END CHAR END DO END SELECT END DO !--- Ensure first input longitudes interval contains first output point boundslon_reg(1,west) dx1=locate(v1,boundslon_reg(1,west))-1 v1=CSHIFT(v1,SHIFT=dx1,DIM=1); v1(nlon_in-dx1+2:)=v1(nlon_in-dx1+2:)+2.*pi !--- Extend input longitudes vector until last interval contains boundslon_reg(nlon_ou,east) dx2=0; DO WHILE(v1(1+dx2)+2.*pi<=boundslon_reg(nlon_ou,east)); dx2=dx2+1; END DO IF (first) THEN first=.FALSE. RETURN ENDIF ENDIF IF (is_mpi_root) THEN !--- Longitudes management: ! * Need to shift data if the origin of input file longitudes /= -pi ! * Need to add some margin in longitude to ensure input interval contains ! all the output intervals => at least one longitudes slice has to be ! duplicated, possibly more for undersampling. IF(l3D) THEN IF (grid_type==unstructured) THEN dx2=0 ELSE !--- Compute input edges longitudes vector (no end point yet) ALLOCATE(v1(nlon_in+1)) v1(1)=(lon_in(nlon_in)+lon_in(1))/2.-pi FORALL(i=2:nlon_in) v1(i)=(lon_in(i-1)+lon_in(i))/2. v1(nlon_in+1)=v1(1)+2.*pi DEALLOCATE(lon_in) !--- Shift input longitudes vector until it contains first output point boundslon_reg(1,west) v1=v1+2*pi*REAL(FLOOR((boundslon_reg(1,west)-v1(1))/(2.*pi))) !--- Ensure first input longitudes interval contains first output point boundslon_reg(1,west) dx1=locate(v1,boundslon_reg(1,west))-1 v1=CSHIFT(v1,SHIFT=dx1,DIM=1) v1(nlon_in-dx1+1:)=v1(nlon_in-dx1+1:)+2.*pi !--- Extend input longitudes vector until last interval contains boundslon_reg(nlat_ou,east) dx2=0; DO WHILE(v1(1+dx2)+2.*pi1) & o3_in2(:j-1,:,l,m)=SPREAD(o3_in2(j,:,l,m),DIM=1,ncopies=j-1) !--- Next to north pole j=nlat_in; DO WHILE(o3_in2(j,1,l,m)==NaN); j=j+1; END DO IF(j