MODULE limit ! !******************************************************************************* ! Author : L. Fairhead, 27/01/94 !------------------------------------------------------------------------------- ! Purpose: Boundary conditions files building for new model using climatologies. ! Both grids have to be regular. !------------------------------------------------------------------------------- ! Note: This routine is designed to work for Earth !------------------------------------------------------------------------------- ! Modification history: ! * 23/03/1994: Z. X. Li ! * 09/1999: L. Fairhead (netcdf reading in LMDZ.3.3) ! * 07/2001: P. Le Van ! * 11/2009: L. Guez (ozone day & night climatos, see etat0_netcdf.F90) ! * 12/2009: D. Cugnet (f77->f90, calendars, files from coupled runs) !------------------------------------------------------------------------------- USE ioipsl, ONLY: flininfo, flinopen, flinget, flinclo, & ioconf_calendar, ioget_calendar, lock_calendar, ioget_mon_len, ioget_year_len USE assert_eq_m, ONLY: assert_eq USE conf_dat_m, ONLY: conf_dat2d, conf_dat3d USE dimphy, ONLY: klon, zmasq USE phys_state_var_mod, ONLY: pctsrf, rlon, rlat USE control_mod USE init_ssrf_m, ONLY: start_init_subsurf CHARACTER(LEN=20), PARAMETER :: & fsst(4)=['amipbc_sst_1x1.nc ','cpl_atm_sst.nc ','histmth_sst.nc '& ,'sstk.nc '] CHARACTER(LEN=20), PARAMETER :: & fsic(4)=['amipbc_sic_1x1.nc ','cpl_atm_sic.nc ','histmth_sic.nc '& ,'ci.nc '] CHARACTER(LEN=10), PARAMETER :: & vsst(4)=['tosbcs ','SISUTESW ','tsol_oce ','sstk '], & vsic(4)=['sicbcs ','SIICECOV ','pourc_sic ','ci '] CHARACTER(LEN=10), PARAMETER :: & frugo='Rugos.nc ', falbe='Albedo.nc ', frelf='Relief.nc ', & vrug='RUGOS ', valb='ALBEDO ', vrel='RELIEF ' CONTAINS !------------------------------------------------------------------------------- ! SUBROUTINE limit_netcdf(masque, phis, extrap) ! !------------------------------------------------------------------------------- ! Author : L. Fairhead, 27/01/94 !------------------------------------------------------------------------------- ! Purpose: Boundary conditions files building for new model using climatologies. ! Both grids have to be regular. !------------------------------------------------------------------------------- ! Note: This routine is designed to work for Earth !------------------------------------------------------------------------------- ! Modification history: ! * 23/03/1994: Z. X. Li ! * 09/1999: L. Fairhead (netcdf reading in LMDZ.3.3) ! * 07/2001: P. Le Van ! * 11/2009: L. Guez (ozone day & night climatos, see etat0_netcdf.F90) ! * 12/2009: D. Cugnet (f77->f90, calendars, files from coupled runs) !------------------------------------------------------------------------------- #ifndef CPP_1D USE indice_sol_mod USE netcdf, ONLY: NF90_OPEN, NF90_CREATE, NF90_CLOSE, & NF90_DEF_DIM, NF90_DEF_VAR, NF90_PUT_VAR, NF90_PUT_ATT, & NF90_NOERR, NF90_NOWRITE, NF90_DOUBLE, NF90_GLOBAL, & NF90_CLOBBER, NF90_ENDDEF, NF90_UNLIMITED, NF90_FLOAT USE inter_barxy_m, ONLY: inter_barxy USE netcdf95, ONLY: nf95_def_var, nf95_put_att, nf95_put_var IMPLICIT NONE !------------------------------------------------------------------------------- ! Arguments: include "iniprint.h" include "dimensions.h" include "paramet.h" REAL, DIMENSION(iip1,jjp1), INTENT(INOUT) :: masque ! land mask REAL, DIMENSION(iip1,jjp1), INTENT(INOUT) :: phis ! ground geopotential LOGICAL, INTENT(IN) :: extrap ! SST extrapolation flag !------------------------------------------------------------------------------- ! Local variables: include "logic.h" include "comgeom2.h" include "comconst.h" !--- INPUT NETCDF FILES NAMES -------------------------------------------------- CHARACTER(LEN=20) :: icefile, sstfile, dumstr, fnam CHARACTER(LEN=10) :: varname !--- OUTPUT VARIABLES FOR NETCDF FILE ------------------------------------------ REAL :: fi_ice(klon), verif(klon) REAL, POINTER :: phy_rug(:,:)=>NULL(), phy_ice(:,:)=>NULL() REAL, POINTER :: phy_sst(:,:)=>NULL(), phy_alb(:,:)=>NULL() REAL, ALLOCATABLE :: phy_bil(:,:), pctsrf_t(:,:,:) INTEGER :: nbad !--- VARIABLES FOR OUTPUT FILE WRITING ----------------------------------------- INTEGER :: ierr, nid, ndim, ntim, k, dims(2), ix_sic, ix_sst INTEGER :: id_tim, id_SST, id_BILS, id_RUG, id_ALB INTEGER :: id_FOCE, id_FSIC, id_FTER, id_FLIC, varid_longitude, varid_latitude INTEGER :: NF90_FORMAT INTEGER :: ndays !--- Depending on the output calendar !--- INITIALIZATIONS ----------------------------------------------------------- #ifdef NC_DOUBLE NF90_FORMAT=NF90_DOUBLE #else NF90_FORMAT=NF90_FLOAT #endif CALL inigeom !--- MASK, GROUND GEOPOT. & SUBSURFACES COMPUTATION (IN CASE ok_etat0==.FALSE.) IF(ALL(masque==-99999.)) THEN CALL start_init_orog0(rlonv,rlatu,phis,masque) CALL gr_dyn_fi(1,iip1,jjp1,klon,masque,zmasq) !--- To physical grid ALLOCATE(rlon(klon),rlat(klon),pctsrf(klon,nbsrf)) CALL start_init_subsurf(.FALSE.) END IF !--- Beware: anneeref (from gcm.def) is used to determine output time sampling ndays=ioget_year_len(anneeref) !--- RUGOSITY TREATMENT -------------------------------------------------------- CALL msg(1,'Traitement de la rugosite') CALL get_2Dfield(frugo,vrug,'RUG',ndays,phy_rug,mask=masque(1:iim,:)) !--- OCEAN TREATMENT ----------------------------------------------------------- CALL msg(1,'Traitement de la glace oceanique') ! Input SIC file selection ! Open file only to test if available DO ix_sic=1,SIZE(fsic) IF ( NF90_OPEN(TRIM(fsic(ix_sic)),NF90_NOWRITE,nid)==NF90_NOERR ) THEN icefile=fsic(ix_sic); varname=vsic(ix_sic); EXIT END IF END DO IF(ix_sic==SIZE(fsic)+1) THEN WRITE(lunout,*) 'ERROR! No sea-ice input file was found.' WRITE(lunout,*) 'One of following files must be available : ' DO k=1,SIZE(fsic); WRITE(lunout,*) TRIM(fsic(k)); END DO CALL abort_physic('limit_netcdf','No sea-ice file was found',1) END IF CALL ncerr(NF90_CLOSE(nid),icefile) CALL msg(-1,'Fichier choisi pour la glace de mer:'//TRIM(icefile)) CALL get_2Dfield(icefile,varname, 'SIC',ndays,phy_ice) ALLOCATE(pctsrf_t(klon,nbsrf,ndays)) DO k=1,ndays fi_ice=phy_ice(:,k) WHERE(fi_ice>=1.0 ) fi_ice=1.0 WHERE(fi_ice=1.0-zmasq) pctsrf_t(:,is_sic,k)=1.0-zmasq pctsrf_t(:,is_oce,k)=0.0 ELSEWHERE pctsrf_t(:,is_oce,k)=1.0-zmasq-pctsrf_t(:,is_sic,k) WHERE(pctsrf_t(:,is_oce,k)0) WRITE(lunout,*) 'pb sous maille pour nb points = ',nbad nbad=COUNT(ABS(SUM(pctsrf_t(:,:,k),DIM=2)-1.0)>EPSFRA) IF(nbad>0) WRITE(lunout,*) 'pb sous surface pour nb points = ',nbad END DO DEALLOCATE(phy_ice) !--- SST TREATMENT ------------------------------------------------------------- CALL msg(1,'Traitement de la sst') ! Input SST file selection ! Open file only to test if available DO ix_sst=1,SIZE(fsst) IF ( NF90_OPEN(TRIM(fsst(ix_sst)),NF90_NOWRITE,nid)==NF90_NOERR ) THEN sstfile=fsst(ix_sst); varname=vsst(ix_sst); EXIT END IF END DO IF(ix_sst==SIZE(fsst)+1) THEN WRITE(lunout,*) 'ERROR! No sst input file was found.' WRITE(lunout,*) 'One of following files must be available : ' DO k=1,SIZE(fsst); WRITE(lunout,*) TRIM(fsst(k)); END DO CALL abort_physic('limit_netcdf','No sst file was found',1) END IF CALL ncerr(NF90_CLOSE(nid),sstfile) CALL msg(-1,'Fichier choisi pour la temperature de mer: '//TRIM(sstfile)) CALL get_2Dfield(sstfile,varname,'SST',ndays,phy_sst,flag=extrap) !--- ALBEDO TREATMENT ---------------------------------------------------------- CALL msg(1,'Traitement de l albedo') CALL get_2Dfield(falbe,valb,'ALB',ndays,phy_alb) !--- REFERENCE GROUND HEAT FLUX TREATMENT -------------------------------------- ALLOCATE(phy_bil(klon,ndays)); phy_bil=0.0 !--- OUTPUT FILE WRITING ------------------------------------------------------- CALL msg(5,'Ecriture du fichier limit : debut') fnam="limit.nc" !--- File creation CALL ncerr(NF90_CREATE(fnam,NF90_CLOBBER,nid),fnam) CALL ncerr(NF90_PUT_ATT(nid,NF90_GLOBAL,"title","Fichier conditions aux limites"),fnam) !--- Dimensions creation CALL ncerr(NF90_DEF_DIM(nid,"points_physiques",klon,ndim),fnam) CALL ncerr(NF90_DEF_DIM(nid,"time",NF90_UNLIMITED,ntim),fnam) dims=[ndim,ntim] !--- Variables creation CALL ncerr(NF90_DEF_VAR(nid,"TEMPS",NF90_FORMAT,[ntim],id_tim),fnam) CALL ncerr(NF90_DEF_VAR(nid,"FOCE", NF90_FORMAT,dims,id_FOCE),fnam) CALL ncerr(NF90_DEF_VAR(nid,"FSIC", NF90_FORMAT,dims,id_FSIC),fnam) CALL ncerr(NF90_DEF_VAR(nid,"FTER", NF90_FORMAT,dims,id_FTER),fnam) CALL ncerr(NF90_DEF_VAR(nid,"FLIC", NF90_FORMAT,dims,id_FLIC),fnam) CALL ncerr(NF90_DEF_VAR(nid,"SST", NF90_FORMAT,dims,id_SST),fnam) CALL ncerr(NF90_DEF_VAR(nid,"BILS", NF90_FORMAT,dims,id_BILS),fnam) CALL ncerr(NF90_DEF_VAR(nid,"ALB", NF90_FORMAT,dims,id_ALB),fnam) CALL ncerr(NF90_DEF_VAR(nid,"RUG", NF90_FORMAT,dims,id_RUG),fnam) call nf95_def_var(nid, "longitude", NF90_FLOAT, ndim, varid_longitude) call nf95_def_var(nid, "latitude", NF90_FLOAT, ndim, varid_latitude) !--- Attributes creation CALL ncerr(NF90_PUT_ATT(nid,id_tim, "title","Jour dans l annee"),fnam) CALL ncerr(NF90_PUT_ATT(nid,id_FOCE,"title","Fraction ocean"),fnam) CALL ncerr(NF90_PUT_ATT(nid,id_FSIC,"title","Fraction glace de mer"),fnam) CALL ncerr(NF90_PUT_ATT(nid,id_FTER,"title","Fraction terre"),fnam) CALL ncerr(NF90_PUT_ATT(nid,id_FLIC,"title","Fraction land ice"),fnam) CALL ncerr(NF90_PUT_ATT(nid,id_SST ,"title","Temperature superficielle de la mer"),fnam) CALL ncerr(NF90_PUT_ATT(nid,id_BILS,"title","Reference flux de chaleur au sol"),fnam) CALL ncerr(NF90_PUT_ATT(nid,id_ALB, "title","Albedo a la surface"),fnam) CALL ncerr(NF90_PUT_ATT(nid,id_RUG, "title","Rugosite"),fnam) call nf95_put_att(nid, varid_longitude, "standard_name", "longitude") call nf95_put_att(nid, varid_longitude, "units", "degrees_east") call nf95_put_att(nid, varid_latitude, "standard_name", "latitude") call nf95_put_att(nid, varid_latitude, "units", "degrees_north") CALL ncerr(NF90_ENDDEF(nid),fnam) !--- Variables saving CALL ncerr(NF90_PUT_VAR(nid,id_tim,[(REAL(k),k=1,ndays)]),fnam) CALL ncerr(NF90_PUT_VAR(nid,id_FOCE,pctsrf_t(:,is_oce,:),[1,1],[klon,ndays]),fnam) CALL ncerr(NF90_PUT_VAR(nid,id_FSIC,pctsrf_t(:,is_sic,:),[1,1],[klon,ndays]),fnam) CALL ncerr(NF90_PUT_VAR(nid,id_FTER,pctsrf_t(:,is_ter,:),[1,1],[klon,ndays]),fnam) CALL ncerr(NF90_PUT_VAR(nid,id_FLIC,pctsrf_t(:,is_lic,:),[1,1],[klon,ndays]),fnam) CALL ncerr(NF90_PUT_VAR(nid,id_SST ,phy_sst(:,:),[1,1],[klon,ndays]),fnam) CALL ncerr(NF90_PUT_VAR(nid,id_BILS,phy_bil(:,:),[1,1],[klon,ndays]),fnam) CALL ncerr(NF90_PUT_VAR(nid,id_ALB ,phy_alb(:,:),[1,1],[klon,ndays]),fnam) CALL ncerr(NF90_PUT_VAR(nid,id_RUG ,phy_rug(:,:),[1,1],[klon,ndays]),fnam) call nf95_put_var(nid, varid_longitude, rlon) call nf95_put_var(nid, varid_latitude, rlat) CALL ncerr(NF90_CLOSE(nid),fnam) CALL msg(6,'Ecriture du fichier limit : fin') DEALLOCATE(pctsrf_t,phy_sst,phy_bil,phy_alb,phy_rug) !=============================================================================== ! CONTAINS ! !=============================================================================== !------------------------------------------------------------------------------- ! SUBROUTINE get_2Dfield(fnam, varname, mode, ndays, champo, flag, mask) ! !----------------------------------------------------------------------------- ! Comments: ! There are two assumptions concerning the NetCDF files, that are satisfied ! with files that are conforming NC convention: ! 1) The last dimension of the variables used is the time record. ! 2) Dimensional variables have the same names as corresponding dimensions. !----------------------------------------------------------------------------- USE netcdf, ONLY: NF90_OPEN, NF90_INQ_VARID, NF90_INQUIRE_VARIABLE, & NF90_CLOSE, NF90_INQ_DIMID, NF90_INQUIRE_DIMENSION, NF90_GET_VAR, & NF90_GET_ATT USE pchsp_95_m, only: pchsp_95 USE pchfe_95_m, only: pchfe_95 USE arth_m, only: arth USE indice_sol_mod IMPLICIT NONE include "dimensions.h" include "paramet.h" include "comgeom2.h" !----------------------------------------------------------------------------- ! Arguments: CHARACTER(LEN=*), INTENT(IN) :: fnam ! NetCDF file name CHARACTER(LEN=10), INTENT(IN) :: varname ! NetCDF variable name CHARACTER(LEN=3), INTENT(IN) :: mode ! RUG, SIC, SST or ALB INTEGER, INTENT(IN) :: ndays ! current year number of days REAL, POINTER, DIMENSION(:, :) :: champo ! output field = f(t) LOGICAL, OPTIONAL, INTENT(IN) :: flag ! extrapol. (SST) old ice (SIC) REAL, OPTIONAL, DIMENSION(iim, jjp1), INTENT(IN) :: mask !------------------------------------------------------------------------------ ! Local variables: !--- NetCDF INTEGER :: ncid, varid ! NetCDF identifiers CHARACTER(LEN=30) :: dnam ! dimension name !--- dimensions INTEGER :: dids(4) ! NetCDF dimensions identifiers REAL, ALLOCATABLE :: dlon_ini(:) ! initial longitudes vector REAL, ALLOCATABLE :: dlat_ini(:) ! initial latitudes vector REAL, POINTER :: dlon(:), dlat(:) ! reordered lon/lat vectors !--- fields INTEGER :: imdep, jmdep, lmdep ! dimensions of 'champ' REAL, ALLOCATABLE :: champ(:,:) ! wanted field on initial grid REAL, ALLOCATABLE :: yder(:), timeyear(:) REAL :: champint(iim,jjp1) ! interpolated field REAL, ALLOCATABLE :: champtime(:,:,:) REAL, ALLOCATABLE :: champan(:,:,:) !--- input files CHARACTER(LEN=20) :: cal_in ! calendar CHARACTER(LEN=20) :: unit_sic ! attribute unit in sea-ice file INTEGER :: ndays_in ! number of days !--- misc INTEGER :: i, j, k, l ! loop counters REAL, ALLOCATABLE :: work(:,:) ! used for extrapolation CHARACTER(LEN=25) :: title ! for messages LOGICAL :: extrp ! flag for extrapolation REAL :: chmin, chmax INTEGER ierr integer n_extrap ! number of extrapolated points logical skip !------------------------------------------------------------------------------ !---Variables depending on keyword 'mode' ------------------------------------- NULLIFY(champo) SELECT CASE(mode) CASE('RUG'); title='Rugosite' CASE('SIC'); title='Sea-ice' CASE('SST'); title='SST' CASE('ALB'); title='Albedo' END SELECT extrp=.FALSE.; IF(PRESENT(flag).AND.mode=='SST') extrp=flag !--- GETTING SOME DIMENSIONAL VARIABLES FROM FILE ----------------------------- CALL msg(5,' Now reading file : '//TRIM(fnam)) CALL ncerr(NF90_OPEN(fnam, NF90_NOWRITE, ncid),fnam) CALL ncerr(NF90_INQ_VARID(ncid, trim(varname), varid),fnam) CALL ncerr(NF90_INQUIRE_VARIABLE(ncid, varid, dimids=dids),fnam) !--- Read unit for sea-ice variable only IF (mode=='SIC') THEN IF(NF90_GET_ATT(ncid, varid, 'units', unit_sic)/=NF90_NOERR) THEN CALL msg(5,'No unit in sea-ice file. Take percentage as default value') unit_sic='X' ELSE CALL msg(5,'Sea-ice cover has unit='//TRIM(unit_sic)) END IF END IF !--- Longitude CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(1), name=dnam, len=imdep),fnam) ALLOCATE(dlon_ini(imdep), dlon(imdep)) CALL ncerr(NF90_INQ_VARID(ncid, dnam, varid), fnam) CALL ncerr(NF90_GET_VAR(ncid, varid, dlon_ini), fnam) CALL msg(5,'variable '//TRIM(dnam)//' dimension ', imdep) !--- Latitude CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(2), name=dnam, len=jmdep),fnam) ALLOCATE(dlat_ini(jmdep), dlat(jmdep)) CALL ncerr(NF90_INQ_VARID(ncid, dnam, varid), fnam) CALL ncerr(NF90_GET_VAR(ncid, varid, dlat_ini), fnam) CALL msg(5,'variable '//TRIM(dnam)//' dimension ', jmdep) !--- Time (variable is not needed - it is rebuilt - but calendar is) CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(3), name=dnam, len=lmdep), fnam) ALLOCATE(timeyear(lmdep)) CALL ncerr(NF90_INQ_VARID(ncid, dnam, varid), fnam) cal_in=' ' IF(NF90_GET_ATT(ncid, varid, 'calendar', cal_in)/=NF90_NOERR) THEN SELECT CASE(mode) CASE('RUG', 'ALB'); cal_in='360d' CASE('SIC', 'SST'); cal_in='gregorian' END SELECT CALL msg(5,'WARNING: missing "calendar" attribute for "time" in '& &//TRIM(fnam)//'. Choosing default value.') END IF CALL msg(5,'var, calendar, dim: '//TRIM(dnam)//' '//TRIM(cal_in), lmdep) !--- CONSTRUCTING THE INPUT TIME VECTOR FOR INTERPOLATION -------------------- !--- Determining input file number of days, depending on calendar ndays_in=year_len(anneeref, cal_in) !--- Time vector reconstruction (time vector from file is not trusted) !--- If input records are not monthly, time sampling has to be constant ! timeyear=mid_months(anneeref, cal_in, lmdep) IF (lmdep /= 12) WRITE(lunout,*) 'Note : les fichiers de ', TRIM(mode), & ' ne comportent pas 12, mais ', lmdep, ' enregistrements.' !--- GETTING THE FIELD AND INTERPOLATING IT ---------------------------------- ALLOCATE(champ(imdep, jmdep), champtime(iim, jjp1, lmdep)) IF(extrp) ALLOCATE(work(imdep, jmdep)) CALL msg(5,'') CALL msg(5,'LECTURE ET INTERPOLATION HORIZ. DE ', lmdep, ' CHAMPS.') CALL ncerr(NF90_INQ_VARID(ncid, varname, varid), fnam) DO l=1, lmdep CALL ncerr(NF90_GET_VAR(ncid,varid,champ,[1,1,l],[imdep,jmdep,1]),fnam) CALL conf_dat2d(title, dlon_ini, dlat_ini, dlon, dlat, champ, .TRUE.) IF(extrp) CALL extrapol(champ,imdep,jmdep,999999.,.TRUE.,.TRUE.,2,work) IF(l==1) THEN CALL msg(5,"----------------------------------------------------------") CALL msg(5,"$$$ Interpolation barycentrique pour "//TRIM(title)//" $$$") CALL msg(5,"----------------------------------------------------------") END IF IF(mode=='RUG') champ=LOG(champ) CALL inter_barxy(dlon,dlat(:jmdep-1),champ,rlonu(:iim),rlatv,champint) IF(mode=='RUG') THEN champint=EXP(champint) WHERE(NINT(mask)/=1) champint=0.001 END IF champtime(:, :, l)=champint END DO CALL ncerr(NF90_CLOSE(ncid), fnam) DEALLOCATE(dlon_ini, dlat_ini, dlon, dlat, champ) IF(extrp) DEALLOCATE(work) !--- TIME INTERPOLATION ------------------------------------------------------ IF(prt_level>5) THEN WRITE(lunout, *) WRITE(lunout, *)'INTERPOLATION TEMPORELLE.' WRITE(lunout, *)' Vecteur temps en entree: ', timeyear WRITE(lunout, *)' Vecteur temps en sortie de 0 a ', ndays END IF ALLOCATE(yder(lmdep), champan(iip1, jjp1, ndays)) skip = .false. n_extrap = 0 DO j=1, jjp1 DO i=1, iim yder = pchsp_95(timeyear, champtime(i, j, :), ibeg=2, iend=2, & vc_beg=0., vc_end=0.) CALL pchfe_95(timeyear, champtime(i, j, :), yder, skip, & arth(0., real(ndays_in) / ndays, ndays), champan(i, j, :), ierr) if (ierr < 0) stop 1 n_extrap = n_extrap + ierr END DO END DO if (n_extrap /= 0) then WRITE(lunout,*) "get_2Dfield pchfe_95: n_extrap = ", n_extrap end if champan(iip1, :, :)=champan(1, :, :) DEALLOCATE(yder, champtime, timeyear) !--- Checking the result DO j=1, jjp1 CALL minmax(iip1, champan(1, j, 10), chmin, chmax) IF (prt_level>5) WRITE(lunout, *)' ',TRIM(title),' au temps 10 ', chmin, chmax, j END DO !--- SPECIAL FILTER FOR SST: SST>271.38 -------------------------------------- IF(mode=='SST') THEN CALL msg(5,'Filtrage de la SST: SST >= 271.38') WHERE(champan<271.38) champan=271.38 END IF !--- SPECIAL FILTER FOR SIC: 0.01.0) champan=1.0 WHERE(champan<0.0) champan=0.0 END IF !--- DYNAMICAL TO PHYSICAL GRID ---------------------------------------------- ALLOCATE(champo(klon, ndays)) DO k=1, ndays CALL gr_dyn_fi(1, iip1, jjp1, klon, champan(1, 1, k), champo(1, k)) END DO DEALLOCATE(champan) END SUBROUTINE get_2Dfield ! !------------------------------------------------------------------------------- !------------------------------------------------------------------------------- ! SUBROUTINE start_init_orog0(lon_in,lat_in,phis,masque) ! !------------------------------------------------------------------------------- IMPLICIT NONE !=============================================================================== ! Purpose: Compute "phis" just like it would be in start_init_orog. !=============================================================================== ! Arguments: REAL, INTENT(IN) :: lon_in(:), lat_in(:) ! dim (iml) (jml) REAL, INTENT(INOUT) :: phis(:,:), masque(:,:) ! dim (iml,jml) !------------------------------------------------------------------------------- ! Local variables: CHARACTER(LEN=256) :: modname="start_init_orog0" INTEGER :: fid, llm_tmp,ttm_tmp, iml,jml, iml_rel,jml_rel, itau(1) REAL :: lev(1), date, dt, deg2rad REAL, ALLOCATABLE :: lon_rad(:), lon_ini(:), lon_rel(:,:), relief_hi(:,:) REAL, ALLOCATABLE :: lat_rad(:), lat_ini(:), lat_rel(:,:) !------------------------------------------------------------------------------- iml=assert_eq(SIZE(lon_in),SIZE(phis,1),SIZE(masque,1),TRIM(modname)//" iml") jml=assert_eq(SIZE(lat_in),SIZE(phis,2),SIZE(masque,2),TRIM(modname)//" jml") IF(iml/=iip1) CALL abort_gcm(TRIM(modname),'iml/=iip1',1) IF(jml/=jjp1) CALL abort_gcm(TRIM(modname),'jml/=jjp1',1) pi=2.0*ASIN(1.0); deg2rad=pi/180.0 IF(ANY(phis/=-99999.)) RETURN !--- phis ALREADY KNOWN !--- HIGH RESOLUTION OROGRAPHY CALL flininfo(frelf, iml_rel, jml_rel, llm_tmp, ttm_tmp, fid) ALLOCATE(lat_rel(iml_rel,jml_rel),lon_rel(iml_rel,jml_rel)) CALL flinopen(frelf, .FALSE., iml_rel, jml_rel, llm_tmp, lon_rel, lat_rel, & lev, ttm_tmp, itau, date, dt, fid) ALLOCATE(relief_hi(iml_rel,jml_rel)) CALL flinget(fid, vrel, iml_rel, jml_rel, llm_tmp, ttm_tmp, 1, 1, relief_hi) CALL flinclo(fid) !--- IF ANGLES ARE IN DEGREES, THEY ARE CONVERTED INTO RADIANS ALLOCATE(lon_ini(iml_rel),lat_ini(jml_rel)) lon_ini(:)=lon_rel(:,1); IF(MAXVAL(lon_rel)>pi) lon_ini=lon_ini*deg2rad lat_ini(:)=lat_rel(1,:); IF(MAXVAL(lat_rel)>pi) lat_ini=lat_ini*deg2rad !--- FIELDS ARE PROCESSED TO BE ON STANDARD ANGULAR DOMAINS ALLOCATE(lon_rad(iml_rel),lat_rad(jml_rel)) CALL conf_dat2d(vrel, lon_ini, lat_ini, lon_rad, lat_rad, relief_hi, .FALSE.) DEALLOCATE(lon_ini,lat_ini) !--- COMPUTING SURFACE GEOPOTENTIAL USING ROUTINE grid_noro0 WRITE(lunout,*) WRITE(lunout,*)'*** Compute surface geopotential ***' !--- CALL OROGRAPHY MODULE (REDUCED VERSION) TO COMPUTE FIELDS CALL grid_noro0(lon_rad, lat_rad, relief_hi, lon_in, lat_in, phis, masque) phis = phis * 9.81 phis(iml,:) = phis(1,:) DEALLOCATE(relief_hi,lon_rad,lat_rad) END SUBROUTINE start_init_orog0 ! !------------------------------------------------------------------------------- !------------------------------------------------------------------------------- ! SUBROUTINE grid_noro0(xd,yd,zd,x,y,zphi,mask) ! !=============================================================================== ! Purpose: Extracted from grid_noro to provide geopotential height for dynamics ! without any call to physics subroutines. !=============================================================================== IMPLICIT NONE !------------------------------------------------------------------------------- ! Arguments: REAL, INTENT(IN) :: xd(:), yd(:) !--- INPUT COORDINATES (imdp) (jmdp) REAL, INTENT(IN) :: zd(:,:) !--- INPUT FIELD (imdp,jmdp) REAL, INTENT(IN) :: x(:), y(:) !--- OUTPUT COORDINATES (imar+1) (jmar) REAL, INTENT(OUT) :: zphi(:,:) !--- GEOPOTENTIAL (imar+1,jmar) REAL, INTENT(INOUT):: mask(:,:) !--- MASK (imar+1,jmar) !------------------------------------------------------------------------------- ! Local variables: CHARACTER(LEN=256) :: modname="grid_noro0" REAL, ALLOCATABLE :: xusn(:), yusn(:) ! dim (imdp+2*iext) (jmdp+2) REAL, ALLOCATABLE :: zusn(:,:) ! dim (imdp+2*iext,jmdp+2) REAL, ALLOCATABLE :: weight(:,:) ! dim (imar+1,jmar) REAL, ALLOCATABLE :: mask_tmp(:,:), zmea(:,:) ! dim (imar+1,jmar) REAL, ALLOCATABLE :: num_tot(:,:), num_lan(:,:) ! dim (imax,jmax) REAL, ALLOCATABLE :: a(:), b(:) ! dim (imax) REAL, ALLOCATABLE :: c(:), d(:) ! dim (jmax) LOGICAL :: masque_lu INTEGER :: i, ii, imdp, imar, iext INTEGER :: j, jj, jmdp, jmar, nn REAL :: xpi, zlenx, weighx, xincr, zbordnor, zmeanor, zweinor, zbordest REAL :: rad, zleny, weighy, masque, zbordsud, zmeasud, zweisud, zbordoue !------------------------------------------------------------------------------- imdp=assert_eq(SIZE(xd),SIZE(zd,1),TRIM(modname)//" imdp") jmdp=assert_eq(SIZE(yd),SIZE(zd,2),TRIM(modname)//" jmdp") imar=assert_eq(SIZE(x),SIZE(zphi,1),SIZE(mask,1),TRIM(modname)//" imar")-1 jmar=assert_eq(SIZE(y),SIZE(zphi,2),SIZE(mask,2),TRIM(modname)//" jmar") IF(imar/=iim) CALL abort_gcm(TRIM(modname),'imar/=iim' ,1) IF(jmar/=jjm+1) CALL abort_gcm(TRIM(modname),'jmar/=jjm+1',1) iext=imdp/10 xpi = ACOS(-1.) rad = 6371229. !--- ARE WE USING A READ MASK ? masque_lu=ANY(mask/=-99999.); IF(.NOT.masque_lu) mask=0.0 WRITE(lunout,*)'Masque lu: ',masque_lu !--- EXTENSION OF THE INPUT DATABASE TO PROCEED COMPUTATIONS AT BOUNDARIES: ALLOCATE(xusn(imdp+2*iext)) xusn(1 +iext:imdp +iext)=xd(:) xusn(1 : iext)=xd(1+imdp-iext:imdp)-2.*xpi xusn(1+imdp+iext:imdp+2*iext)=xd(1 :iext)+2.*xpi ALLOCATE(yusn(jmdp+2)) yusn(1 )=yd(1) +(yd(1) -yd(2)) yusn(2:jmdp+1)=yd(:) yusn( jmdp+2)=yd(jmdp)+(yd(jmdp)-yd(jmdp-1)) ALLOCATE(zusn(imdp+2*iext,jmdp+2)) zusn(1 +iext:imdp +iext,2:jmdp+1)=zd (: , :) zusn(1 : iext,2:jmdp+1)=zd (imdp-iext+1:imdp , :) zusn(1+imdp +iext:imdp+2*iext,2:jmdp+1)=zd (1:iext , :) zusn(1 :imdp/2+iext, 1)=zusn(1+imdp/2:imdp +iext, 2) zusn(1+imdp/2+iext:imdp+2*iext, 1)=zusn(1 :imdp/2+iext, 2) zusn(1 :imdp/2+iext, jmdp+2)=zusn(1+imdp/2:imdp +iext,jmdp+1) zusn(1+imdp/2+iext:imdp+2*iext, jmdp+2)=zusn(1 :imdp/2+iext,jmdp+1) !--- COMPUTE LIMITS OF MODEL GRIDPOINT AREA (REGULAR GRID) ALLOCATE(a(imar+1),b(imar+1)) b(1:imar)=(x(1:imar )+ x(2:imar+1))/2.0 b(imar+1)= x( imar+1)+(x( imar+1)-x(imar))/2.0 a(1)=x(1)-(x(2)-x(1))/2.0 a(2:imar+1)= b(1:imar) ALLOCATE(c(jmar),d(jmar)) d(1:jmar-1)=(y(1:jmar-1)+ y(2:jmar))/2.0 d( jmar )= y( jmar )+(y( jmar)-y(jmar-1))/2.0 c(1)=y(1)-(y(2)-y(1))/2.0 c(2:jmar)=d(1:jmar-1) !--- INITIALIZATIONS: ALLOCATE(weight(imar+1,jmar)); weight(:,:)= 0.0 ALLOCATE(zmea (imar+1,jmar)); zmea (:,:)= 0.0 !--- SUMMATION OVER GRIDPOINT AREA zleny=xpi/REAL(jmdp)*rad xincr=xpi/REAL(jmdp)/2. ALLOCATE(num_tot(imar+1,jmar)); num_tot(:,:)=0. ALLOCATE(num_lan(imar+1,jmar)); num_lan(:,:)=0. DO ii = 1, imar+1 DO jj = 1, jmar DO j = 2,jmdp+1 zlenx =zleny *COS(yusn(j)) zbordnor=(xincr+c(jj)-yusn(j))*rad zbordsud=(xincr-d(jj)+yusn(j))*rad weighy=AMAX1(0.,AMIN1(zbordnor,zbordsud,zleny)) IF(weighy/=0) THEN DO i = 2, imdp+2*iext-1 zbordest=(xusn(i)-a(ii)+xincr)*rad*COS(yusn(j)) zbordoue=(b(ii)+xincr-xusn(i))*rad*COS(yusn(j)) weighx=AMAX1(0.,AMIN1(zbordest,zbordoue,zlenx)) IF(weighx/=0)THEN num_tot(ii,jj)=num_tot(ii,jj)+1.0 IF(zusn(i,j)>=1.)num_lan(ii,jj)=num_lan(ii,jj)+1.0 weight(ii,jj)=weight(ii,jj)+weighx*weighy zmea (ii,jj)=zmea (ii,jj)+zusn(i,j)*weighx*weighy !--- MEAN END IF END DO END IF END DO END DO END DO !--- COMPUTE PARAMETERS NEEDED BY LOTT & MILLER (1997) AND LOTT (1999) SSO SCHEME IF(.NOT.masque_lu) THEN WHERE(weight(:,1:jmar-1)/=0.0) mask=num_lan(:,:)/num_tot(:,:) END IF nn=COUNT(weight(:,1:jmar-1)==0.0) IF(nn/=0) WRITE(lunout,*)'Problem with weight ; vanishing occurrences: ',nn WHERE(weight/=0.0) zmea(:,:)=zmea(:,:)/weight(:,:) !--- MASK BASED ON GROUND MAXIMUM, 10% THRESHOLD (<10%: SURF PARAMS MEANINGLESS) ALLOCATE(mask_tmp(imar+1,jmar)); mask_tmp(:,:)=0.0 WHERE(mask>=0.1) mask_tmp = 1. WHERE(weight(:,:)/=0.0) zphi(:,:)=mask_tmp(:,:)*zmea(:,:) zmea(:,:)=mask_tmp(:,:)*zmea(:,:) END WHERE WRITE(lunout,*)' MEAN ORO:' ,MAXVAL(zmea) !--- Values at poles zphi(imar+1,:)=zphi(1,:) zweinor=SUM(weight(1:imar, 1),DIM=1) zweisud=SUM(weight(1:imar,jmar),DIM=1) zmeanor=SUM(weight(1:imar, 1)*zmea(1:imar, 1),DIM=1) zmeasud=SUM(weight(1:imar,jmar)*zmea(1:imar,jmar),DIM=1) zphi(:,1)=zmeanor/zweinor; zphi(:,jmar)=zmeasud/zweisud END SUBROUTINE grid_noro0 ! !------------------------------------------------------------------------------- !------------------------------------------------------------------------------- ! FUNCTION year_len(y,cal_in) ! !------------------------------------------------------------------------------- IMPLICIT NONE !------------------------------------------------------------------------------- ! Arguments: INTEGER :: year_len INTEGER, INTENT(IN) :: y CHARACTER(LEN=*), INTENT(IN) :: cal_in !------------------------------------------------------------------------------- ! Local variables: CHARACTER(LEN=20) :: cal_out ! calendar (for outputs) !------------------------------------------------------------------------------- !--- Getting the input calendar to reset at the end of the function CALL ioget_calendar(cal_out) !--- Unlocking calendar and setting it to wanted one CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_in)) !--- Getting the number of days in this year year_len=ioget_year_len(y) !--- Back to original calendar CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_out)) END FUNCTION year_len ! !------------------------------------------------------------------------------- !------------------------------------------------------------------------------- ! FUNCTION mid_months(y,cal_in,nm) ! !------------------------------------------------------------------------------- IMPLICIT NONE !------------------------------------------------------------------------------- ! Arguments: INTEGER, INTENT(IN) :: y ! year CHARACTER(LEN=*), INTENT(IN) :: cal_in ! calendar INTEGER, INTENT(IN) :: nm ! months/year number REAL, DIMENSION(nm) :: mid_months ! mid-month times !------------------------------------------------------------------------------- ! Local variables: CHARACTER(LEN=99) :: mess ! error message CHARACTER(LEN=20) :: cal_out ! calendar (for outputs) INTEGER, DIMENSION(nm) :: mnth ! months lengths (days) INTEGER :: m ! months counter INTEGER :: nd ! number of days !------------------------------------------------------------------------------- nd=year_len(y,cal_in) IF(nm==12) THEN !--- Getting the input calendar to reset at the end of the function CALL ioget_calendar(cal_out) !--- Unlocking calendar and setting it to wanted one CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_in)) !--- Getting the length of each month DO m=1,nm; mnth(m)=ioget_mon_len(y,m); END DO !--- Back to original calendar CALL lock_calendar(.FALSE.); CALL ioconf_calendar(TRIM(cal_out)) ELSE IF(MODULO(nd,nm)/=0) THEN WRITE(mess,'(a,i3,a,i3,a)')'Unconsistent calendar: ',nd,' days/year, but ',& nm,' months/year. Months number should divide days number.' CALL abort_physic('mid_months',TRIM(mess),1) ELSE mnth=[(m,m=1,nm,nd/nm)] END IF !--- Mid-months times mid_months(1)=0.5*REAL(mnth(1)) DO k=2,nm mid_months(k)=mid_months(k-1)+0.5*REAL(mnth(k-1)+mnth(k)) END DO END FUNCTION mid_months ! !------------------------------------------------------------------------------- !------------------------------------------------------------------------------- ! SUBROUTINE msg(lev,str1,i,str2) ! !------------------------------------------------------------------------------- ! Arguments: INTEGER, INTENT(IN) :: lev CHARACTER(LEN=*), INTENT(IN) :: str1 INTEGER, OPTIONAL, INTENT(IN) :: i CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: str2 !------------------------------------------------------------------------------- IF(prt_level>lev) THEN IF(PRESENT(str2)) THEN WRITE(lunout,*) TRIM(str1), i, TRIM(str2) ELSE IF(PRESENT(i)) THEN WRITE(lunout,*) TRIM(str1), i ELSE WRITE(lunout,*) TRIM(str1) END IF END IF END SUBROUTINE msg ! !------------------------------------------------------------------------------- !------------------------------------------------------------------------------- ! SUBROUTINE ncerr(ncres,fnam) ! !------------------------------------------------------------------------------- ! Purpose: NetCDF errors handling. !------------------------------------------------------------------------------- USE netcdf, ONLY : NF90_NOERR, NF90_STRERROR IMPLICIT NONE !------------------------------------------------------------------------------- ! Arguments: INTEGER, INTENT(IN) :: ncres CHARACTER(LEN=*), INTENT(IN) :: fnam !------------------------------------------------------------------------------- IF(ncres/=NF90_NOERR) THEN WRITE(lunout,*)'Problem with file '//TRIM(fnam)//' in routine limit_netcdf.' CALL abort_physic('limit_netcdf',NF90_STRERROR(ncres),1) END IF END SUBROUTINE ncerr ! !------------------------------------------------------------------------------- #endif ! of #ifndef CPP_1D END SUBROUTINE limit_netcdf END MODULE limit ! !*******************************************************************************