[2336] | 1 | MODULE limit |
---|
[1687] | 2 | ! |
---|
[2336] | 3 | !******************************************************************************* |
---|
| 4 | ! Author : L. Fairhead, 27/01/94 |
---|
[1687] | 5 | !------------------------------------------------------------------------------- |
---|
[2336] | 6 | ! Purpose: Boundary conditions files building for new model using climatologies. |
---|
| 7 | ! Both grids have to be regular. |
---|
| 8 | !------------------------------------------------------------------------------- |
---|
| 9 | ! Note: This routine is designed to work for Earth |
---|
| 10 | !------------------------------------------------------------------------------- |
---|
| 11 | ! Modification history: |
---|
| 12 | ! * 23/03/1994: Z. X. Li |
---|
| 13 | ! * 09/1999: L. Fairhead (netcdf reading in LMDZ.3.3) |
---|
| 14 | ! * 07/2001: P. Le Van |
---|
| 15 | ! * 11/2009: L. Guez (ozone day & night climatos, see etat0_netcdf.F90) |
---|
| 16 | ! * 12/2009: D. Cugnet (f77->f90, calendars, files from coupled runs) |
---|
| 17 | !------------------------------------------------------------------------------- |
---|
| 18 | |
---|
[2766] | 19 | USE ioipsl, ONLY: flininfo, flinopen, flinget, flinclo |
---|
[2336] | 20 | USE assert_eq_m, ONLY: assert_eq |
---|
[2766] | 21 | USE cal_tools_m, ONLY: year_len, mid_month |
---|
[2336] | 22 | USE conf_dat_m, ONLY: conf_dat2d, conf_dat3d |
---|
| 23 | USE dimphy, ONLY: klon, zmasq |
---|
[2766] | 24 | USE geometry_mod, ONLY: longitude_deg, latitude_deg |
---|
[2399] | 25 | USE phys_state_var_mod, ONLY: pctsrf |
---|
[2766] | 26 | USE control_mod, ONLY: anneeref |
---|
[2336] | 27 | USE init_ssrf_m, ONLY: start_init_subsurf |
---|
| 28 | |
---|
[3163] | 29 | INTEGER, PARAMETER :: ns=256 |
---|
| 30 | CHARACTER(LEN=ns), PARAMETER :: & |
---|
[2889] | 31 | fsst(5)=['amipbc_sst_1x1.nc ','amip_sst_1x1.nc ','cpl_atm_sst.nc '& |
---|
[3163] | 32 | ,'histmth_sst.nc ','sstk.nc '], & |
---|
[2889] | 33 | fsic(5)=['amipbc_sic_1x1.nc ','amip_sic_1x1.nc ','cpl_atm_sic.nc '& |
---|
[3163] | 34 | ,'histmth_sic.nc ','ci.nc '], & |
---|
[2889] | 35 | vsst(5)=['tosbcs ','tos ','SISUTESW ','tsol_oce ','sstk '], & |
---|
[3163] | 36 | vsic(5)=['sicbcs ','sic ','SIICECOV ','pourc_sic ','ci '], & |
---|
| 37 | frugo='Rugos.nc ', falbe='Albedo.nc ', frelf='Relief.nc ', & |
---|
| 38 | vrug='RUGOS ', valb='ALBEDO ', vrel='RELIEF ', & |
---|
| 39 | DegK(11)=['degK ','degree_K ','degreeK ','deg_K '& |
---|
| 40 | ,'degsK ','degrees_K ','degreesK ','degs_K '& |
---|
| 41 | ,'degree_kelvin ','degrees_kelvin','K '], & |
---|
| 42 | DegC(10)=['degC ','degree_C ','degreeC ','deg_C '& |
---|
| 43 | ,'degsC ','degrees_C ','degreesC ','degs_C '& |
---|
| 44 | ,'degree_Celsius','celsius '], & |
---|
[3168] | 45 | Perc(2) =['% ','percent '], & |
---|
| 46 | Frac(2) =['1.0 ','1 '] |
---|
[2336] | 47 | |
---|
| 48 | CONTAINS |
---|
| 49 | |
---|
| 50 | !------------------------------------------------------------------------------- |
---|
| 51 | ! |
---|
| 52 | SUBROUTINE limit_netcdf(masque, phis, extrap) |
---|
| 53 | ! |
---|
| 54 | !------------------------------------------------------------------------------- |
---|
[1687] | 55 | ! Author : L. Fairhead, 27/01/94 |
---|
| 56 | !------------------------------------------------------------------------------- |
---|
| 57 | ! Purpose: Boundary conditions files building for new model using climatologies. |
---|
| 58 | ! Both grids have to be regular. |
---|
| 59 | !------------------------------------------------------------------------------- |
---|
| 60 | ! Note: This routine is designed to work for Earth |
---|
| 61 | !------------------------------------------------------------------------------- |
---|
| 62 | ! Modification history: |
---|
| 63 | ! * 23/03/1994: Z. X. Li |
---|
| 64 | ! * 09/1999: L. Fairhead (netcdf reading in LMDZ.3.3) |
---|
| 65 | ! * 07/2001: P. Le Van |
---|
| 66 | ! * 11/2009: L. Guez (ozone day & night climatos, see etat0_netcdf.F90) |
---|
| 67 | ! * 12/2009: D. Cugnet (f77->f90, calendars, files from coupled runs) |
---|
[2540] | 68 | ! * 04/2016: D. Cugnet (12/14 recs SST/SIC files: cyclic/interannual runs) |
---|
[2889] | 69 | ! * 05/2017: D. Cugnet (linear time interpolation for BCS files) |
---|
[1687] | 70 | !------------------------------------------------------------------------------- |
---|
[2293] | 71 | #ifndef CPP_1D |
---|
[5282] | 72 | USE iniprint_mod_h |
---|
[5281] | 73 | USE comgeom2_mod_h |
---|
[1785] | 74 | USE indice_sol_mod |
---|
[5084] | 75 | USE netcdf, ONLY: NF90_OPEN, NF90_CREATE, NF90_CLOSE, & |
---|
[2336] | 76 | NF90_DEF_DIM, NF90_DEF_VAR, NF90_PUT_VAR, NF90_PUT_ATT, & |
---|
[5084] | 77 | NF90_NOERR, NF90_NOWRITE, NF90_DOUBLE, NF90_GLOBAL, & |
---|
[3803] | 78 | NF90_CLOBBER, NF90_ENDDEF, NF90_UNLIMITED, NF90_FLOAT, & |
---|
[5084] | 79 | NF90_64BIT_OFFSET |
---|
[2336] | 80 | USE inter_barxy_m, ONLY: inter_barxy |
---|
| 81 | USE netcdf95, ONLY: nf95_def_var, nf95_put_att, nf95_put_var |
---|
[2597] | 82 | USE comconst_mod, ONLY: pi |
---|
[2665] | 83 | USE phys_cal_mod, ONLY: calend |
---|
[5249] | 84 | USE lmdz_cppkeys_wrapper, ONLY: nf90_format |
---|
[5271] | 85 | USE dimensions_mod, ONLY: iim, jjm, llm, ndm |
---|
[5285] | 86 | USE paramet_mod_h |
---|
[5281] | 87 | USE comgeom2_mod_h |
---|
| 88 | IMPLICIT NONE |
---|
[1687] | 89 | !------------------------------------------------------------------------------- |
---|
| 90 | ! Arguments: |
---|
[5271] | 91 | |
---|
[5272] | 92 | |
---|
[2336] | 93 | REAL, DIMENSION(iip1,jjp1), INTENT(INOUT) :: masque ! land mask |
---|
| 94 | REAL, DIMENSION(iip1,jjp1), INTENT(INOUT) :: phis ! ground geopotential |
---|
| 95 | LOGICAL, INTENT(IN) :: extrap ! SST extrapolation flag |
---|
[1687] | 96 | |
---|
[3168] | 97 | !--- INPUT NETCDF FILES AND VARIABLES NAMES ------------------------------------ |
---|
| 98 | CHARACTER(LEN=ns) :: icefile, sstfile, fnam, varname |
---|
[2336] | 99 | |
---|
[1687] | 100 | !--- OUTPUT VARIABLES FOR NETCDF FILE ------------------------------------------ |
---|
[3168] | 101 | REAL :: fi_ice(klon) |
---|
[2293] | 102 | REAL, POINTER :: phy_rug(:,:)=>NULL(), phy_ice(:,:)=>NULL() |
---|
| 103 | REAL, POINTER :: phy_sst(:,:)=>NULL(), phy_alb(:,:)=>NULL() |
---|
| 104 | REAL, ALLOCATABLE :: phy_bil(:,:), pctsrf_t(:,:,:) |
---|
| 105 | INTEGER :: nbad |
---|
[1687] | 106 | |
---|
| 107 | !--- VARIABLES FOR OUTPUT FILE WRITING ----------------------------------------- |
---|
[3168] | 108 | INTEGER :: nid, ndim, ntim, k, dims(2), ix_sic, ix_sst |
---|
[1687] | 109 | INTEGER :: id_tim, id_SST, id_BILS, id_RUG, id_ALB |
---|
[2128] | 110 | INTEGER :: id_FOCE, id_FSIC, id_FTER, id_FLIC, varid_longitude, varid_latitude |
---|
[1687] | 111 | INTEGER :: ndays !--- Depending on the output calendar |
---|
[3168] | 112 | CHARACTER(LEN=ns) :: str |
---|
[1687] | 113 | |
---|
| 114 | !--- INITIALIZATIONS ----------------------------------------------------------- |
---|
| 115 | CALL inigeom |
---|
| 116 | |
---|
[2336] | 117 | !--- MASK, GROUND GEOPOT. & SUBSURFACES COMPUTATION (IN CASE ok_etat0==.FALSE.) |
---|
| 118 | IF(ALL(masque==-99999.)) THEN |
---|
| 119 | CALL start_init_orog0(rlonv,rlatu,phis,masque) |
---|
| 120 | CALL gr_dyn_fi(1,iip1,jjp1,klon,masque,zmasq) !--- To physical grid |
---|
[2399] | 121 | ALLOCATE(pctsrf(klon,nbsrf)) |
---|
[2336] | 122 | CALL start_init_subsurf(.FALSE.) |
---|
[2576] | 123 | !--- TO MATCH EXACTLY WHAT WOULD BE DONE IN etat0phys_netcdf |
---|
| 124 | WHERE( masque(:,:)<EPSFRA) masque(:,:)=0. |
---|
| 125 | WHERE(1.-masque(:,:)<EPSFRA) masque(:,:)=1. |
---|
[2336] | 126 | END IF |
---|
| 127 | |
---|
[1687] | 128 | !--- Beware: anneeref (from gcm.def) is used to determine output time sampling |
---|
[2766] | 129 | ndays=year_len(anneeref) |
---|
[1687] | 130 | |
---|
| 131 | !--- RUGOSITY TREATMENT -------------------------------------------------------- |
---|
[2947] | 132 | CALL msg(0,""); CALL msg(0," *** TRAITEMENT DE LA RUGOSITE ***") |
---|
[2336] | 133 | CALL get_2Dfield(frugo,vrug,'RUG',ndays,phy_rug,mask=masque(1:iim,:)) |
---|
[1687] | 134 | |
---|
| 135 | !--- OCEAN TREATMENT ----------------------------------------------------------- |
---|
[2947] | 136 | CALL msg(0,""); CALL msg(0," *** TRAITEMENT DE LA GLACE OCEANIQUE ***") |
---|
[1687] | 137 | |
---|
| 138 | ! Input SIC file selection |
---|
| 139 | ! Open file only to test if available |
---|
[2293] | 140 | DO ix_sic=1,SIZE(fsic) |
---|
| 141 | IF ( NF90_OPEN(TRIM(fsic(ix_sic)),NF90_NOWRITE,nid)==NF90_NOERR ) THEN |
---|
| 142 | icefile=fsic(ix_sic); varname=vsic(ix_sic); EXIT |
---|
| 143 | END IF |
---|
| 144 | END DO |
---|
| 145 | IF(ix_sic==SIZE(fsic)+1) THEN |
---|
[1687] | 146 | WRITE(lunout,*) 'ERROR! No sea-ice input file was found.' |
---|
[2293] | 147 | WRITE(lunout,*) 'One of following files must be available : ' |
---|
| 148 | DO k=1,SIZE(fsic); WRITE(lunout,*) TRIM(fsic(k)); END DO |
---|
[2311] | 149 | CALL abort_physic('limit_netcdf','No sea-ice file was found',1) |
---|
[1687] | 150 | END IF |
---|
[2293] | 151 | CALL ncerr(NF90_CLOSE(nid),icefile) |
---|
[2947] | 152 | CALL msg(0,'Fichier choisi pour la glace de mer:'//TRIM(icefile)) |
---|
[1687] | 153 | |
---|
[2336] | 154 | CALL get_2Dfield(icefile,varname, 'SIC',ndays,phy_ice) |
---|
[1687] | 155 | |
---|
| 156 | ALLOCATE(pctsrf_t(klon,nbsrf,ndays)) |
---|
| 157 | DO k=1,ndays |
---|
| 158 | fi_ice=phy_ice(:,k) |
---|
| 159 | WHERE(fi_ice>=1.0 ) fi_ice=1.0 |
---|
| 160 | WHERE(fi_ice<EPSFRA) fi_ice=0.0 |
---|
| 161 | pctsrf_t(:,is_ter,k)=pctsrf(:,is_ter) ! land soil |
---|
| 162 | pctsrf_t(:,is_lic,k)=pctsrf(:,is_lic) ! land ice |
---|
[2293] | 163 | SELECT CASE(ix_sic) |
---|
[2889] | 164 | CASE(3) ! SIC=pICE*(1-LIC-TER) (CPL) |
---|
[1687] | 165 | pctsrf_t(:,is_sic,k)=fi_ice(:)*(1.-pctsrf(:,is_lic)-pctsrf(:,is_ter)) |
---|
[2889] | 166 | CASE(4) ! SIC=pICE (HIST) |
---|
[1687] | 167 | pctsrf_t(:,is_sic,k)=fi_ice(:) |
---|
[2293] | 168 | CASE DEFAULT ! SIC=pICE-LIC (AMIP,ERAI) |
---|
[1687] | 169 | pctsrf_t(:,is_sic,k)=fi_ice-pctsrf_t(:,is_lic,k) |
---|
[2293] | 170 | END SELECT |
---|
[1687] | 171 | WHERE(pctsrf_t(:,is_sic,k)<=0) pctsrf_t(:,is_sic,k)=0. |
---|
| 172 | WHERE(1.0-zmasq<EPSFRA) |
---|
| 173 | pctsrf_t(:,is_sic,k)=0.0 |
---|
| 174 | pctsrf_t(:,is_oce,k)=0.0 |
---|
| 175 | ELSEWHERE |
---|
| 176 | WHERE(pctsrf_t(:,is_sic,k)>=1.0-zmasq) |
---|
| 177 | pctsrf_t(:,is_sic,k)=1.0-zmasq |
---|
| 178 | pctsrf_t(:,is_oce,k)=0.0 |
---|
| 179 | ELSEWHERE |
---|
| 180 | pctsrf_t(:,is_oce,k)=1.0-zmasq-pctsrf_t(:,is_sic,k) |
---|
| 181 | WHERE(pctsrf_t(:,is_oce,k)<EPSFRA) |
---|
| 182 | pctsrf_t(:,is_oce,k)=0.0 |
---|
| 183 | pctsrf_t(:,is_sic,k)=1.0-zmasq |
---|
| 184 | END WHERE |
---|
| 185 | END WHERE |
---|
| 186 | END WHERE |
---|
| 187 | nbad=COUNT(pctsrf_t(:,is_oce,k)<0.0) |
---|
[2293] | 188 | IF(nbad>0) WRITE(lunout,*) 'pb sous maille pour nb points = ',nbad |
---|
| 189 | nbad=COUNT(ABS(SUM(pctsrf_t(:,:,k),DIM=2)-1.0)>EPSFRA) |
---|
[1687] | 190 | IF(nbad>0) WRITE(lunout,*) 'pb sous surface pour nb points = ',nbad |
---|
| 191 | END DO |
---|
| 192 | DEALLOCATE(phy_ice) |
---|
| 193 | |
---|
| 194 | !--- SST TREATMENT ------------------------------------------------------------- |
---|
[2947] | 195 | CALL msg(0,""); CALL msg(0," *** TRAITEMENT DE LA SST ***") |
---|
[1687] | 196 | |
---|
| 197 | ! Input SST file selection |
---|
| 198 | ! Open file only to test if available |
---|
[2293] | 199 | DO ix_sst=1,SIZE(fsst) |
---|
| 200 | IF ( NF90_OPEN(TRIM(fsst(ix_sst)),NF90_NOWRITE,nid)==NF90_NOERR ) THEN |
---|
| 201 | sstfile=fsst(ix_sst); varname=vsst(ix_sst); EXIT |
---|
| 202 | END IF |
---|
| 203 | END DO |
---|
| 204 | IF(ix_sst==SIZE(fsst)+1) THEN |
---|
[1687] | 205 | WRITE(lunout,*) 'ERROR! No sst input file was found.' |
---|
[2293] | 206 | WRITE(lunout,*) 'One of following files must be available : ' |
---|
| 207 | DO k=1,SIZE(fsst); WRITE(lunout,*) TRIM(fsst(k)); END DO |
---|
[2311] | 208 | CALL abort_physic('limit_netcdf','No sst file was found',1) |
---|
[1687] | 209 | END IF |
---|
[2293] | 210 | CALL ncerr(NF90_CLOSE(nid),sstfile) |
---|
[2947] | 211 | CALL msg(0,'Fichier choisi pour la temperature de mer: '//TRIM(sstfile)) |
---|
[1687] | 212 | |
---|
[2336] | 213 | CALL get_2Dfield(sstfile,varname,'SST',ndays,phy_sst,flag=extrap) |
---|
[1687] | 214 | |
---|
| 215 | !--- ALBEDO TREATMENT ---------------------------------------------------------- |
---|
[2947] | 216 | CALL msg(0,""); CALL msg(0," *** TRAITEMENT DE L'ALBEDO ***") |
---|
[2336] | 217 | CALL get_2Dfield(falbe,valb,'ALB',ndays,phy_alb) |
---|
[1687] | 218 | |
---|
| 219 | !--- REFERENCE GROUND HEAT FLUX TREATMENT -------------------------------------- |
---|
| 220 | ALLOCATE(phy_bil(klon,ndays)); phy_bil=0.0 |
---|
| 221 | |
---|
| 222 | !--- OUTPUT FILE WRITING ------------------------------------------------------- |
---|
[2947] | 223 | CALL msg(0,""); CALL msg(0,' *** Ecriture du fichier limit : debut ***') |
---|
[2293] | 224 | fnam="limit.nc" |
---|
[1687] | 225 | |
---|
| 226 | !--- File creation |
---|
[3803] | 227 | CALL ncerr(NF90_CREATE(fnam,IOR(NF90_CLOBBER,NF90_64BIT_OFFSET),nid),fnam) |
---|
[2293] | 228 | CALL ncerr(NF90_PUT_ATT(nid,NF90_GLOBAL,"title","Fichier conditions aux limites"),fnam) |
---|
[2889] | 229 | str='File produced using ce0l executable.' |
---|
| 230 | str=TRIM(str)//NEW_LINE(' ')//'Sea Ice Concentration built from' |
---|
| 231 | SELECT CASE(ix_sic) |
---|
| 232 | CASE(1); str=TRIM(str)//' Amip mid-month boundary condition (BCS).' |
---|
| 233 | CASE(2); str=TRIM(str)//' Amip monthly mean observations.' |
---|
| 234 | CASE(3); str=TRIM(str)//' IPSL coupled model outputs.' |
---|
| 235 | CASE(4); str=TRIM(str)//' LMDZ model outputs.' |
---|
| 236 | CASE(5); str=TRIM(str)//' ci.nc file.' |
---|
| 237 | END SELECT |
---|
| 238 | str=TRIM(str)//NEW_LINE(' ')//'Sea Surface Temperature built from' |
---|
| 239 | SELECT CASE(ix_sst) |
---|
| 240 | CASE(1); str=TRIM(str)//' Amip mid-month boundary condition (BCS).' |
---|
| 241 | CASE(2); str=TRIM(str)//' Amip monthly mean observations.' |
---|
| 242 | CASE(3); str=TRIM(str)//' IPSL coupled model outputs.' |
---|
| 243 | CASE(4); str=TRIM(str)//' LMDZ model outputs.' |
---|
| 244 | CASE(5); str=TRIM(str)//' sstk.nc file.' |
---|
| 245 | END SELECT |
---|
| 246 | CALL ncerr(NF90_PUT_ATT(nid,NF90_GLOBAL,"history",TRIM(str)),fnam) |
---|
[1687] | 247 | |
---|
| 248 | !--- Dimensions creation |
---|
[2293] | 249 | CALL ncerr(NF90_DEF_DIM(nid,"points_physiques",klon,ndim),fnam) |
---|
| 250 | CALL ncerr(NF90_DEF_DIM(nid,"time",NF90_UNLIMITED,ntim),fnam) |
---|
[1687] | 251 | |
---|
[2293] | 252 | dims=[ndim,ntim] |
---|
[1687] | 253 | |
---|
| 254 | !--- Variables creation |
---|
[2293] | 255 | CALL ncerr(NF90_DEF_VAR(nid,"TEMPS",NF90_FORMAT,[ntim],id_tim),fnam) |
---|
| 256 | CALL ncerr(NF90_DEF_VAR(nid,"FOCE", NF90_FORMAT,dims,id_FOCE),fnam) |
---|
| 257 | CALL ncerr(NF90_DEF_VAR(nid,"FSIC", NF90_FORMAT,dims,id_FSIC),fnam) |
---|
| 258 | CALL ncerr(NF90_DEF_VAR(nid,"FTER", NF90_FORMAT,dims,id_FTER),fnam) |
---|
| 259 | CALL ncerr(NF90_DEF_VAR(nid,"FLIC", NF90_FORMAT,dims,id_FLIC),fnam) |
---|
| 260 | CALL ncerr(NF90_DEF_VAR(nid,"SST", NF90_FORMAT,dims,id_SST),fnam) |
---|
| 261 | CALL ncerr(NF90_DEF_VAR(nid,"BILS", NF90_FORMAT,dims,id_BILS),fnam) |
---|
| 262 | CALL ncerr(NF90_DEF_VAR(nid,"ALB", NF90_FORMAT,dims,id_ALB),fnam) |
---|
| 263 | CALL ncerr(NF90_DEF_VAR(nid,"RUG", NF90_FORMAT,dims,id_RUG),fnam) |
---|
[2128] | 264 | call nf95_def_var(nid, "longitude", NF90_FLOAT, ndim, varid_longitude) |
---|
[2293] | 265 | call nf95_def_var(nid, "latitude", NF90_FLOAT, ndim, varid_latitude) |
---|
[1687] | 266 | |
---|
| 267 | !--- Attributes creation |
---|
[2293] | 268 | CALL ncerr(NF90_PUT_ATT(nid,id_tim, "title","Jour dans l annee"),fnam) |
---|
[2665] | 269 | CALL ncerr(NF90_PUT_ATT(nid,id_tim, "calendar",calend),fnam) |
---|
[2293] | 270 | CALL ncerr(NF90_PUT_ATT(nid,id_FOCE,"title","Fraction ocean"),fnam) |
---|
| 271 | CALL ncerr(NF90_PUT_ATT(nid,id_FSIC,"title","Fraction glace de mer"),fnam) |
---|
| 272 | CALL ncerr(NF90_PUT_ATT(nid,id_FTER,"title","Fraction terre"),fnam) |
---|
| 273 | CALL ncerr(NF90_PUT_ATT(nid,id_FLIC,"title","Fraction land ice"),fnam) |
---|
| 274 | CALL ncerr(NF90_PUT_ATT(nid,id_SST ,"title","Temperature superficielle de la mer"),fnam) |
---|
| 275 | CALL ncerr(NF90_PUT_ATT(nid,id_BILS,"title","Reference flux de chaleur au sol"),fnam) |
---|
| 276 | CALL ncerr(NF90_PUT_ATT(nid,id_ALB, "title","Albedo a la surface"),fnam) |
---|
| 277 | CALL ncerr(NF90_PUT_ATT(nid,id_RUG, "title","Rugosite"),fnam) |
---|
[1687] | 278 | |
---|
[2128] | 279 | call nf95_put_att(nid, varid_longitude, "standard_name", "longitude") |
---|
| 280 | call nf95_put_att(nid, varid_longitude, "units", "degrees_east") |
---|
| 281 | |
---|
| 282 | call nf95_put_att(nid, varid_latitude, "standard_name", "latitude") |
---|
| 283 | call nf95_put_att(nid, varid_latitude, "units", "degrees_north") |
---|
| 284 | |
---|
[2293] | 285 | CALL ncerr(NF90_ENDDEF(nid),fnam) |
---|
[1687] | 286 | |
---|
| 287 | !--- Variables saving |
---|
[2293] | 288 | CALL ncerr(NF90_PUT_VAR(nid,id_tim,[(REAL(k),k=1,ndays)]),fnam) |
---|
| 289 | CALL ncerr(NF90_PUT_VAR(nid,id_FOCE,pctsrf_t(:,is_oce,:),[1,1],[klon,ndays]),fnam) |
---|
| 290 | CALL ncerr(NF90_PUT_VAR(nid,id_FSIC,pctsrf_t(:,is_sic,:),[1,1],[klon,ndays]),fnam) |
---|
| 291 | CALL ncerr(NF90_PUT_VAR(nid,id_FTER,pctsrf_t(:,is_ter,:),[1,1],[klon,ndays]),fnam) |
---|
| 292 | CALL ncerr(NF90_PUT_VAR(nid,id_FLIC,pctsrf_t(:,is_lic,:),[1,1],[klon,ndays]),fnam) |
---|
| 293 | CALL ncerr(NF90_PUT_VAR(nid,id_SST ,phy_sst(:,:),[1,1],[klon,ndays]),fnam) |
---|
| 294 | CALL ncerr(NF90_PUT_VAR(nid,id_BILS,phy_bil(:,:),[1,1],[klon,ndays]),fnam) |
---|
| 295 | CALL ncerr(NF90_PUT_VAR(nid,id_ALB ,phy_alb(:,:),[1,1],[klon,ndays]),fnam) |
---|
| 296 | CALL ncerr(NF90_PUT_VAR(nid,id_RUG ,phy_rug(:,:),[1,1],[klon,ndays]),fnam) |
---|
[2399] | 297 | call nf95_put_var(nid, varid_longitude, longitude_deg) |
---|
| 298 | call nf95_put_var(nid, varid_latitude, latitude_deg) |
---|
[1687] | 299 | |
---|
[2293] | 300 | CALL ncerr(NF90_CLOSE(nid),fnam) |
---|
[1687] | 301 | |
---|
[2947] | 302 | CALL msg(0,""); CALL msg(0,' *** Ecriture du fichier limit : fin ***') |
---|
[1687] | 303 | |
---|
| 304 | DEALLOCATE(pctsrf_t,phy_sst,phy_bil,phy_alb,phy_rug) |
---|
| 305 | |
---|
| 306 | |
---|
| 307 | !=============================================================================== |
---|
| 308 | ! |
---|
| 309 | CONTAINS |
---|
| 310 | ! |
---|
| 311 | !=============================================================================== |
---|
| 312 | |
---|
| 313 | |
---|
| 314 | !------------------------------------------------------------------------------- |
---|
| 315 | ! |
---|
[2336] | 316 | SUBROUTINE get_2Dfield(fnam, varname, mode, ndays, champo, flag, mask) |
---|
[1687] | 317 | ! |
---|
| 318 | !----------------------------------------------------------------------------- |
---|
| 319 | ! Comments: |
---|
| 320 | ! There are two assumptions concerning the NetCDF files, that are satisfied |
---|
| 321 | ! with files that are conforming NC convention: |
---|
| 322 | ! 1) The last dimension of the variables used is the time record. |
---|
| 323 | ! 2) Dimensional variables have the same names as corresponding dimensions. |
---|
| 324 | !----------------------------------------------------------------------------- |
---|
[5084] | 325 | USE netcdf, ONLY: NF90_OPEN, NF90_INQ_VARID, NF90_INQUIRE_VARIABLE, & |
---|
[1687] | 326 | NF90_CLOSE, NF90_INQ_DIMID, NF90_INQUIRE_DIMENSION, NF90_GET_VAR, & |
---|
| 327 | NF90_GET_ATT |
---|
[1785] | 328 | USE pchsp_95_m, only: pchsp_95 |
---|
| 329 | USE pchfe_95_m, only: pchfe_95 |
---|
| 330 | USE arth_m, only: arth |
---|
| 331 | USE indice_sol_mod |
---|
[1687] | 332 | |
---|
[5271] | 333 | USE dimensions_mod, ONLY: iim, jjm, llm, ndm |
---|
[5285] | 334 | USE paramet_mod_h |
---|
[5271] | 335 | IMPLICIT NONE |
---|
| 336 | |
---|
[5272] | 337 | |
---|
[1687] | 338 | !----------------------------------------------------------------------------- |
---|
| 339 | ! Arguments: |
---|
| 340 | CHARACTER(LEN=*), INTENT(IN) :: fnam ! NetCDF file name |
---|
[3168] | 341 | CHARACTER(LEN=*), INTENT(IN) :: varname ! NetCDF variable name |
---|
| 342 | CHARACTER(LEN=*), INTENT(IN) :: mode ! RUG, SIC, SST or ALB |
---|
[1687] | 343 | INTEGER, INTENT(IN) :: ndays ! current year number of days |
---|
| 344 | REAL, POINTER, DIMENSION(:, :) :: champo ! output field = f(t) |
---|
| 345 | LOGICAL, OPTIONAL, INTENT(IN) :: flag ! extrapol. (SST) old ice (SIC) |
---|
| 346 | REAL, OPTIONAL, DIMENSION(iim, jjp1), INTENT(IN) :: mask |
---|
| 347 | !------------------------------------------------------------------------------ |
---|
| 348 | ! Local variables: |
---|
| 349 | !--- NetCDF |
---|
[2293] | 350 | INTEGER :: ncid, varid ! NetCDF identifiers |
---|
[3168] | 351 | CHARACTER(LEN=ns) :: dnam ! dimension name |
---|
[1687] | 352 | !--- dimensions |
---|
[2293] | 353 | INTEGER :: dids(4) ! NetCDF dimensions identifiers |
---|
| 354 | REAL, ALLOCATABLE :: dlon_ini(:) ! initial longitudes vector |
---|
| 355 | REAL, ALLOCATABLE :: dlat_ini(:) ! initial latitudes vector |
---|
| 356 | REAL, POINTER :: dlon(:), dlat(:) ! reordered lon/lat vectors |
---|
[1687] | 357 | !--- fields |
---|
[2293] | 358 | INTEGER :: imdep, jmdep, lmdep ! dimensions of 'champ' |
---|
| 359 | REAL, ALLOCATABLE :: champ(:,:) ! wanted field on initial grid |
---|
| 360 | REAL, ALLOCATABLE :: yder(:), timeyear(:) |
---|
| 361 | REAL :: champint(iim,jjp1) ! interpolated field |
---|
| 362 | REAL, ALLOCATABLE :: champtime(:,:,:) |
---|
| 363 | REAL, ALLOCATABLE :: champan(:,:,:) |
---|
[1687] | 364 | !--- input files |
---|
[3168] | 365 | CHARACTER(LEN=ns) :: fnam_m, fnam_p ! previous/next files names |
---|
| 366 | CHARACTER(LEN=ns) :: cal_in ! calendar |
---|
| 367 | CHARACTER(LEN=ns) :: units ! attribute "units" in sic/sst file |
---|
[2293] | 368 | INTEGER :: ndays_in ! number of days |
---|
[3163] | 369 | REAL :: value ! mean/max value near equator |
---|
[1687] | 370 | !--- misc |
---|
[3168] | 371 | INTEGER :: i, j, k, l ! loop counters |
---|
[2293] | 372 | REAL, ALLOCATABLE :: work(:,:) ! used for extrapolation |
---|
[3168] | 373 | CHARACTER(LEN=ns) :: title, mess ! for messages |
---|
[2889] | 374 | LOGICAL :: is_bcs ! flag for BCS data |
---|
[2293] | 375 | LOGICAL :: extrp ! flag for extrapolation |
---|
[3168] | 376 | LOGICAL :: ll |
---|
[2889] | 377 | REAL :: chmin, chmax, timeday, al |
---|
[2576] | 378 | INTEGER ierr, idx |
---|
[1687] | 379 | integer n_extrap ! number of extrapolated points |
---|
| 380 | logical skip |
---|
| 381 | |
---|
| 382 | !------------------------------------------------------------------------------ |
---|
| 383 | !---Variables depending on keyword 'mode' ------------------------------------- |
---|
| 384 | NULLIFY(champo) |
---|
| 385 | |
---|
| 386 | SELECT CASE(mode) |
---|
| 387 | CASE('RUG'); title='Rugosite' |
---|
| 388 | CASE('SIC'); title='Sea-ice' |
---|
| 389 | CASE('SST'); title='SST' |
---|
| 390 | CASE('ALB'); title='Albedo' |
---|
| 391 | END SELECT |
---|
[2336] | 392 | extrp=.FALSE.; IF(PRESENT(flag).AND.mode=='SST') extrp=flag |
---|
[2889] | 393 | is_bcs=(mode=='SIC'.AND.ix_sic==1).OR.(mode=='SST'.AND.ix_sst==1) |
---|
[2576] | 394 | idx=INDEX(fnam,'.nc')-1 |
---|
[1687] | 395 | |
---|
| 396 | !--- GETTING SOME DIMENSIONAL VARIABLES FROM FILE ----------------------------- |
---|
[2293] | 397 | CALL msg(5,' Now reading file : '//TRIM(fnam)) |
---|
| 398 | CALL ncerr(NF90_OPEN(fnam, NF90_NOWRITE, ncid),fnam) |
---|
| 399 | CALL ncerr(NF90_INQ_VARID(ncid, trim(varname), varid),fnam) |
---|
| 400 | CALL ncerr(NF90_INQUIRE_VARIABLE(ncid, varid, dimids=dids),fnam) |
---|
[1687] | 401 | |
---|
| 402 | !--- Longitude |
---|
[2293] | 403 | CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(1), name=dnam, len=imdep),fnam) |
---|
| 404 | ALLOCATE(dlon_ini(imdep), dlon(imdep)) |
---|
| 405 | CALL ncerr(NF90_INQ_VARID(ncid, dnam, varid), fnam) |
---|
| 406 | CALL ncerr(NF90_GET_VAR(ncid, varid, dlon_ini), fnam) |
---|
| 407 | CALL msg(5,'variable '//TRIM(dnam)//' dimension ', imdep) |
---|
[1687] | 408 | |
---|
| 409 | !--- Latitude |
---|
[2293] | 410 | CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(2), name=dnam, len=jmdep),fnam) |
---|
| 411 | ALLOCATE(dlat_ini(jmdep), dlat(jmdep)) |
---|
| 412 | CALL ncerr(NF90_INQ_VARID(ncid, dnam, varid), fnam) |
---|
| 413 | CALL ncerr(NF90_GET_VAR(ncid, varid, dlat_ini), fnam) |
---|
| 414 | CALL msg(5,'variable '//TRIM(dnam)//' dimension ', jmdep) |
---|
[1687] | 415 | |
---|
| 416 | !--- Time (variable is not needed - it is rebuilt - but calendar is) |
---|
[2293] | 417 | CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(3), name=dnam, len=lmdep), fnam) |
---|
[2576] | 418 | ALLOCATE(timeyear(lmdep+2)) |
---|
[2293] | 419 | CALL ncerr(NF90_INQ_VARID(ncid, dnam, varid), fnam) |
---|
[1687] | 420 | cal_in=' ' |
---|
[2293] | 421 | IF(NF90_GET_ATT(ncid, varid, 'calendar', cal_in)/=NF90_NOERR) THEN |
---|
[1687] | 422 | SELECT CASE(mode) |
---|
[4361] | 423 | CASE('RUG', 'ALB'); cal_in='360_day' |
---|
[1687] | 424 | CASE('SIC', 'SST'); cal_in='gregorian' |
---|
| 425 | END SELECT |
---|
[2947] | 426 | CALL msg(0,'WARNING: missing "calendar" attribute for "time" in '& |
---|
[2293] | 427 | &//TRIM(fnam)//'. Choosing default value.') |
---|
[1687] | 428 | END IF |
---|
[2665] | 429 | CALL strclean(cal_in) !--- REMOVE (WEIRD) NULL CHARACTERS |
---|
[2947] | 430 | CALL msg(0,'var, calendar, dim: '//TRIM(dnam)//' '//TRIM(cal_in), lmdep) |
---|
[1687] | 431 | |
---|
[2540] | 432 | !--- Determining input file number of days, depending on calendar |
---|
[1687] | 433 | ndays_in=year_len(anneeref, cal_in) |
---|
| 434 | |
---|
[2540] | 435 | !--- Rebuilding input time vector (field from input file might be unreliable) |
---|
| 436 | IF(lmdep==12) THEN |
---|
[2766] | 437 | timeyear=mid_month(anneeref, cal_in) |
---|
[2576] | 438 | CALL msg(0,'Monthly input file(s) for '//TRIM(title)//'.') |
---|
[2540] | 439 | ELSE IF(lmdep==ndays_in) THEN |
---|
[2576] | 440 | timeyear=[(REAL(k)-0.5,k=0,ndays_in+1)] |
---|
[2540] | 441 | CALL msg(0,'Daily input file (no time interpolation).') |
---|
| 442 | ELSE |
---|
| 443 | WRITE(mess,'(a,i3,a,i3,a)')'Mismatching input file: found',lmdep, & |
---|
[2576] | 444 | ' records, 12/',ndays_in,' (monthly/daily needed).' |
---|
[2766] | 445 | CALL abort_physic('mid_month',TRIM(mess),1) |
---|
[2540] | 446 | END IF |
---|
[1687] | 447 | |
---|
| 448 | !--- GETTING THE FIELD AND INTERPOLATING IT ---------------------------------- |
---|
[2576] | 449 | ALLOCATE(champ(imdep, jmdep), champtime(iim, jjp1, lmdep+2)) |
---|
[1687] | 450 | IF(extrp) ALLOCATE(work(imdep, jmdep)) |
---|
[2293] | 451 | CALL msg(5,'') |
---|
[2766] | 452 | CALL msg(5,'READ AND INTERPOLATE HORIZONTALLY ', lmdep, ' FIELDS.') |
---|
[2293] | 453 | CALL ncerr(NF90_INQ_VARID(ncid, varname, varid), fnam) |
---|
[1687] | 454 | DO l=1, lmdep |
---|
[2293] | 455 | CALL ncerr(NF90_GET_VAR(ncid,varid,champ,[1,1,l],[imdep,jmdep,1]),fnam) |
---|
[3163] | 456 | CALL conf_dat2d(title, dlon_ini, dlat_ini, dlon, dlat, champ, .TRUE.) |
---|
| 457 | |
---|
| 458 | !--- FOR SIC/SST FIELDS ONLY |
---|
| 459 | IF(l==1.AND.is_in(mode,['SIC','SST'])) THEN |
---|
| 460 | |
---|
| 461 | !--- DETERMINE THE UNIT: READ FROM FILE OR ASSUMED USING FIELD VALUES |
---|
| 462 | ierr=NF90_GET_ATT(ncid, varid, 'units', units) |
---|
| 463 | IF(ierr==NF90_NOERR) THEN !--- ATTRIBUTE "units" FOUND IN THE FILE |
---|
| 464 | CALL strclean(units) |
---|
| 465 | IF(mode=='SIC'.AND.is_in(units,Perc)) units="%" |
---|
| 466 | IF(mode=='SIC'.AND.is_in(units,Frac)) units="1" |
---|
| 467 | IF(mode=='SST'.AND.is_in(units,DegC)) units="C" |
---|
| 468 | IF(mode=='SST'.AND.is_in(units,DegK)) units="K" |
---|
| 469 | ELSE !--- CHECK THE FIELD VALUES |
---|
| 470 | IF(mode=='SIC') value=MAXVAL(champ(:,:)) |
---|
| 471 | IF(mode=='SST') value= SUM(champ(:,jmdep/2),DIM=1)/REAL(imdep) |
---|
| 472 | IF(mode=='SIC') THEN; units="1"; IF(value>= 10.) units="%"; END IF |
---|
| 473 | IF(mode=='SST') THEN; units="C"; IF(value>=100.) units="K"; END IF |
---|
[2766] | 474 | END IF |
---|
[3163] | 475 | CALL msg(0,'INPUT FILE '//TRIM(title)//' UNIT IS: "'//TRIM(units)//'".') |
---|
| 476 | IF(ierr/=NF90_NOERR) CALL msg(0,'WARNING ! UNIT TO BE CHECKED ! ' & |
---|
| 477 | //'No "units" attribute, so only based on the fields values.') |
---|
| 478 | |
---|
| 479 | !--- CHECK VALUES ARE IN THE EXPECTED RANGE |
---|
| 480 | SELECT CASE(units) |
---|
| 481 | CASE('%'); ll=ANY(champ>100.0+EPSFRA); str='percentages > 100.' |
---|
| 482 | CASE('1'); ll=ANY(champ> 1.0+EPSFRA); str='fractions > 1.' |
---|
| 483 | CASE('C'); ll=ANY(champ<-100.).OR.ANY(champ> 60.); str='<-100 or >60 DegC' |
---|
| 484 | CASE('K'); ll=ANY(champ< 180.).OR.ANY(champ>330.); str='<180 or >330 DegK' |
---|
| 485 | CASE DEFAULT; CALL abort_physic(mode, 'Unrecognized '//TRIM(title) & |
---|
| 486 | //' unit: '//TRIM(units),1) |
---|
| 487 | END SELECT |
---|
| 488 | |
---|
| 489 | !--- DROPPED FOR BCS DATA (FRACTIONS CAN BE HIGHER THAN 1) |
---|
| 490 | IF(ll.AND.ix_sic/=1.AND.mode=='SIC') & |
---|
[4257] | 491 | CALL abort_physic(mode,'unrealistic '//TRIM(mode)//' found: '//TRIM(str), 1) |
---|
[3163] | 492 | |
---|
[2766] | 493 | END IF |
---|
[3163] | 494 | |
---|
[2293] | 495 | IF(extrp) CALL extrapol(champ,imdep,jmdep,999999.,.TRUE.,.TRUE.,2,work) |
---|
[2336] | 496 | IF(l==1) THEN |
---|
[3163] | 497 | CALL msg(5,"--------------------------------------------------------") |
---|
| 498 | CALL msg(5,"$$$ Barycentric interpolation for "//TRIM(title)//" $$$") |
---|
| 499 | CALL msg(5,"--------------------------------------------------------") |
---|
[1687] | 500 | END IF |
---|
[2336] | 501 | IF(mode=='RUG') champ=LOG(champ) |
---|
| 502 | CALL inter_barxy(dlon,dlat(:jmdep-1),champ,rlonu(:iim),rlatv,champint) |
---|
| 503 | IF(mode=='RUG') THEN |
---|
| 504 | champint=EXP(champint) |
---|
| 505 | WHERE(NINT(mask)/=1) champint=0.001 |
---|
| 506 | END IF |
---|
[2576] | 507 | champtime(:, :, l+1)=champint |
---|
[1687] | 508 | END DO |
---|
[2293] | 509 | CALL ncerr(NF90_CLOSE(ncid), fnam) |
---|
[1687] | 510 | |
---|
[2576] | 511 | !--- FIRST RECORD: LAST ONE OF PREVIOUS YEAR (CURRENT YEAR IF UNAVAILABLE) |
---|
| 512 | fnam_m=fnam(1:idx)//'_m.nc' |
---|
| 513 | IF(NF90_OPEN(fnam_m,NF90_NOWRITE,ncid)==NF90_NOERR) THEN |
---|
| 514 | CALL msg(0,'Reading previous year file ("'//TRIM(fnam_m)//'") last record for '//TRIM(title)) |
---|
| 515 | CALL ncerr(NF90_INQ_VARID(ncid, varname, varid),fnam_m) |
---|
| 516 | CALL ncerr(NF90_INQUIRE_VARIABLE(ncid, varid, dimids=dids),fnam_m) |
---|
| 517 | CALL ncerr(NF90_INQUIRE_DIMENSION(ncid, dids(3), len=l), fnam_m) |
---|
| 518 | CALL ncerr(NF90_GET_VAR(ncid,varid,champ,[1,1,l],[imdep,jmdep,1]),fnam_m) |
---|
| 519 | CALL ncerr(NF90_CLOSE(ncid), fnam_m) |
---|
| 520 | CALL conf_dat2d(title, dlon_ini, dlat_ini, dlon, dlat, champ, .TRUE.) |
---|
| 521 | IF(extrp) CALL extrapol(champ,imdep,jmdep,999999.,.TRUE.,.TRUE.,2,work) |
---|
| 522 | IF(mode=='RUG') champ=LOG(champ) |
---|
| 523 | CALL inter_barxy(dlon,dlat(:jmdep-1),champ,rlonu(:iim),rlatv,champint) |
---|
| 524 | IF(mode=='RUG') THEN |
---|
| 525 | champint=EXP(champint) |
---|
| 526 | WHERE(NINT(mask)/=1) champint=0.001 |
---|
| 527 | END IF |
---|
| 528 | champtime(:, :, 1)=champint |
---|
| 529 | ELSE |
---|
| 530 | CALL msg(0,'Using current year file ("'//TRIM(fnam)//'") last record for '//TRIM(title)) |
---|
| 531 | champtime(:, :, 1)=champtime(:, :, lmdep+1) |
---|
| 532 | END IF |
---|
| 533 | |
---|
| 534 | !--- LAST RECORD: FIRST ONE OF NEXT YEAR (CURRENT YEAR IF UNAVAILABLE) |
---|
| 535 | fnam_p=fnam(1:idx)//'_p.nc' |
---|
| 536 | IF(NF90_OPEN(fnam_p,NF90_NOWRITE,ncid)==NF90_NOERR) THEN |
---|
[2665] | 537 | CALL msg(0,'Reading next year file ("'//TRIM(fnam_p)//'") first record for '//TRIM(title)) |
---|
[2576] | 538 | CALL ncerr(NF90_INQ_VARID(ncid, varname, varid),fnam_p) |
---|
| 539 | CALL ncerr(NF90_GET_VAR(ncid,varid,champ,[1,1,1],[imdep,jmdep,1]),fnam_p) |
---|
| 540 | CALL ncerr(NF90_CLOSE(ncid), fnam_p) |
---|
| 541 | CALL conf_dat2d(title, dlon_ini, dlat_ini, dlon, dlat, champ, .TRUE.) |
---|
| 542 | IF(extrp) CALL extrapol(champ,imdep,jmdep,999999.,.TRUE.,.TRUE.,2,work) |
---|
| 543 | IF(mode=='RUG') champ=LOG(champ) |
---|
| 544 | CALL inter_barxy(dlon,dlat(:jmdep-1),champ,rlonu(:iim),rlatv,champint) |
---|
| 545 | IF(mode=='RUG') THEN |
---|
| 546 | champint=EXP(champint) |
---|
| 547 | WHERE(NINT(mask)/=1) champint=0.001 |
---|
| 548 | END IF |
---|
| 549 | champtime(:, :, lmdep+2)=champint |
---|
| 550 | ELSE |
---|
| 551 | CALL msg(0,'Using current year file ("'//TRIM(fnam)//'") first record for '//TRIM(title)) |
---|
| 552 | champtime(:, :, lmdep+2)=champtime(:, :, 2) |
---|
| 553 | END IF |
---|
[1687] | 554 | DEALLOCATE(dlon_ini, dlat_ini, dlon, dlat, champ) |
---|
| 555 | IF(extrp) DEALLOCATE(work) |
---|
| 556 | |
---|
| 557 | !--- TIME INTERPOLATION ------------------------------------------------------ |
---|
[2540] | 558 | IF(prt_level>0) THEN |
---|
| 559 | IF(ndays/=ndays_in) THEN |
---|
[2766] | 560 | WRITE(lunout,*)'DIFFERENT YEAR LENGTHS:' |
---|
[2947] | 561 | WRITE(lunout,*)' In the input file: ',ndays_in |
---|
[2766] | 562 | WRITE(lunout,*)' In the output file: ',ndays |
---|
[2540] | 563 | END IF |
---|
| 564 | IF(lmdep==ndays_in) THEN |
---|
[2766] | 565 | WRITE(lunout, *)'NO TIME INTERPOLATION.' |
---|
| 566 | WRITE(lunout, *)' Daily input file.' |
---|
[2540] | 567 | ELSE |
---|
[2889] | 568 | IF( is_bcs) WRITE(lunout, *)'LINEAR TIME INTERPOLATION.' |
---|
| 569 | IF(.NOT.is_bcs) WRITE(lunout, *)'SPLINES TIME INTERPOLATION.' |
---|
[2766] | 570 | WRITE(lunout, *)' Input time vector: ', timeyear |
---|
[3163] | 571 | WRITE(lunout, *)' Output time vector: from 0.5 to ', ndays-0.5 |
---|
[2540] | 572 | END IF |
---|
[1687] | 573 | END IF |
---|
[2889] | 574 | ALLOCATE(champan(iip1, jjp1, ndays)) |
---|
| 575 | |
---|
| 576 | IF(lmdep==ndays_in) THEN !--- DAILY DATA: NO TIME INTERPOLATION |
---|
[3380] | 577 | DO l=1,lmdep |
---|
| 578 | champan(1:iim,:,l)=champtime(:,:,l+1) |
---|
| 579 | END DO |
---|
[2889] | 580 | ELSE IF(is_bcs) THEN !--- BCS DATA: LINEAR TIME INTERPOLATION |
---|
| 581 | l=1 |
---|
| 582 | DO k=1, ndays |
---|
[2893] | 583 | timeday = (REAL(k)-0.5)*REAL(ndays_in)/ndays |
---|
[2889] | 584 | IF(timeyear(l+1)<timeday) l=l+1 |
---|
| 585 | al=(timeday-timeyear(l))/(timeyear(l+1)-timeyear(l)) |
---|
[3380] | 586 | champan(1:iim,:,k) = champtime(1:iim,:,l)+al*(champtime(1:iim,:,l+1)-champtime(1:iim,:,l)) |
---|
[2889] | 587 | END DO |
---|
| 588 | ELSE !--- AVE DATA: SPLINE TIME INTERPOLATION |
---|
[2540] | 589 | skip = .false. |
---|
| 590 | n_extrap = 0 |
---|
[2889] | 591 | ALLOCATE(yder(lmdep+2)) |
---|
[2540] | 592 | DO j=1, jjp1 |
---|
| 593 | DO i=1, iim |
---|
| 594 | yder = pchsp_95(timeyear, champtime(i, j, :), ibeg=2, iend=2, & |
---|
| 595 | vc_beg=0., vc_end=0.) |
---|
| 596 | CALL pchfe_95(timeyear, champtime(i, j, :), yder, skip, & |
---|
[2893] | 597 | arth(0.5, real(ndays_in) / ndays, ndays), champan(i, j, :), ierr) |
---|
[4467] | 598 | if (ierr < 0) call abort_physic("get_2Dfield", "", 1) |
---|
[2540] | 599 | n_extrap = n_extrap + ierr |
---|
| 600 | END DO |
---|
| 601 | END DO |
---|
| 602 | IF(n_extrap /= 0) WRITE(lunout,*) "get_2Dfield pchfe_95: n_extrap = ", n_extrap |
---|
[2889] | 603 | DEALLOCATE(yder) |
---|
[2540] | 604 | END IF |
---|
[1687] | 605 | champan(iip1, :, :)=champan(1, :, :) |
---|
[2889] | 606 | DEALLOCATE(champtime, timeyear) |
---|
[1687] | 607 | |
---|
| 608 | !--- Checking the result |
---|
| 609 | DO j=1, jjp1 |
---|
| 610 | CALL minmax(iip1, champan(1, j, 10), chmin, chmax) |
---|
[2766] | 611 | IF (prt_level>5) WRITE(lunout, *)' ',TRIM(title),' at time 10 ', chmin, chmax, j |
---|
[1687] | 612 | END DO |
---|
| 613 | |
---|
| 614 | !--- SPECIAL FILTER FOR SST: SST>271.38 -------------------------------------- |
---|
| 615 | IF(mode=='SST') THEN |
---|
[3163] | 616 | SELECT CASE(units) |
---|
| 617 | CASE("K"); CALL msg(0,'SST field is already in kelvins.') |
---|
| 618 | CASE("C"); CALL msg(0,'SST field converted from celcius degrees to kelvins.') |
---|
| 619 | champan(:, :, :)=champan(:, :, :)+273.15 |
---|
| 620 | END SELECT |
---|
| 621 | CALL msg(0,'Filtering SST: Sea Surface Temperature >= 271.38') |
---|
[1687] | 622 | WHERE(champan<271.38) champan=271.38 |
---|
| 623 | END IF |
---|
| 624 | |
---|
| 625 | !--- SPECIAL FILTER FOR SIC: 0.0<SIC<1.0 ------------------------------------- |
---|
| 626 | IF(mode=='SIC') THEN |
---|
[3163] | 627 | SELECT CASE(units) |
---|
| 628 | CASE("1"); CALL msg(0,'SIC field already in fraction of 1') |
---|
| 629 | CASE("%"); CALL msg(0,'SIC field converted from percentage to fraction of 1.') |
---|
[1687] | 630 | champan(:, :, :)=champan(:, :, :)/100. |
---|
[3163] | 631 | END SELECT |
---|
| 632 | CALL msg(0,'Filtering SIC: 0.0 <= Sea-ice <=1.0') |
---|
[1687] | 633 | WHERE(champan>1.0) champan=1.0 |
---|
| 634 | WHERE(champan<0.0) champan=0.0 |
---|
| 635 | END IF |
---|
| 636 | |
---|
| 637 | !--- DYNAMICAL TO PHYSICAL GRID ---------------------------------------------- |
---|
| 638 | ALLOCATE(champo(klon, ndays)) |
---|
| 639 | DO k=1, ndays |
---|
[4261] | 640 | CALL gr_dyn_fi(1, iip1, jjp1, klon, champan(1, 1, k), champo(:, k)) |
---|
[1687] | 641 | END DO |
---|
| 642 | DEALLOCATE(champan) |
---|
| 643 | |
---|
| 644 | END SUBROUTINE get_2Dfield |
---|
| 645 | ! |
---|
| 646 | !------------------------------------------------------------------------------- |
---|
| 647 | |
---|
| 648 | |
---|
| 649 | !------------------------------------------------------------------------------- |
---|
| 650 | ! |
---|
[2336] | 651 | SUBROUTINE start_init_orog0(lon_in,lat_in,phis,masque) |
---|
| 652 | ! |
---|
| 653 | !------------------------------------------------------------------------------- |
---|
[2576] | 654 | USE grid_noro_m, ONLY: grid_noro0 |
---|
[2336] | 655 | IMPLICIT NONE |
---|
| 656 | !=============================================================================== |
---|
| 657 | ! Purpose: Compute "phis" just like it would be in start_init_orog. |
---|
| 658 | !=============================================================================== |
---|
| 659 | ! Arguments: |
---|
| 660 | REAL, INTENT(IN) :: lon_in(:), lat_in(:) ! dim (iml) (jml) |
---|
| 661 | REAL, INTENT(INOUT) :: phis(:,:), masque(:,:) ! dim (iml,jml) |
---|
| 662 | !------------------------------------------------------------------------------- |
---|
| 663 | ! Local variables: |
---|
[3168] | 664 | CHARACTER(LEN=ns) :: modname="start_init_orog0" |
---|
[2336] | 665 | INTEGER :: fid, llm_tmp,ttm_tmp, iml,jml, iml_rel,jml_rel, itau(1) |
---|
| 666 | REAL :: lev(1), date, dt, deg2rad |
---|
| 667 | REAL, ALLOCATABLE :: lon_rad(:), lon_ini(:), lon_rel(:,:), relief_hi(:,:) |
---|
| 668 | REAL, ALLOCATABLE :: lat_rad(:), lat_ini(:), lat_rel(:,:) |
---|
| 669 | !------------------------------------------------------------------------------- |
---|
| 670 | iml=assert_eq(SIZE(lon_in),SIZE(phis,1),SIZE(masque,1),TRIM(modname)//" iml") |
---|
| 671 | jml=assert_eq(SIZE(lat_in),SIZE(phis,2),SIZE(masque,2),TRIM(modname)//" jml") |
---|
| 672 | IF(iml/=iip1) CALL abort_gcm(TRIM(modname),'iml/=iip1',1) |
---|
| 673 | IF(jml/=jjp1) CALL abort_gcm(TRIM(modname),'jml/=jjp1',1) |
---|
| 674 | pi=2.0*ASIN(1.0); deg2rad=pi/180.0 |
---|
| 675 | IF(ANY(phis/=-99999.)) RETURN !--- phis ALREADY KNOWN |
---|
| 676 | |
---|
| 677 | !--- HIGH RESOLUTION OROGRAPHY |
---|
| 678 | CALL flininfo(frelf, iml_rel, jml_rel, llm_tmp, ttm_tmp, fid) |
---|
| 679 | |
---|
| 680 | ALLOCATE(lat_rel(iml_rel,jml_rel),lon_rel(iml_rel,jml_rel)) |
---|
| 681 | CALL flinopen(frelf, .FALSE., iml_rel, jml_rel, llm_tmp, lon_rel, lat_rel, & |
---|
| 682 | lev, ttm_tmp, itau, date, dt, fid) |
---|
| 683 | ALLOCATE(relief_hi(iml_rel,jml_rel)) |
---|
| 684 | CALL flinget(fid, vrel, iml_rel, jml_rel, llm_tmp, ttm_tmp, 1, 1, relief_hi) |
---|
| 685 | CALL flinclo(fid) |
---|
| 686 | |
---|
| 687 | !--- IF ANGLES ARE IN DEGREES, THEY ARE CONVERTED INTO RADIANS |
---|
| 688 | ALLOCATE(lon_ini(iml_rel),lat_ini(jml_rel)) |
---|
| 689 | lon_ini(:)=lon_rel(:,1); IF(MAXVAL(lon_rel)>pi) lon_ini=lon_ini*deg2rad |
---|
| 690 | lat_ini(:)=lat_rel(1,:); IF(MAXVAL(lat_rel)>pi) lat_ini=lat_ini*deg2rad |
---|
| 691 | |
---|
| 692 | !--- FIELDS ARE PROCESSED TO BE ON STANDARD ANGULAR DOMAINS |
---|
| 693 | ALLOCATE(lon_rad(iml_rel),lat_rad(jml_rel)) |
---|
| 694 | CALL conf_dat2d(vrel, lon_ini, lat_ini, lon_rad, lat_rad, relief_hi, .FALSE.) |
---|
| 695 | DEALLOCATE(lon_ini,lat_ini) |
---|
| 696 | |
---|
| 697 | !--- COMPUTING SURFACE GEOPOTENTIAL USING ROUTINE grid_noro0 |
---|
| 698 | WRITE(lunout,*) |
---|
| 699 | WRITE(lunout,*)'*** Compute surface geopotential ***' |
---|
| 700 | |
---|
| 701 | !--- CALL OROGRAPHY MODULE (REDUCED VERSION) TO COMPUTE FIELDS |
---|
| 702 | CALL grid_noro0(lon_rad, lat_rad, relief_hi, lon_in, lat_in, phis, masque) |
---|
| 703 | phis = phis * 9.81 |
---|
| 704 | phis(iml,:) = phis(1,:) |
---|
| 705 | DEALLOCATE(relief_hi,lon_rad,lat_rad) |
---|
| 706 | |
---|
| 707 | END SUBROUTINE start_init_orog0 |
---|
| 708 | ! |
---|
| 709 | !------------------------------------------------------------------------------- |
---|
| 710 | |
---|
| 711 | |
---|
| 712 | !------------------------------------------------------------------------------- |
---|
| 713 | ! |
---|
[2293] | 714 | SUBROUTINE msg(lev,str1,i,str2) |
---|
| 715 | ! |
---|
| 716 | !------------------------------------------------------------------------------- |
---|
| 717 | ! Arguments: |
---|
| 718 | INTEGER, INTENT(IN) :: lev |
---|
| 719 | CHARACTER(LEN=*), INTENT(IN) :: str1 |
---|
| 720 | INTEGER, OPTIONAL, INTENT(IN) :: i |
---|
| 721 | CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: str2 |
---|
| 722 | !------------------------------------------------------------------------------- |
---|
[2947] | 723 | IF(prt_level>=lev) THEN |
---|
[2293] | 724 | IF(PRESENT(str2)) THEN |
---|
| 725 | WRITE(lunout,*) TRIM(str1), i, TRIM(str2) |
---|
| 726 | ELSE IF(PRESENT(i)) THEN |
---|
| 727 | WRITE(lunout,*) TRIM(str1), i |
---|
| 728 | ELSE |
---|
| 729 | WRITE(lunout,*) TRIM(str1) |
---|
| 730 | END IF |
---|
| 731 | END IF |
---|
| 732 | |
---|
| 733 | END SUBROUTINE msg |
---|
| 734 | ! |
---|
| 735 | !------------------------------------------------------------------------------- |
---|
| 736 | |
---|
| 737 | |
---|
| 738 | !------------------------------------------------------------------------------- |
---|
| 739 | ! |
---|
[1687] | 740 | SUBROUTINE ncerr(ncres,fnam) |
---|
| 741 | ! |
---|
| 742 | !------------------------------------------------------------------------------- |
---|
| 743 | ! Purpose: NetCDF errors handling. |
---|
| 744 | !------------------------------------------------------------------------------- |
---|
[5084] | 745 | USE netcdf, ONLY : NF90_NOERR, NF90_STRERROR |
---|
[1687] | 746 | IMPLICIT NONE |
---|
| 747 | !------------------------------------------------------------------------------- |
---|
| 748 | ! Arguments: |
---|
| 749 | INTEGER, INTENT(IN) :: ncres |
---|
| 750 | CHARACTER(LEN=*), INTENT(IN) :: fnam |
---|
| 751 | !------------------------------------------------------------------------------- |
---|
| 752 | IF(ncres/=NF90_NOERR) THEN |
---|
| 753 | WRITE(lunout,*)'Problem with file '//TRIM(fnam)//' in routine limit_netcdf.' |
---|
[2311] | 754 | CALL abort_physic('limit_netcdf',NF90_STRERROR(ncres),1) |
---|
[1687] | 755 | END IF |
---|
| 756 | |
---|
| 757 | END SUBROUTINE ncerr |
---|
| 758 | ! |
---|
| 759 | !------------------------------------------------------------------------------- |
---|
| 760 | |
---|
[2665] | 761 | |
---|
| 762 | !------------------------------------------------------------------------------- |
---|
| 763 | ! |
---|
| 764 | SUBROUTINE strclean(s) |
---|
| 765 | ! |
---|
| 766 | !------------------------------------------------------------------------------- |
---|
| 767 | IMPLICIT NONE |
---|
| 768 | !------------------------------------------------------------------------------- |
---|
| 769 | ! Purpose: Remove tail null characters from the input string. |
---|
| 770 | !------------------------------------------------------------------------------- |
---|
| 771 | ! Parameters: |
---|
| 772 | CHARACTER(LEN=*), INTENT(INOUT) :: s |
---|
| 773 | !------------------------------------------------------------------------------- |
---|
| 774 | ! Local variable: |
---|
| 775 | INTEGER :: k |
---|
| 776 | !------------------------------------------------------------------------------- |
---|
| 777 | k=LEN_TRIM(s); DO WHILE(ICHAR(s(k:k))==0); s(k:k)=' '; k=LEN_TRIM(s); END DO |
---|
| 778 | |
---|
| 779 | END SUBROUTINE strclean |
---|
| 780 | ! |
---|
| 781 | !------------------------------------------------------------------------------- |
---|
| 782 | |
---|
[3163] | 783 | |
---|
| 784 | !------------------------------------------------------------------------------- |
---|
| 785 | ! |
---|
| 786 | FUNCTION is_in(s1,s2) RESULT(res) |
---|
| 787 | ! |
---|
| 788 | !------------------------------------------------------------------------------- |
---|
| 789 | IMPLICIT NONE |
---|
| 790 | !------------------------------------------------------------------------------- |
---|
| 791 | ! Purpose: Check wether s1 is present in the s2(:) list (case insensitive). |
---|
| 792 | !------------------------------------------------------------------------------- |
---|
| 793 | ! Arguments: |
---|
| 794 | CHARACTER(LEN=*), INTENT(IN) :: s1, s2(:) |
---|
| 795 | LOGICAL :: res |
---|
| 796 | !------------------------------------------------------------------------------- |
---|
| 797 | res=.FALSE.; DO k=1,SIZE(s2); res=res.OR.strLow(s1)==strLow(s2(k)); END DO |
---|
| 798 | |
---|
| 799 | END FUNCTION is_in |
---|
| 800 | ! |
---|
| 801 | !------------------------------------------------------------------------------- |
---|
| 802 | |
---|
| 803 | |
---|
| 804 | !------------------------------------------------------------------------------- |
---|
| 805 | ! |
---|
| 806 | ELEMENTAL FUNCTION strLow(s) RESULT(res) |
---|
| 807 | ! |
---|
| 808 | !------------------------------------------------------------------------------- |
---|
| 809 | IMPLICIT NONE |
---|
| 810 | !------------------------------------------------------------------------------- |
---|
| 811 | ! Purpose: Lower case conversion. |
---|
| 812 | !------------------------------------------------------------------------------- |
---|
| 813 | ! Arguments: |
---|
| 814 | CHARACTER(LEN=*), INTENT(IN) :: s |
---|
| 815 | CHARACTER(LEN=ns) :: res |
---|
| 816 | !------------------------------------------------------------------------------- |
---|
| 817 | ! Local variable: |
---|
| 818 | INTEGER :: k, ix |
---|
| 819 | !------------------------------------------------------------------------------- |
---|
| 820 | res=s |
---|
| 821 | DO k=1,LEN(s); ix=IACHAR(s(k:k)) |
---|
| 822 | IF(64<ix.AND.ix<91) res(k:k)=ACHAR(ix+32) |
---|
| 823 | END DO |
---|
| 824 | |
---|
| 825 | END FUNCTION strLow |
---|
| 826 | ! |
---|
| 827 | !------------------------------------------------------------------------------- |
---|
| 828 | |
---|
[1687] | 829 | #endif |
---|
[2239] | 830 | ! of #ifndef CPP_1D |
---|
[1687] | 831 | END SUBROUTINE limit_netcdf |
---|
[2336] | 832 | |
---|
| 833 | END MODULE limit |
---|
| 834 | ! |
---|
| 835 | !******************************************************************************* |
---|
[2540] | 836 | |
---|