!-------------------------------------- LICENCE BEGIN ------------------------------------
!Environment Canada - Atmospheric Science and Technology License/Disclaimer,
! version 3; Last Modified: May 7, 2008.
!This is free but copyrighted software; you can use/redistribute/modify it under the terms
!of the Environment Canada - Atmospheric Science and Technology License/Disclaimer
!version 3 or (at your option) any later version that should be found at:
!http://collaboration.cmc.ec.gc.ca/science/rpn.comm/license.html
!
!This software is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
!without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
!See the above mentioned License/Disclaimer for more details.
!You should have received a copy of the License/Disclaimer along with this software;
!if not, you can write to: EC-RPN COMM Group, 2121 TransCanada, suite 500, Dorval (Quebec),
!CANADA, H9P 1J3; or send e-mail to service.rpn@ec.gc.ca
!-------------------------------------- LICENCE END --------------------------------------
***s/r itf_phy_inincr - define physics surface forcing increments
*
#include "model_macros_f.h"
*
subroutine itf_phy_inincr 1,3
*
implicit none
*
*author
* Bernard Dugas - july 97
*
*revision
* v2_21 - Bernard Dugas - Adaption to GEM DM v2.2 and to physics v3.67:
* 1) change all includes and variables to DM specs;
* 2) each PE does the interpolations for its own area;
* 3) INCRTS is only used for water points (TWATER);
* 4) INCRTG is only used for glacier points (TGLACIER(2));
* 5) INCRTP is only used for deep soil points when the
* force-restore surface scheme is selected (TSOIL(2));
* 6) INCRICD is added for sea ice depht increments (ICEDP);
* 7) INCRGL is now used with GLSEAS0, as GLSEA can be
* modified to account for lake ice evolution;
* 8) INCRAL is no longer defined, since CALCALB
* is used instead (called by CLIMPHS2).
* v2_31 - Desgagne M. - remove stkmemw
* v2_31 - Bernard Dugas - adapt to new climatological file descriptors
* v2_32 - Bernard Dugas - correct exit condition codes
* v3_00 - Desgagne & Lee - Lam configuration
* v3_01 - Bernard Dugas - add call to INTOZON
* v3_02 - Lubos Spacek - remplacer l_ni et l_nj par p_ni et p_nj
* v3_02 - Bernard Dugas - remettre le code de controle d'interpolation O3
* - toujours initialiser incrtp (independamment de soli_L)
* v3_10 - Lee V. - RPN_bcastc for bcast on MPI_CHARACTER
* v3_11 - Bernard Dugas - use glacier AND (maybe) sea ice SNODP to define INCRNE
* - Account for P_pbl_iceme_L = .true.
* v3_22 - Bernard Dugas - add check for Done_inincr and correct xncibl(5) = 12
* v3_22 - Katja Winger - make sure that BC are read twice in first month
* (at Lctl_step=2 and in the middle of the month)
* v3_30 - Dugas B. - new itf_phy interface
*
*object
*
* Initialize the climatological increments field that
* are in the permanent physics bus. Care has to be taken
* with the active physics configuration, especially with
* regards to the active surface parameterization.
*
* The increments are those required by the routine CLIMPHS2
* to simulate the evolution of some of the surface forcing
* terms accounted for in the RPN/CMC physics package.
*
* *** It is important to note that the physics package is
* *** supposed to have been initialized by a previous call
* *** to phy_init.
*
*implicits
#include "dcst.cdk"
#include "lun.cdk"
#include "mem.cdk"
#include "path.cdk"
#include "itf_phy_buses.cdk"
#include "itf_phy_busind.cdk"
#include "ptopo.cdk"
#include "lctl.cdk"
#include "rstr.cdk"
#include "out3.cdk"
#include "itf_phy_config.cdk"
#include "cstv.cdk"
#include "glb_ld.cdk"
*
*modules
INTEGER newdate,fnom,fclos,idatmg2,wkoffit,fstopc
EXTERNAL newdate,fnom,fclos,idatmg2,wkoffit,fstopc
INTEGER fstouv,fstinf,fstsui,fstlir,fstprm,fstfrm
EXTERNAL fstouv,fstinf,fstsui,fstlir,fstprm,fstfrm
**
*
REAL*8 heur_8,scal_8,rad_8
INTEGER xndate(14), xnday1 , xncibl(14), xnday2
EQUIVALENCE (xndate(14), xnday1),(xncibl(14), xnday2)
*
INTEGER xncfld
PARAMETER (xncfld=6)
CHARACTER gtyp_S(xncfld),gtypz_S(xncfld)
INTEGER xnzig1(xncfld),xnzig2(xncfld),xnzig3(xncfld)
INTEGER xnig1(xncfld), xnig2(xncfld), xnig3(xncfld)
INTEGER xnig4(xncfld), xnzig4(xncfld)
INTEGER xnilu(xncfld), xnjlu(xncfld)
CHARACTER*4 clim_S(xncfld)
*
LOGICAL acli_L,IntSol_L,IntIce_L
INTEGER xndate1,xndate2,xniz,xnjz
INTEGER Lun_clim,Lun_acli,Lun_activ,ppjour,aujour
INTEGER xnerr,xncle,xnhold,xnpos,xnzpos,xncoupe
INTEGER i,xni,xnj,xnk,xnim,xnjm,xnzm,xnmxfld
CHARACTER typvar_S,hold_S*12
CHARACTER*256 clima_S,anclima_S
CHARACTER*8 soli_S,soli2_S
*
integer Done_Inincr
save Done_Inincr
*
real xrwrk3, neige
pointer (iwork3,xrwrk3(1))
real xrlon(p_ni), xrlat(p_ni), xrwrk1(p_ni), xrwrk2(p_ni)
real, dimension (:), allocatable :: xrclim
real, dimension (:), allocatable :: xrtic,xrtac
real busper, busper2(max(1,p_bper_siz))
pointer (pabusper,busper(*))
*
DATA clim_S /'SD','HS','I8','LG','TP','TM'/
DATA Done_Inincr / -1 /
*
*
** -------------------------------------------------
** 0. Check the current date and time as the increments
** are only updated at the beginning of a run * OR *
** close to the 15th of each month at 00Z
** -------------------------------------------------
*
Lun_clim = 0
Lun_acli = 0
*
heur_8 = Lctl_step * (Cstv_dt_8 / 3600. )
call incdatr( xnday1,Out3_date,heur_8 )
call datmgp2( xndate )
*
** Update the climatological ozone field every 24 hours
** (and note that it is first done at timestep 1)
*
if (Cstv_dt_8 .lt. 43200. ) then
*
ppjour = nint( 86400. / Cstv_dt_8 )
aujour = mod( Lctl_step, ppjour )
*
if (Rstri_sdon .gt. 0 .and.
% aujour .eq. 1 )
% call intozon
( xndate(3), xndate(2), lun_out.gt.0 )
*
endif
*
** Check to see if the current model date and time
** are "close" to the 15th of the month at 12Z
*
xncibl(1) = 0
xncibl(2) = xndate(2)
xncibl(3) = 15
xncibl(4) = xndate(4)
xncibl(5) = 12
xncibl(6) = 00
*
xnday2 = idatmg2( xncibl )
call difdatr( xnday2,xnday1,heur_8 )
heur_8 = heur_8 * (3600. / Cstv_dt_8 )
*
if (((Done_Inincr.ne.xnday2) .and.
% (heur_8 .ge. -0.5 ) .and.
% (heur_8 .lt. +0.5 )) .or.
% ( Lctl_step .eq. 2 )) then
*
* Make sure we only do this once per target date
* in the middel of each month and at Lctl_step=2
if ( Lctl_step .ne. 2 ) Done_Inincr = xnday2
*
xnerr = newdate( xnday1, xndate1,xndate2, -3 )
if (Lun_out.gt.0) write(Lun_out,1000) xndate1,xndate2
*
if ( Ptopo_myproc .eq. 0 ) then
*
** 1. Processor 0 opens the climatological file
** and finds the maximum size of a grid
** -----------------------------------------
*
xnerr = 0
xnmxfld = 0
*
xnim = 0
xnjm = 0
xnzm = 0
*
clima_S = trim(Path_input_S)//'/imclima'
*
xnerr = wkoffit( clima_S )
*
if ( xnerr .eq. 1 .or. xnerr .eq. 33 ) then
*
xnerr = fstopc( 'MSGLVL','INFORM',.false. )
xnerr = fnom( Lun_clim, clima_S , 'STD+RND+OLD', 0 )
xnerr = fstouv( Lun_clim, 'RND' )
*
if ( xnerr .ge. 0 ) then
*
if (Lun_out.gt.0) write( Lun_out,* )
% 'Opening file: ', trim( clima_S )
*
xnerr = fstinf( Lun_clim, xni,xnj,xnk,
% -1, ' ', -1, -1, -1, ' ', ' ')
*
0010 if ( xnerr .ge. 0 ) then
xnim = max( xnim,xni )
xnjm = max( xnjm,xnj )
xnmxfld = xnim*xnjm
xnerr = fstsui( Lun_clim, xni,xnj,xnk )
goto 0010
endif
*
endif
*
endif
*
if (( xnmxfld .le. 0 ) .or.
% ( Lun_clim .eq. 0 )) then
if (Lun_out.gt.0) write( Lun_out,1001 ) trim( clima_S )
goto 0030
endif
*
** 2. Try opening the analysed/climatological file
** --------------------------------------------
*
xnerr = 0
acli_L = .false.
*
anclima_S = trim(Path_input_S)//'/imanclima'
*
xnerr = wkoffit( anclima_S )
*
if ( xnerr .eq. 1 .or. xnerr .eq. 33 ) then
*
xnerr = fnom( Lun_acli, anclima_S , 'STD+RND+OLD', 0 )
xnerr = fstouv( Lun_acli, 'RND' )
*
if ( xnerr .ge. 0 ) then
*
acli_L = .true.
*
if (Lun_out.gt.0) write( Lun_out,* )
% 'Opening file: ', trim( anclima_S )
*
xnerr = fstinf( Lun_acli, xni,xnj,xnk,
% -1, ' ', -1, -1, -1, ' ', ' ')
*
0020 if ( xnerr .ge. 0 ) then
xnim = max( xnim,xni )
xnjm = max( xnjm,xnj )
xnmxfld = xnim*xnjm
xnerr = fstsui( Lun_acli, xni,xnj,xnk )
goto 0020
endif
*
endif
*
endif
*
if ( .not. acli_L ) then
if (Lun_out.gt.0) write( Lun_out,1002 )
else
if (Lun_out.gt.0) write( Lun_out,1003 ) trim( anclima_S )
endif
*
xnzm = max( xnim,xnjm )
*
endif
*
0030 call RPN_COMM_bcast( Lun_clim ,1,"MPI_INTEGER",0,"grid",xnerr )
call RPN_COMM_bcast( xnmxfld ,1,"MPI_INTEGER",0,"grid",xnerr )
call RPN_COMM_bcast( xnzm ,1,"MPI_INTEGER",0,"grid",xnerr )
*
if (( xnmxfld .le. 0 ) .or. ( Lun_clim .eq. 0 ))
% call gem_stop
( 'itf_phy_inincr',-1 )
*
** Every one allocates the xncfld 2D climatology fields
*
allocate (xrtic(xnzm*xncfld),xrtac(xnzm*xncfld))
allocate (xrclim(xnmxfld*xncfld))
*
** What (Month/Year) is the target for the increments ?
*
xncibl(1) = 0
*
** Same month or next one ?
*
if ( xndate(3) .ge. 15 ) then
xncibl(2) = mod( xndate(2),12 )+1
else
xncibl(2) = xndate(2)
endif
*
xncibl(3) = 15
*
** Same year or next one ?
*
xncibl(4) = xndate(4)
if ( xncibl(2) .lt. xndate(2) )
% xncibl(4) = xncibl(4) + 1
*
xncibl(5) = 0
xncibl(6) = 0
*
xnday2 = idatmg2( xncibl )
call difdatr( xnday2,xnday1,heur_8 )
scal_8 = 1./(heur_8*(3600./Cstv_dt_8))
*
xnpos = 1
*
C print *,'Dans itf_phy_inincr: scal_8,xncibl=',scal_8,xncibl
C call flush( Lun_out )
*
rad_8 = 180./acos( -1.0D0 )
*
** Determine whether the surface land scheme is fully
** interactive or not. If it is, the TP increments will
** only be used for TGLACIER(2)
*
call low2up( P_pbl_schsl_S,soli_S )
soli2_S = soli_S
if (P_pset_second_L) call low2up( P_pbl2_schsl_S,soli2_S )
if ((soli_S .eq.'CLASS' .or. soli_S .eq.'ISBA') .and.
% (soli2_S .eq.'CLASS' .or. soli2_S .eq.'ISBA')) then
IntSol_L = .true.
else
IntSol_L = .false.
endif
*
** Check for interactive snow (in glacier & sea ice
** modules) and ice thichness (in sea ice module)
*
IntIce_L = P_pbl_iceme_L
if (P_pset_second_L) IntIce_L = (IntIce_L .and. P_pbl2_iceme_L)
*
C print *,'Dans itf_phy_inincr: soli_S ',soli_S
C call flush( Lun_out )
*
if ( Ptopo_myproc .eq. 0 ) then
*
** 3. Processor 0 reads the target month's climatic fields
** ----------------------------------------------------
*
do i=1,xncfld
*
if ((clim_S(i).eq.'SD' .and. IntIce_L .and. IntSol_L)
% .or. (clim_S(i).eq.'I8' .and. IntIce_L)
% .or. (clim_S(i).eq.'HS' .and. IntSol_L)) then
*
if (Lun_out.gt.0) write( Lun_out,1004 ) clim_S(i)
*
else
*
xncle = -1
*
** Always try reading from Lun_acli first, if possible
*
if (acli_L)
% xncle = fstlir( xrclim(xnpos), Lun_acli,
% xnilu(i),xnjlu(i), xnk, -1,' ',-1,
% xncibl(2),xncibl(4),'A',clim_S(i) )
*
** Otherwise read from Lun_clim
*
if ( xncle .lt. 0 )
% xncle = fstlir( xrclim(xnpos), Lun_clim,
% xnilu(i),xnjlu(i), xnk, -1,' ',
% -1,xncibl(2),-1,'C',clim_S(i) )
*
if ( xncle .lt. 0 ) then
if (Lun_out.gt.0) write( Lun_out,1005 ) clim_S(i)
goto 0040
endif
*
** Save the grid description parameters
*
xnerr = fstprm( xncle, xnhold,xnhold,xnhold,
% xnilu(i),xnjlu(i),xnk,
% xnhold,xnhold,xnhold,xnhold,xnhold,
% typvar_S,clim_S(i),hold_S,
% gtyp_S(i),xnig1(i),xnig2(i),
% xnig3(i),xnig4(i),
% xnhold,xnhold,xnhold,xnhold,
% xnhold,xnhold,xnhold )
*
if (gtyp_S(i).eq.'Z') then
*
** Retreive necessary Z-grid descriptors
*
if (typvar_S.eq.'A') then
Lun_activ = Lun_acli
else if (typvar_S.eq.'C') then
Lun_activ = Lun_clim
endif
*
xnzpos = (i-1)*xnzm+1
*
xncle = fstlir( xrtic(xnzpos), Lun_activ,
% xniz,xnjz,xnhold, -1,' ',
% xnig1(i),xnig2(i),xnig3(i),
% typvar_S,'>>')
*
if ( xncle .lt. 0 .or.
% xniz .ne. xnilu(i) ) then
if (Lun_out.gt.0) write( Lun_out,1006 )
xncle = -1
goto 0040
endif
*
xncle = fstlir( xrtac(xnzpos), Lun_activ,
% xniz,xnjz,xnhold, -1,' ',
% xnig1(i),xnig2(i),xnig3(i),
% typvar_S,'^^')
*
if ( xncle .lt. 0 .or.
% xnjz .ne. xnjlu(i) ) then
if (Lun_out.gt.0) write( Lun_out,1007 )
xncle = -1
goto 0040
endif
*
xnerr = fstprm( xncle, xnhold,xnhold,xnhold,
% xnhold,xnhold,xnhold,
% xnhold,xnhold,xnhold,xnhold,xnhold,
% hold_S,hold_S,hold_S,
% gtypz_S(i),xnzig1(i),xnzig2(i),
% xnzig3(i),xnzig4(i),
% xnhold,xnhold,xnhold,xnhold,
% xnhold,xnhold,xnhold )
endif
*
C print *,'Dans itf_phy_inincr: lu ',clim_S(i)
C call flush( Lun_out )
*
endif
*
xnpos = xnpos+xnmxfld
*
enddo
*
xnerr = fstfrm( Lun_clim )
xnerr = fclos ( Lun_clim )
*
if ( acli_L ) then
xnerr = fstfrm( Lun_acli )
xnerr = fclos ( Lun_acli )
endif
*
xncle = 0
endif
*
0040 call RPN_COMM_bcast( xncle ,1,"MPI_INTEGER",0,"grid",xnerr )
if (xncle.lt.0) call gem_stop
( 'itf_phy_inincr',-1 )
*
call RPN_COMM_bcast (xrtic ,xncfld*xnzm ,"MPI_REAL" ,0,"grid",xnerr)
call RPN_COMM_bcast (xrtac ,xncfld*xnzm ,"MPI_REAL" ,0,"grid",xnerr)
call RPN_COMM_bcast (xrclim ,xncfld*xnmxfld,"MPI_REAL" ,0,"grid",xnerr)
call RPN_COMM_bcast (xnilu ,xncfld ,"MPI_INTEGER" ,0,"grid",xnerr)
call RPN_COMM_bcast (xnjlu ,xncfld ,"MPI_INTEGER" ,0,"grid",xnerr)
call RPN_COMM_bcast (xnig1 ,xncfld ,"MPI_INTEGER" ,0,"grid",xnerr)
call RPN_COMM_bcast (xnig2 ,xncfld ,"MPI_INTEGER" ,0,"grid",xnerr)
call RPN_COMM_bcast (xnig3 ,xncfld ,"MPI_INTEGER" ,0,"grid",xnerr)
call RPN_COMM_bcast (xnig4 ,xncfld ,"MPI_INTEGER" ,0,"grid",xnerr)
call RPN_COMM_bcastc(gtyp_S ,xncfld ,"MPI_CHARACTER",0,"grid",xnerr)
call RPN_COMM_bcast (xnzig1 ,xncfld ,"MPI_INTEGER" ,0,"grid",xnerr)
call RPN_COMM_bcast (xnzig2 ,xncfld ,"MPI_INTEGER" ,0,"grid",xnerr)
call RPN_COMM_bcast (xnzig3 ,xncfld ,"MPI_INTEGER" ,0,"grid",xnerr)
call RPN_COMM_bcast (xnzig4 ,xncfld ,"MPI_INTEGER" ,0,"grid",xnerr)
call RPN_COMM_bcastc(gtypz_S,xncfld ,"MPI_CHARACTER",0,"grid",xnerr)
call RPN_COMM_bcastc(clim_S ,xncfld*4 ,"MPI_CHARACTER",0,"grid",xnerr)
*
** 4. Loop on the p_nj physics slabs.
** The increments are updated in place
** -----------------------------------
*
do xncoupe=1,p_nj
*
Mem_pslic = xncoupe
*
pabusper = loc (Phy_busper3D((xncoupe-1)*p_bper_siz+1))
*
** Convert Busp_idlat,Busp_idlon from radians to degrees
*
do i=0,p_ni-1
xrlat(i+1) = busper(dlat+i)*rad_8
xrlon(i+1) = busper(dlon+i)*rad_8
enddo
*
** Interpolate the climatological
** fields to the current slab
*
call ez_rgoptc( 'INTERP','VOISIN',.true. )
*
** Start with the snow depth. It is used over glaciers and sea-ice
** when ICEMELT is false. Soil uses it only with the Force-restore
** scheme. Either the maximum GLACIER and GLSEA snow depth (which
** are currently in the second and fourth rows) or the soil snow
** depth is used here.
*
xnzpos = 1
xnpos = 1
iwork3 = loc( xrclim(xnpos) )
*
C xnpos = 0
*
if (.not.IntIce_L .or. .not.IntSol_L) then
*
if (gtyp_S(1).ne.'Z') then
call ez_rgscint( xrwrk1,p_ni,1,xrlat,xrlon,
% xrwrk3,xnilu(1),xnjlu(1),
% gtyp_S(1),xnig1(1),xnig2(1),
% xnig3(1),xnig4(1),.TRUE. )
c call ez_rgfree
else
call ez_igscint( xrwrk1,p_ni,1,xrlat,xrlon,
% xrwrk3,xnilu(1),xnjlu(1),'Z',
% gtypz_S(1),xnzig1(1),xnzig2(1),
% xnzig3(1),xnzig4(1),.TRUE.,
% xrtic(xnzpos),xrtac(xnzpos) )
endif
*
if (.not.IntIce_L) then
do i=0,p_ni-1
neige = max( busper(snodp+ p_ni+i) ,
% busper(snodp+3*p_ni+i) )
busper(incrne+i) = (xrwrk1(i+1) - neige) * scal_8
enddo
else
do i=0,p_ni-1
busper(incrne+i) = (xrwrk1(i+1) -busper(snodp+i))
% * scal_8
enddo
endif
*
C print *,'Dans itf_phy_inincr: Busp_snodp interpole/initialise'
C call flush( Lun_out )
*
endif
*
** Continue with the soil moisture (again, only for Force-restore) ...
*
xnzpos = xnzpos+xnzm
xnpos = xnpos+xnmxfld
iwork3 = loc( xrclim(xnpos) )
*
if (.not.IntSol_L) then
*
if (gtyp_S(2).ne.'Z') then
call ez_rgscint( xrwrk1,p_ni,1,xrlat,xrlon,
% xrwrk3,xnilu(2),xnjlu(2),
% gtyp_S(2),xnig1(2),xnig2(2),
% xnig3(2),xnig4(2),.TRUE. )
c call ez_rgfree
else
call ez_igscint( xrwrk1,p_ni,1,xrlat,xrlon,
% xrwrk3,xnilu(2),xnjlu(2),'Z',
% gtypz_S(2),xnzig1(2),xnzig2(2),
% xnzig3(2),xnzig4(2),.TRUE.,
% xrtic(xnzpos),xrtac(xnzpos) )
endif
*
do i=0,p_ni-1
busper(incrhs+i) = (xrwrk1(i+1) - busper(wsoil+i))
% * scal_8
enddo
*
C print *,'Dans itf_phy_inincr: Busp_wsoil interpole/initialise'
C call flush( Lun_out )
*
endif
*
** ... With the sea ice thickness (if ICEMELT is false) ...
*
xnzpos = xnzpos+xnzm
xnpos = xnpos+xnmxfld
iwork3 = loc( xrclim(xnpos) )
*
if (.not.IntIce_L) then
*
if (gtyp_S(3).ne.'Z') then
call ez_rgscint( xrwrk1,p_ni,1,xrlat,xrlon,
% xrwrk3,xnilu(3),xnjlu(3),
% gtyp_S(3),xnig1(3),xnig2(3),
% xnig3(3),xnig4(3),.TRUE. )
c call ez_rgfree
else
call ez_igscint( xrwrk1,p_ni,1,xrlat,xrlon,
% xrwrk3,xnilu(3),xnjlu(3),'Z',
% gtypz_S(3),xnzig1(3),xnzig2(3),
% xnzig3(3),xnzig4(3),.TRUE.,
% xrtic(xnzpos),xrtac(xnzpos) )
endif
*
do i=0,p_ni-1
busper(incricd+i) = (xrwrk1(i+1) - busper(icedp+i))
% * scal_8
enddo
*
C print *,'Dans itf_phy_inincr: Busp_icedp interpole/initialise'
C call flush( Lun_out )
*
endif
*
** ... And with the sea ice mask
*
xnzpos = xnzpos+xnzm
xnpos = xnpos+xnmxfld
iwork3 = loc( xrclim(xnpos) )
*
if (gtyp_S(4).ne.'Z') then
call ez_rgscint( xrwrk1,p_ni,1,xrlat,xrlon,
% xrwrk3,xnilu(4),xnjlu(4),
% gtyp_S(4),xnig1(4),xnig2(4),
% xnig3(4),xnig4(4),.TRUE. )
c call ez_rgfree
else
call ez_igscint( xrwrk1,p_ni,1,xrlat,xrlon,
% xrwrk3,xnilu(4),xnjlu(4),'Z',
% gtypz_S(4),xnzig1(4),xnzig2(4),
% xnzig3(4),xnzig4(4),.TRUE.,
% xrtic(xnzpos),xrtac(xnzpos) )
endif
*
do i=0,p_ni-1
busper(incrgl+i) = (xrwrk1(i+1) - busper(glsea0+i))
% * scal_8
enddo
*
C print *,'Dans itf_phy_inincr: Busp_glsea0 interpole/initialise'
C call flush( Lun_out )
*
call ez_rgoptc( 'INTERP','CUBIQUE',.true. )
*
** Define the two temperature increments
*
** Start by interpolating the deep soil temperature
** TP into xrwrk2. TP will be used by Force-restore
** and by the glacier module.
*
xnzpos = xnzpos+xnzm
xnpos = xnpos+xnmxfld
iwork3 = loc( xrclim(xnpos) )
*
if (gtyp_S(5).ne.'Z') then
call ez_rgscint( xrwrk2,p_ni,1,xrlat,xrlon,
% xrwrk3,xnilu(5),xnjlu(5),
% gtyp_S(5),xnig1(5),xnig2(5),
% xnig3(5),xnig4(5),.TRUE. )
c call ez_rgfree
else
call ez_igscint( xrwrk2,p_ni,1,xrlat,xrlon,
% xrwrk3,xnilu(5),xnjlu(5),'Z',
% gtypz_S(5),xnzig1(5),xnzig2(5),
% xnzig3(5),xnzig4(5),.TRUE.,
% xrtic(xnzpos),xrtac(xnzpos) )
endif
*
** The sea surface temperature TM goes into xrwrk1
*
xnzpos = xnzpos+xnzm
xnpos = xnpos+xnmxfld
iwork3 = loc( xrclim(xnpos) )
*
if (gtyp_S(6).ne.'Z') then
call ez_rgscint( xrwrk1,p_ni,1,xrlat,xrlon,
% xrwrk3,xnilu(6),xnjlu(6),
% gtyp_S(6),xnig1(6),xnig2(6),
% xnig3(6),xnig4(6),.TRUE. )
c call ez_rgfree
else
call ez_igscint( xrwrk1,p_ni,1,xrlat,xrlon,
% xrwrk3,xnilu(6),xnjlu(6),'Z',
% gtypz_S(6),xnzig1(6),xnzig2(6),
% xnzig3(6),xnzig4(6),.TRUE.,
% xrtic(xnzpos),xrtac(xnzpos) )
endif
*
** Convert the new temperatures from Celsius to Kelvin
*
do i=1,p_ni
xrwrk1(i) = xrwrk1(i) + Dcst_tcdk_8
xrwrk2(i) = xrwrk2(i) + Dcst_tcdk_8
enddo
*
** Busp_incrts is used for SSTs, Busp_incrtg is used in the
** second layer Continental Ice temperature and Busp_incrtp
** is used with the Force-Restore.
*
do i=0,p_ni-1
busper(incrts+i) = (xrwrk1(i+1) - busper(twater+i))
% * scal_8
busper(incrtg+i) = (min( Dcst_trpl_8,
% 1.0d0*xrwrk2(i+1)) - busper(tglacier+p_ni+i)) * scal_8
busper(incrtp+i) = (xrwrk2(i+1) - busper(tsoil+p_ni+i))
% * scal_8
enddo
*
C print *,'Dans itf_phy_inincr: TMP/TS interpole/initialise'
C call flush( Lun_out )
*
** 5. Write the modified slice
** ------------------------
*
enddo
*
** De-allocate working memory
*
deallocate (xrclim)
*
endif
*
1000 format(
+//,'INITIALIZATION OF PHYSICS INCREMENTS (S/R ITF_PHY_ININCR)',
+ /,'=================================================',
+//,' The current date is ',I8,'/',I8.8,'(YYYYMMDD/HHMMSShh)',
+//)
*
1001 format(' Unable to open climatological file ',A)
1002 format(' Analysed/climatological file not opened')
1003 format(' Using analysed/climatological file ',A)
1004 format(' Climatological ',A4,' will not be read')
1005 format(' Unable to read climatological ',A4)
1006 format(' Unable to read >> record')
1007 format(' Unable to read ^^ record')
1008 format(' Physics not initialized before itf_phy_inincr')
return
end