! ! $Id: limit_netcdf.F90 1508 2011-04-15 13:05:34Z jghattas $ !------------------------------------------------------------------------------- ! SUBROUTINE limit_netcdf(interbar, extrap, oldice, masque) ! !------------------------------------------------------------------------------- ! 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 control_mod USE indice_sol_mod #ifdef CPP_EARTH USE dimphy USE ioipsl, ONLY : ioget_year_len USE phys_state_var_mod, ONLY : pctsrf 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 #endif IMPLICIT NONE !------------------------------------------------------------------------------- ! Arguments: #include "dimensions.h" #include "paramet.h" #include "iniprint.h" LOGICAL, INTENT(IN) :: interbar ! barycentric interpolation LOGICAL, INTENT(IN) :: extrap ! SST extrapolation flag LOGICAL, INTENT(IN) :: oldice ! old way ice computation REAL, DIMENSION(iip1,jjp1), INTENT(IN) :: masque ! land mask #ifndef CPP_EARTH CALL abort_gcm('limit_netcdf','Earth-specific routine, needs Earth physics',1) #else !------------------------------------------------------------------------------- ! Local variables: #include "logic.h" #include "comvert.h" #include "comgeom2.h" #include "comconst.h" !--- INPUT NETCDF FILES NAMES -------------------------------------------------- CHARACTER(LEN=25) :: icefile, sstfile, dumstr CHARACTER(LEN=25), PARAMETER :: famipsst='amipbc_sst_1x1.nc ', & famipsic='amipbc_sic_1x1.nc ', & fcpldsst='cpl_atm_sst.nc ', & fcpldsic='cpl_atm_sic.nc ', & fhistsst='histmth_sst.nc ', & fhistsic='histmth_sic.nc ', & frugo ='Rugos.nc ', & falbe ='Albedo.nc ' CHARACTER(LEN=10) :: varname !--- OUTPUT VARIABLES FOR NETCDF FILE ------------------------------------------ REAL, DIMENSION(klon) :: fi_ice, verif REAL, DIMENSION(:,:), POINTER :: phy_rug=>NULL(), phy_ice=>NULL() REAL, DIMENSION(:,:), POINTER :: phy_sst=>NULL(), phy_alb=>NULL() REAL, DIMENSION(:,:), ALLOCATABLE :: phy_bil REAL, DIMENSION(:,:,:), ALLOCATABLE :: pctsrf_t INTEGER :: nbad !--- VARIABLES FOR OUTPUT FILE WRITING ----------------------------------------- INTEGER :: ierr, nid, ndim, ntim, k INTEGER, DIMENSION(2) :: dims INTEGER :: id_tim, id_SST, id_BILS, id_RUG, id_ALB INTEGER :: id_FOCE, id_FSIC, id_FTER, id_FLIC INTEGER :: NF90_FORMAT INTEGER :: ndays !--- Depending on the output calendar !--- INITIALIZATIONS ----------------------------------------------------------- #ifdef NC_DOUBLE NF90_FORMAT=NF90_DOUBLE #else NF90_FORMAT=NF90_FLOAT #endif pi = 4.*ATAN(1.) rad = 6371229. daysec= 86400. omeg = 2.*pi/daysec g = 9.8 kappa = 0.2857143 cpp = 1004.70885 dtvr = daysec/REAL(day_step) CALL inigeom !--- Beware: anneeref (from gcm.def) is used to determine output time sampling ndays=ioget_year_len(anneeref) !--- RUGOSITY TREATMENT -------------------------------------------------------- IF (prt_level>1) WRITE(lunout,*) 'Traitement de la rugosite' varname='RUGOS' CALL get_2Dfield(frugo,varname,'RUG',interbar,ndays,phy_rug,mask=masque(1:iim,:)) !--- OCEAN TREATMENT ----------------------------------------------------------- IF (prt_level>1) WRITE(lunout,*) 'Traitement de la glace oceanique' ! Input SIC file selection ! Open file only to test if available IF ( NF90_OPEN(TRIM(famipsic),NF90_NOWRITE,nid)==NF90_NOERR ) THEN icefile=TRIM(famipsic) varname='sicbcs' ELSE IF( NF90_OPEN(TRIM(fcpldsic),NF90_NOWRITE,nid)==NF90_NOERR ) THEN icefile=TRIM(fcpldsic) varname='SIICECOV' ELSE IF ( NF90_OPEN(TRIM(fhistsic),NF90_NOWRITE,nid)==NF90_NOERR ) THEN icefile=TRIM(fhistsic) varname='pourc_sic' ELSE WRITE(lunout,*) 'ERROR! No sea-ice input file was found.' WRITE(lunout,*) 'One of following files must be availible : ',trim(famipsic),', ',trim(fcpldsic),', ',trim(fhistsic) CALL abort_gcm('limit_netcdf','No sea-ice file was found',1) END IF ierr=NF90_CLOSE(nid) IF (prt_level>=0) WRITE(lunout,*)'Pour la glace de mer a ete choisi le fichier '//TRIM(icefile) CALL get_2Dfield(icefile,varname, 'SIC',interbar,ndays,phy_ice,flag=oldice) 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 point = ',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 ------------------------------------------------------------- IF (prt_level>1) WRITE(lunout,*) 'Traitement de la sst' ! Input SST file selection ! Open file only to test if available IF ( NF90_OPEN(TRIM(famipsst),NF90_NOWRITE,nid)==NF90_NOERR ) THEN sstfile=TRIM(famipsst) varname='tosbcs' ELSE IF ( NF90_OPEN(TRIM(fcpldsst),NF90_NOWRITE,nid)==NF90_NOERR ) THEN sstfile=TRIM(fcpldsst) varname='SISUTESW' ELSE IF ( NF90_OPEN(TRIM(fhistsst),NF90_NOWRITE,nid)==NF90_NOERR ) THEN sstfile=TRIM(fhistsst) varname='tsol_oce' ELSE WRITE(lunout,*) 'ERROR! No sst input file was found.' WRITE(lunout,*) 'One of following files must be availible : ',trim(famipsst),trim(fcpldsst),trim(fhistsst) CALL abort_gcm('limit_netcdf','No sst file was found',1) END IF ierr=NF90_CLOSE(nid) IF (prt_level>=0) WRITE(lunout,*)'Pour la temperature de mer a ete choisi un fichier '//TRIM(sstfile) CALL get_2Dfield(sstfile,varname,'SST',interbar,ndays,phy_sst,flag=extrap) !--- ALBEDO TREATMENT ---------------------------------------------------------- IF (prt_level>1) WRITE(lunout,*) 'Traitement de l albedo' varname='ALBEDO' CALL get_2Dfield(falbe,varname,'ALB',interbar,ndays,phy_alb) !--- REFERENCE GROUND HEAT FLUX TREATMENT -------------------------------------- ALLOCATE(phy_bil(klon,ndays)); phy_bil=0.0 !--- OUTPUT FILE WRITING ------------------------------------------------------- IF (prt_level>5) WRITE(lunout,*) 'Ecriture du fichier limit : debut' !--- File creation ierr=NF90_CREATE("limit.nc",NF90_CLOBBER,nid) ierr=NF90_PUT_ATT(nid,NF90_GLOBAL,"title","Fichier conditions aux limites") !--- Dimensions creation ierr=NF90_DEF_DIM(nid,"points_physiques",klon,ndim) ierr=NF90_DEF_DIM(nid,"time",NF90_UNLIMITED,ntim) dims=(/ndim,ntim/) !--- Variables creation ierr=NF90_DEF_VAR(nid,"TEMPS",NF90_FORMAT,(/ntim/),id_tim) ierr=NF90_DEF_VAR(nid,"FOCE", NF90_FORMAT,dims,id_FOCE) ierr=NF90_DEF_VAR(nid,"FSIC", NF90_FORMAT,dims,id_FSIC) ierr=NF90_DEF_VAR(nid,"FTER", NF90_FORMAT,dims,id_FTER) ierr=NF90_DEF_VAR(nid,"FLIC", NF90_FORMAT,dims,id_FLIC) ierr=NF90_DEF_VAR(nid,"SST", NF90_FORMAT,dims,id_SST) ierr=NF90_DEF_VAR(nid,"BILS", NF90_FORMAT,dims,id_BILS) ierr=NF90_DEF_VAR(nid,"ALB", NF90_FORMAT,dims,id_ALB) ierr=NF90_DEF_VAR(nid,"RUG", NF90_FORMAT,dims,id_RUG) !--- Attributes creation ierr=NF90_PUT_ATT(nid,id_tim, "title","Jour dans l annee") ierr=NF90_PUT_ATT(nid,id_FOCE,"title","Fraction ocean") ierr=NF90_PUT_ATT(nid,id_FSIC,"title","Fraction glace de mer") ierr=NF90_PUT_ATT(nid,id_FTER,"title","Fraction terre") ierr=NF90_PUT_ATT(nid,id_FLIC,"title","Fraction land ice") ierr=NF90_PUT_ATT(nid,id_SST ,"title","Temperature superficielle de la mer") ierr=NF90_PUT_ATT(nid,id_BILS,"title","Reference flux de chaleur au sol") ierr=NF90_PUT_ATT(nid,id_ALB, "title","Albedo a la surface") ierr=NF90_PUT_ATT(nid,id_RUG, "title","Rugosite") ierr=NF90_ENDDEF(nid) !--- Variables saving ierr=NF90_PUT_VAR(nid,id_tim,(/(REAL(k),k=1,ndays)/)) ierr=NF90_PUT_VAR(nid,id_FOCE,pctsrf_t(:,is_oce,:),(/1,1/),(/klon,ndays/)) ierr=NF90_PUT_VAR(nid,id_FSIC,pctsrf_t(:,is_sic,:),(/1,1/),(/klon,ndays/)) ierr=NF90_PUT_VAR(nid,id_FTER,pctsrf_t(:,is_ter,:),(/1,1/),(/klon,ndays/)) ierr=NF90_PUT_VAR(nid,id_FLIC,pctsrf_t(:,is_lic,:),(/1,1/),(/klon,ndays/)) ierr=NF90_PUT_VAR(nid,id_SST ,phy_sst(:,:),(/1,1/),(/klon,ndays/)) ierr=NF90_PUT_VAR(nid,id_BILS,phy_bil(:,:),(/1,1/),(/klon,ndays/)) ierr=NF90_PUT_VAR(nid,id_ALB ,phy_alb(:,:),(/1,1/),(/klon,ndays/)) ierr=NF90_PUT_VAR(nid,id_RUG ,phy_rug(:,:),(/1,1/),(/klon,ndays/)) ierr=NF90_CLOSE(nid) IF (prt_level>5) WRITE(lunout,*) 'Ecriture du fichier limit : fin' DEALLOCATE(pctsrf_t,phy_sst,phy_bil,phy_alb,phy_rug) !=============================================================================== ! CONTAINS ! !=============================================================================== !------------------------------------------------------------------------------- ! SUBROUTINE get_2Dfield(fnam, varname, mode, ibar, 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 dimphy, ONLY : klon USE phys_state_var_mod, ONLY : pctsrf USE control_mod 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" #include "iniprint.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 LOGICAL, INTENT(IN) :: ibar ! interp on pressure levels 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, DIMENSION(4) :: dids ! NetCDF dimensions identifiers REAL, ALLOCATABLE, DIMENSION(:) :: dlon_ini ! initial longitudes vector REAL, ALLOCATABLE, DIMENSION(:) :: dlat_ini ! initial latitudes vector REAL, POINTER, DIMENSION(:) :: dlon, dlat ! reordered lon/lat vectors !--- fields INTEGER :: imdep, jmdep, lmdep ! dimensions of 'champ' REAL, ALLOCATABLE, DIMENSION(:, :) :: champ ! wanted field on initial grid REAL, ALLOCATABLE, DIMENSION(:) :: yder, timeyear REAL, DIMENSION(iim, jjp1) :: champint ! interpolated field REAL, ALLOCATABLE, DIMENSION(:, :, :) :: champtime REAL, ALLOCATABLE, DIMENSION(:, :, :) :: 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, DIMENSION(:, :) :: work ! used for extrapolation CHARACTER(LEN=25) :: title ! for messages LOGICAL :: extrp ! flag for extrapolation LOGICAL :: oldice ! flag for old way ice computation 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. oldice=.FALSE. IF ( PRESENT(flag) ) THEN IF ( flag .AND. mode=='SST' ) extrp=.TRUE. IF ( flag .AND. mode=='SIC' ) oldice=.TRUE. END IF !--- GETTING SOME DIMENSIONAL VARIABLES FROM FILE ----------------------------- IF (prt_level>5) WRITE(lunout,*) ' Now reading file : ',fnam ierr=NF90_OPEN(fnam, NF90_NOWRITE, ncid); CALL ncerr(ierr, fnam) ierr=NF90_INQ_VARID(ncid, trim(varname), varid); CALL ncerr(ierr, fnam) ierr=NF90_INQUIRE_VARIABLE(ncid, varid, dimids=dids); CALL ncerr(ierr, fnam) !--- Read unit for sea-ice variable only IF (mode=='SIC') THEN ierr=NF90_GET_ATT(ncid, varid, 'units', unit_sic) IF(ierr/=NF90_NOERR) THEN IF (prt_level>5) WRITE(lunout,*) 'No unit was given in sea-ice file. Take percentage as default value' unit_sic='X' ELSE IF (prt_level>5) WRITE(lunout,*) ' Sea-ice cover has unit=',unit_sic END IF END IF !--- Longitude ierr=NF90_INQUIRE_DIMENSION(ncid, dids(1), name=dnam, len=imdep) CALL ncerr(ierr, fnam); ALLOCATE(dlon_ini(imdep), dlon(imdep)) ierr=NF90_INQ_VARID(ncid, dnam, varid); CALL ncerr(ierr, fnam) ierr=NF90_GET_VAR(ncid, varid, dlon_ini); CALL ncerr(ierr, fnam) IF (prt_level>5) WRITE(lunout, *) 'variable ', dnam, 'dimension ', imdep !--- Latitude ierr=NF90_INQUIRE_DIMENSION(ncid, dids(2), name=dnam, len=jmdep) CALL ncerr(ierr, fnam); ALLOCATE(dlat_ini(jmdep), dlat(jmdep)) ierr=NF90_INQ_VARID(ncid, dnam, varid); CALL ncerr(ierr, fnam) ierr=NF90_GET_VAR(ncid, varid, dlat_ini); CALL ncerr(ierr, fnam) IF (prt_level>5) WRITE(lunout, *) 'variable ', dnam, 'dimension ', jmdep !--- Time (variable is not needed - it is rebuilt - but calendar is) ierr=NF90_INQUIRE_DIMENSION(ncid, dids(3), name=dnam, len=lmdep) CALL ncerr(ierr, fnam); ALLOCATE(timeyear(lmdep)) ierr=NF90_INQ_VARID(ncid, dnam, varid); CALL ncerr(ierr, fnam) cal_in=' ' ierr=NF90_GET_ATT(ncid, varid, 'calendar', cal_in) IF(ierr/=NF90_NOERR) THEN SELECT CASE(mode) CASE('RUG', 'ALB'); cal_in='360d' CASE('SIC', 'SST'); cal_in='gregorian' END SELECT IF (prt_level>5) WRITE(lunout, *)'ATTENTION: variable "time" sans attribut "calendrier" ' & // 'dans '//TRIM(fnam)//'. On choisit la valeur par defaut.' END IF IF (prt_level>5) WRITE(lunout, *) 'variable ', dnam, 'dimension ', lmdep, 'calendrier ', & cal_in !--- 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)) IF (prt_level>5) WRITE(lunout, *) IF (prt_level>5) WRITE(lunout,*)'LECTURE ET INTERPOLATION HORIZ. DE ', lmdep, ' CHAMPS.' ierr=NF90_INQ_VARID(ncid, varname, varid); CALL ncerr(ierr, fnam) DO l=1, lmdep ierr=NF90_GET_VAR(ncid, varid, champ, (/1, 1, l/), (/imdep, jmdep, 1/)) CALL ncerr(ierr, fnam) CALL conf_dat2d(title, imdep, jmdep, dlon_ini, dlat_ini, dlon, dlat, & champ, ibar) IF (extrp) CALL extrapol(champ, imdep, jmdep, 999999., .TRUE., .TRUE., 2, & work) IF(ibar .AND. .NOT.oldice) THEN IF(l==1 .AND. prt_level>5) THEN WRITE(lunout, *) '-------------------------------------------------------------------------' WRITE(lunout, *) 'Utilisation de l''interpolation barycentrique pour ',TRIM(title),' $$$' WRITE(lunout, *) '-------------------------------------------------------------------------' 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 ELSE SELECT CASE(mode) CASE('RUG'); CALL rugosite(imdep, jmdep, dlon, dlat, champ, & iim, jjp1, rlonv, rlatu, champint, mask) CASE('SIC'); CALL sea_ice (imdep, jmdep, dlon, dlat, champ, & iim, jjp1, rlonv, rlatu, champint) CASE('SST', 'ALB'); CALL grille_m(imdep, jmdep, dlon, dlat, champ, & iim, jjp1, rlonv, rlatu, champint) END SELECT END IF champtime(:, :, l)=champint END DO ierr=NF90_CLOSE(ncid); CALL ncerr(ierr, 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 IF (prt_level>5) WRITE(lunout, *) 'Filtrage de la SST: SST >= 271.38' WHERE(champan<271.38) champan=271.38 END IF !--- SPECIAL FILTER FOR SIC: 0.05) WRITE(lunout, *) 'Filtrage de la SIC: 0.0 < Sea-ice < 1.0' IF (unit_sic=='1') THEN ! Nothing to be done for sea-ice field is already in fraction of 1 ! This is the case for sea-ice in file cpl_atm_sic.nc IF (prt_level>5) WRITE(lunout,*) 'Sea-ice field already in fraction of 1' ELSE ! Convert sea ice from percentage to fraction of 1 IF (prt_level>5) WRITE(lunout,*) 'Transformt sea-ice field from percentage to fraction of 1.' champan(:, :, :)=champan(:, :, :)/100. END IF champan(iip1, :, :)=champan(1, :, :) WHERE(champan>1.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 ! !------------------------------------------------------------------------------- !------------------------------------------------------------------------------- ! FUNCTION year_len(y,cal_in) ! !------------------------------------------------------------------------------- USE ioipsl, ONLY : ioget_calendar,ioconf_calendar,lock_calendar,ioget_year_len 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) ! !------------------------------------------------------------------------------- USE ioipsl, ONLY : ioget_calendar,ioconf_calendar,lock_calendar,ioget_mon_len 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_gcm('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 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 !------------------------------------------------------------------------------- #include "iniprint.h" IF(ncres/=NF90_NOERR) THEN WRITE(lunout,*)'Problem with file '//TRIM(fnam)//' in routine limit_netcdf.' CALL abort_gcm('limit_netcdf',NF90_STRERROR(ncres),1) END IF END SUBROUTINE ncerr ! !------------------------------------------------------------------------------- #endif ! of #ifdef CPP_EARTH END SUBROUTINE limit_netcdf