Browse Source

Receive and send fluxes from and to the flux calculator.

distribute-radiation-on-surface-types
Sven Karsten 2 years ago
parent
commit
9be54516b0
  1. 54
      src/oasis_interface/oas_define.F90
  2. 156
      src/oasis_interface/oas_exchange_fields.F90
  3. 24
      src/oasis_interface/oas_receive_field.F90
  4. 4
      src/oasis_interface/oas_vardef.F90

54
src/oasis_interface/oas_define.F90

@ -215,10 +215,63 @@ SUBROUTINE oas_define(local_comm)
ssnd(10)%var_name = MOM5_instance_letter//'STSUR04'
ssnd(11)%var_name = MOM5_instance_letter//'STSUR05'
ssnd(12)%var_name = MOM5_instance_letter//'STSUR06'
ssnd(13)%var_name = MOM5_instance_letter//'SALBE01' ! shortwave albedo of this ice class (1=open sea, 2..6=ice)
ssnd(14)%var_name = MOM5_instance_letter//'SALBE02'
ssnd(15)%var_name = MOM5_instance_letter//'SALBE03'
ssnd(16)%var_name = MOM5_instance_letter//'SALBE04'
ssnd(17)%var_name = MOM5_instance_letter//'SALBE05'
ssnd(18)%var_name = MOM5_instance_letter//'SALBE06'
#endif
! Store info which fields to receive
ALLOCATE( srcv(nfld_rcv_tot), stat = ierror )
#ifdef OASIS_IOW_ESM
srcv( 1)%var_name = 'MRMRAI01' ! precipitation
srcv( 2)%var_name = 'MRMEVA01' ! evaporation
srcv( 3)%var_name = 'MRMEVA02' ! evaporation
srcv( 4)%var_name = 'MRMEVA03' ! evaporation
srcv( 5)%var_name = 'MRMEVA04' ! evaporation
srcv( 6)%var_name = 'MRMEVA05' ! evaporation
srcv( 7)%var_name = 'MRMEVA06' ! evaporation
srcv( 8)%var_name = 'MRMSNO01' ! snow
srcv( 9)%var_name = 'MRPSUR01' ! sea level pressure
srcv(10)%var_name = 'MRU10M01' ! u velocity 10m
srcv(11)%var_name = 'MRV10M01' ! v velocity 10m
srcv(12)%var_name = 'MRUMOM01' ! u wind stress
srcv(13)%var_name = 'MRUMOM02' ! u wind stress
srcv(14)%var_name = 'MRUMOM03' ! u wind stress
srcv(15)%var_name = 'MRUMOM04' ! u wind stress
srcv(16)%var_name = 'MRUMOM05' ! u wind stress
srcv(17)%var_name = 'MRUMOM06' ! u wind stress
srcv(18)%var_name = 'MRVMOM01' ! v wind stress
srcv(19)%var_name = 'MRVMOM02' ! v wind stress
srcv(20)%var_name = 'MRVMOM03' ! v wind stress
srcv(21)%var_name = 'MRVMOM04' ! v wind stress
srcv(22)%var_name = 'MRVMOM05' ! v wind stress
srcv(23)%var_name = 'MRVMOM06' ! v wind stress
srcv(24)%var_name = 'MRRBBR01' ! longwave radiation upward
srcv(25)%var_name = 'MRRBBR02' ! longwave radiation upward
srcv(26)%var_name = 'MRRBBR03' ! longwave radiation upward
srcv(27)%var_name = 'MRRBBR04' ! longwave radiation upward
srcv(28)%var_name = 'MRRBBR05' ! longwave radiation upward
srcv(29)%var_name = 'MRRBBR06' ! longwave radiation upward
srcv(30)%var_name = 'MRRLWD01' ! longwave radiation downward
srcv(31)%var_name = 'MRRSDD01' ! shortwave radiation downward direct
srcv(32)%var_name = 'MRRSIN01' ! shortwave radiation downward diffusive
srcv(33)%var_name = 'MRHLAT01' ! latent heat flux
srcv(34)%var_name = 'MRHLAT02' ! latent heat flux
srcv(35)%var_name = 'MRHLAT03' ! latent heat flux
srcv(36)%var_name = 'MRHLAT04' ! latent heat flux
srcv(37)%var_name = 'MRHLAT05' ! latent heat flux
srcv(38)%var_name = 'MRHLAT06' ! latent heat flux
srcv(39)%var_name = 'MRHSEN01' ! sensible heat flux
srcv(40)%var_name = 'MRHSEN02' ! sensible heat flux
srcv(41)%var_name = 'MRHSEN03' ! sensible heat flux
srcv(42)%var_name = 'MRHSEN04' ! sensible heat flux
srcv(43)%var_name = 'MRHSEN05' ! sensible heat flux
srcv(44)%var_name = 'MRHSEN06' ! sensible heat flux
#else
srcv( 1)%var_name = 'OCEREPRE' ! precipitation
srcv( 2)%var_name = 'OCEREEVA' ! evaporation
srcv( 3)%var_name = 'OCERESNO' ! snow
@ -233,6 +286,7 @@ SUBROUTINE oas_define(local_comm)
srcv(12)%var_name = 'OCERESWF' ! shortwave radiation downward diffusive
srcv(13)%var_name = 'OCERELHF' ! latent heat flux
srcv(14)%var_name = 'OCERESHF' ! sensible heat flux
#endif
! Allocate field for temporary storage of coupler data
ALLOCATE( exfld1(isc:iec, jsc:jec), stat = ierror )

156
src/oasis_interface/oas_exchange_fields.F90

@ -16,19 +16,27 @@ USE ice_model_mod, only: ice_data_type,atmos_ice_boundary_type
USE ocean_model_mod, only: ice_ocean_boundary_type
USE time_manager_mod, only : get_time,get_time,operator(-),time_type
USE netcdf
IMPLICIT NONE
!
! local parameters, variables and arrays
!
INTEGER :: &
maskt(iec-isc+1,jec-jsc+1), & ! mask array
maskc(iec-isc+1,jec-jsc+1) ! mask array
INTEGER :: &
isec, &
kinfo, &
jn, &
sc,dy,dummy
ncfileid, ncvarid, & ! NetCDF IDs
sc,dy,k,dummy, &
nrcvinfo(nfld_rcv_tot) ! OASIS info argument
REAL (KIND=8) :: &
ztmp1 (iec,jec)
ztmp1 (iec,jec)
TYPE (ice_data_type), INTENT(IN) :: Ice
TYPE (atmos_ice_boundary_type), INTENT(INOUT) :: Ice_boundary
@ -36,6 +44,9 @@ type(ice_ocean_boundary_type), intent(inout) :: Ice_Ocean_Boundary
type(time_type), intent(in) :: Time_start,Timet ! for coupling sandra
integer, parameter :: dt_cpld=600
ztmp1 (:,:) = 0
nrcvinfo (:) = OASIS_idle
call get_time(Timet-Time_start, sc, dy)
isec = (864e2*dy+sc)
@ -58,6 +69,147 @@ integer, parameter :: dt_cpld=600
ztmp1(isc:iec,jsc:jec) =Ice%t_surf(isc:iec,jsc:jec,jn)
CALL oas_send (jn+6, isec,ztmp1, kinfo ) ! send MSTSUR01..MSTSUR06
ENDDO
write(*,*) 'Sending MSALBE01..MSALBE06 at isec=',isec
DO jn=1,6
ztmp1(isc:iec,jsc:jec) =Ice%albedo_vis_dir(isc:iec,jsc:jec,jn)
CALL oas_send (jn+12, isec,ztmp1, kinfo ) ! send MSALBE01..MSALBE06
ENDDO
!----------------------------------------------------------------------------
! STEP 2: Receive all fluxes from flux calculator
!----------------------------------------------------------------------------
! 2.1 Receive fields via coupler from flux calculator
DO jn = 1, nfld_rcv_tot
! IF( srcv(jn)%laction ) THEN
CALL oas_recieve( jn, isec, ztmp1(:,:), nrcvinfo(jn) )
!write(*,*) " nrcvinfo(jn) ", nrcvinfo(jn), jn, OASIS_Rcv, isec
IF( nrcvinfo(jn) == OASIS_Rcv ) frcv(:,:,jn)=ztmp1(:,:)
! ENDIF
ENDDO
!write(*,*)" oasis receive laction: " , srcv(jn)%laction, isec
!IF (ltime) CALL get_timings (i_cpl_get, ntstep, dt, izerror)
! 2.2 Obtain coupling mask for fileds
istatus=nf90_open('masks.nc', NF90_NOWRITE, ncfileid)
istatus=nf90_inq_varid(ncfileid, 'tmom.msk' , ncvarid)
istatus=nf90_get_var(ncfileid, ncvarid, maskt, &
(/ isc, jsc /), (/ iec-isc+1,jec-jsc+1 /))
istatus=nf90_inq_varid(ncfileid, 'cmom.msk' , ncvarid)
istatus=nf90_get_var(ncfileid, ncvarid, maskc, &
(/ isc, jsc /), (/ iec-isc+1,jec-jsc+1 /))
istatus=nf90_close(ncfileid)
! 2.3 Store the received fluxes
jn = 1
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! rain 1)
do k = 1, size(Ice_boundary%lprec,3)
WHERE (maskt == 0) Ice_boundary%lprec(isc:iec,jsc:jec,k) = -frcv(isc:iec,jsc:jec,jn)
enddo
!ENDIF
jn = 2
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! evaporation 2..7)
do k = 1, size(Ice_boundary%q_flux,3)
WHERE (maskt == 0) Ice_boundary%q_flux(isc:iec,jsc:jec,k) = frcv(isc:iec,jsc:jec,jn)
jn = jn + 1 ! count up for different surface types
enddo
!ENDIF
jn = 8
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! snow 8)
do k = 1, size(Ice_boundary%fprec,3)
WHERE (maskt == 0) Ice_boundary%fprec(isc:iec,jsc:jec,k) = -frcv(isc:iec,jsc:jec,jn)
enddo
!ENDIF
jn = 9
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! sea level pressure 9)
do k = 1, size(Ice_boundary%p,3)
WHERE (maskt == 0) Ice_boundary%p(isc:iec,jsc:jec,k) = frcv(isc:iec,jsc:jec,jn)
enddo
!ENDIF
jn = 10
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! u velocity 10)
WHERE (maskc == 0) Ice_Ocean_Boundary%u_wind(isc:iec,jsc:jec) = frcv(isc:iec,jsc:jec,jn)
!ENDIF
jn = 11
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! v velocity 11)
WHERE (maskc == 0) Ice_Ocean_Boundary%v_wind(isc:iec,jsc:jec) = frcv(isc:iec,jsc:jec,jn)
!ENDIF
jn = 12
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! u wind stress 12..17)
do k = 1, size(Ice_boundary%u_flux,3)
WHERE (maskc == 0) Ice_boundary%u_flux(isc:iec,jsc:jec,k) = frcv(isc:iec,jsc:jec,jn)
jn = jn + 1
enddo
!ENDIF
jn = 18
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! v wind stress 18..23)
do k = 1, size(Ice_boundary%v_flux,3)
WHERE (maskc == 0) Ice_boundary%v_flux(isc:iec,jsc:jec,k) = frcv(isc:iec,jsc:jec,jn)
jn = jn + 1
enddo
!ENDIF
jn = 24
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! longwave radiation upward 24..29)
do k = 1, size(Ice_boundary%lw_flux,3)
WHERE (maskt == 0) Ice_boundary%lw_flux(isc:iec,jsc:jec,k) = frcv(isc:iec,jsc:jec,jn)
jn = jn + 1
enddo
!ENDIF
jn = 30
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! longwave radiation downward 30)
do k = 1, size(Ice_boundary%lw_flux,3)
WHERE (maskt == 0) Ice_boundary%lw_flux(isc:iec,jsc:jec,k) = -Ice_boundary%lw_flux(isc:iec,jsc:jec,k) - frcv(isc:iec,jsc:jec,jn)
enddo
!ENDIF
jn = 31
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! shortwave radiation direct 31)
do k = 1, size(Ice_boundary%sw_flux_vis_dir,3)
WHERE (maskt == 0) Ice_boundary%sw_flux_vis_dir(isc:iec,jsc:jec,k) = -frcv(isc:iec,jsc:jec,jn)
enddo
!ENDIF
jn = 32
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! shortware radiation diffusive 32)
do k = 1, size(Ice_boundary%sw_flux_vis_dif,3)
WHERE (maskt == 0) Ice_boundary%sw_flux_vis_dif(isc:iec,jsc:jec,k) = -frcv(isc:iec,jsc:jec,jn)
enddo
!ENDIF
jn = 33
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! latent heat flux 33..38)
do k = 1, size(Ice_boundary%lh_flux,3)
WHERE (maskt == 0) Ice_boundary%lh_flux(isc:iec,jsc:jec,k) = frcv(isc:iec,jsc:jec,jn)
jn = jn + 1
enddo
!ENDIF
jn = 39
!IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! sensible heat flux 39..44)
do k = 1, size(Ice_boundary%t_flux,3)
WHERE (maskt == 0) Ice_boundary%t_flux(isc:iec,jsc:jec,k) = frcv(isc:iec,jsc:jec,jn)
jn = jn + 1
enddo
!ENDIF
write(*,*) 'oas_exchange_fields finished.'
END SUBROUTINE oas_exchange_fields
#ENDIF

24
src/oasis_interface/oas_receive_field.F90

@ -53,20 +53,20 @@ REAL(kind=8) :: &
!-------------------------------------------------------------------------------
#ifdef OASIS_IOW_ESM
! Do not read here, use instead oas_exchange_fields module
#else
DO jn = 1, nfld_rcv_tot
! IF( srcv(jn)%laction ) THEN
CALL oas_recieve( jn, isec, ztmp1(:,:), nrcvinfo(jn) )
! write(*,*) " nrcvinfo(jn) ", nrcvinfo(jn), jn, OASIS_Rcv, isec
IF( nrcvinfo(jn) == OASIS_Rcv ) frcv(:,:,jn)=ztmp1(:,:)
! ENDIF
ENDDO
! IF( srcv(jn)%laction ) THEN
CALL oas_recieve( jn, isec, ztmp1(:,:), nrcvinfo(jn) )
! write(*,*) " nrcvinfo(jn) ", nrcvinfo(jn), jn, OASIS_Rcv, isec
IF( nrcvinfo(jn) == OASIS_Rcv ) frcv(:,:,jn)=ztmp1(:,:)
! ENDIF
ENDDO
!write(*,*)" oasis receive laction: " , srcv(jn)%laction, isec
!IF (ltime) CALL get_timings (i_cpl_get, ntstep, dt, izerror)
istatus=nf90_open('masks.nc', NF90_NOWRITE, ncfileid)
istatus=nf90_inq_varid(ncfileid, 'tmom.msk' , ncvarid)
istatus=nf90_get_var(ncfileid, ncvarid, maskt, &
@ -86,8 +86,8 @@ REAL(kind=8) :: &
jn = jn + 1
! IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! evaporation 2)
do k = 1, size(Ice_boundary%q_flux,3)
WHERE (maskt == 0) Ice_boundary%q_flux(isc:iec,jsc:jec,k) = -frcv(isc:iec,jsc:jec,jn)
enddo
WHERE (maskt == 0) Ice_boundary%q_flux(isc:iec,jsc:jec,k) = -frcv(isc:iec,jsc:jec,jn)
enddo
! ENDIF
jn = jn + 1
! IF( nrcvinfo(jn) == OASIS_Rcv ) THEN ! snow 3)
@ -157,7 +157,7 @@ REAL(kind=8) :: &
WHERE (maskt == 0) Ice_boundary%t_flux(isc:iec,jsc:jec,k) = -frcv(isc:iec,jsc:jec,jn)
enddo
! ENDIF
#endif
END SUBROUTINE oas_receive_field
#ENDIF

4
src/oasis_interface/oas_vardef.F90

@ -40,8 +40,8 @@ IMPLICIT NONE
END TYPE CPL_FLD
#IFDEF OASIS_IOW_ESM
INTEGER, PARAMETER :: nfld_snd_tot=12
INTEGER, PARAMETER :: nfld_rcv_tot=14
INTEGER, PARAMETER :: nfld_snd_tot=18
INTEGER, PARAMETER :: nfld_rcv_tot=44
#ELSE
INTEGER, PARAMETER :: nfld_snd_tot=2
INTEGER, PARAMETER :: nfld_rcv_tot=14

Loading…
Cancel
Save