Skip to content
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
22 changes: 18 additions & 4 deletions mediator/med.F90
Original file line number Diff line number Diff line change
Expand Up @@ -37,6 +37,8 @@ module MED
use med_methods_mod , only : FB_diagnose => med_methods_FB_diagnose
use med_methods_mod , only : FB_getFieldN => med_methods_FB_getFieldN
use med_methods_mod , only : clock_timeprint => med_methods_clock_timeprint
use med_field_info_mod , only : med_field_info_type
use med_field_info_mod , only : med_field_info_array_from_names_wtracers_ungridded, med_field_info_array_from_state
use med_utils_mod , only : memcheck => med_memcheck
use med_internalstate_mod , only : InternalState, med_internalstate_init, med_internalstate_coupling
use med_internalstate_mod , only : med_internalstate_defaultmasks, logunit, maintask
Expand Down Expand Up @@ -1636,6 +1638,7 @@ subroutine DataInitialize(gcomp, rc)

! local variables
type(InternalState) :: is_local
type(med_field_info_type), allocatable :: field_info_array(:)
type(ESMF_Clock) :: clock
type(ESMF_State) :: importState, exportState
type(ESMF_Time) :: time
Expand Down Expand Up @@ -1749,19 +1752,25 @@ subroutine DataInitialize(gcomp, rc)
trim(compname(n1))//'_'//trim(compname(n2))
end if

call med_field_info_array_from_state( &
state = is_local%wrap%NStateImp(n1), &
field_info_array = field_info_array, &
rc = rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return

! Check import FB, if there is no field in it then use export FB
! to provide mesh information
call State_GetNumFields(is_local%wrap%NStateImp(n2), fieldCount, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
if (fieldCount == 0) then
call FB_init(is_local%wrap%FBImp(n1,n2), is_local%wrap%flds_scalar_name, &
field_info_array=field_info_array, &
STgeom=is_local%wrap%NStateExp(n2), &
STflds=is_local%wrap%NStateImp(n1), &
name='FBImp'//trim(compname(n1))//'_'//trim(compname(n2)), rc=rc)
else
call FB_init(is_local%wrap%FBImp(n1,n2), is_local%wrap%flds_scalar_name, &
field_info_array=field_info_array, &
STgeom=is_local%wrap%NStateImp(n2), &
STflds=is_local%wrap%NStateImp(n1), &
name='FBImp'//trim(compname(n1))//'_'//trim(compname(n2)), rc=rc)
end if
if (ChkErr(rc,__LINE__,u_FILE_u)) return
Expand Down Expand Up @@ -1789,14 +1798,19 @@ subroutine DataInitialize(gcomp, rc)
allocate(fldnames(fieldCount))
call med_fldList_getfldnames(fldListMed_ocnalb%fields, fldnames, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call med_field_info_array_from_names_wtracers_ungridded( &
field_names = fldnames, &
field_info_array = field_info_array, &
rc = rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call FB_init(is_local%wrap%FBMed_ocnalb_a, is_local%wrap%flds_scalar_name, &
STgeom=is_local%wrap%NStateImp(compatm), fieldnamelist=fldnames, name='FBMed_ocnalb_a', rc=rc)
field_info_array=field_info_array, STgeom=is_local%wrap%NStateImp(compatm), name='FBMed_ocnalb_a', rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
if (maintask) then
write(logunit,'(a)') trim(subname)//' initializing FB FBMed_ocnalb_a'
end if
call FB_init(is_local%wrap%FBMed_ocnalb_o, is_local%wrap%flds_scalar_name, &
STgeom=is_local%wrap%NStateImp(compocn), fieldnamelist=fldnames, name='FBMed_ocnalb_o', rc=rc)
field_info_array = field_info_array, STgeom=is_local%wrap%NStateImp(compocn), name='FBMed_ocnalb_o', rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
if (maintask) then
write(logunit,'(a)') trim(subname)//' initializing FB FBMed_ocnalb_o'
Expand Down
222 changes: 222 additions & 0 deletions mediator/med_field_info_mod.F90
Original file line number Diff line number Diff line change
@@ -0,0 +1,222 @@
module med_field_info_mod

!-----------------------------------------------------------------------------
! Defines a type and related operations for storing metadata about fields that can be
! used to create an ESMF FieldBundle.
!-----------------------------------------------------------------------------

use ESMF, only : ESMF_MAXSTR, ESMF_SUCCESS
use ESMF, only : ESMF_Field, ESMF_State, ESMF_AttributeGet, ESMF_StateGet
use med_utils_mod, only : ChkErr => med_utils_ChkErr
use shr_log_mod, only : shr_log_error
use shr_string_mod, only : shr_string_withoutSuffix
use shr_wtracers_mod, only : WTRACERS_SUFFIX

implicit none
private

!-----------------------------------------------
! Public methods
!-----------------------------------------------

public :: med_field_info_create ! Create a single field
public :: med_field_info_array_from_names_wtracers_ungridded ! Create an array of field_info objects based on an array of names, where water tracers are given an ungridded dimension
public :: med_field_info_array_from_state ! Create an array of field_info objects based on the fields in an ESMF State

!-----------------------------------------------
! Types
!-----------------------------------------------

type, public :: med_field_info_type
character(ESMF_MAXSTR) :: name
integer :: n_ungridded ! number of ungridded dimensions

! These arrays will be allocated to be of size ungridded_count
integer, allocatable :: ungridded_lbound(:)
integer, allocatable :: ungridded_ubound(:)
end type med_field_info_type

character(len=*),parameter :: u_FILE_u = &
__FILE__

!================================================================================
contains
!================================================================================

function med_field_info_create(name, ungridded_lbound, ungridded_ubound, rc) result(field_info)
! Create a single field

! input/output variables
character(len=*), intent(in) :: name

! ungridded_lbound and ungridded_ubound must either both be present or both be absent;
! if present, they must be the same size
integer, intent(in), optional :: ungridded_lbound(:)
integer, intent(in), optional :: ungridded_ubound(:)

integer, intent(out) :: rc
type(med_field_info_type) :: field_info ! function result

! local variables
integer :: n_ungridded
character(len=*), parameter :: subname = '(med_field_info_create)'
! ----------------------------------------------

rc = ESMF_SUCCESS

if (present(ungridded_lbound) .neqv. present(ungridded_ubound)) then
call shr_log_error( &
subname//": ERROR: ungridded_lbound and ungridded_ubound must both be present or both absent.", &
line=__LINE__, file=u_FILE_u, rc=rc)
return
end if

field_info%name = name

if (present(ungridded_lbound)) then
n_ungridded = size(ungridded_lbound)
if (size(ungridded_ubound) /= n_ungridded) then
call shr_log_error( &
subname//": ERROR: ungridded_lbound and ungridded_ubound must have the same size.", &
line=__LINE__, file=u_FILE_u, rc=rc)
return
end if
field_info%n_ungridded = n_ungridded
allocate(field_info%ungridded_lbound(n_ungridded))
allocate(field_info%ungridded_ubound(n_ungridded))
field_info%ungridded_lbound = ungridded_lbound
field_info%ungridded_ubound = ungridded_ubound
else
field_info%n_ungridded = 0
end if

end function med_field_info_create

!-----------------------------------------------------------------------------

subroutine med_field_info_array_from_names_wtracers_ungridded(field_names, field_info_array, rc)
! Create an array of field_info objects based on an array of names, where water
! tracers are given an ungridded dimension.
!
! It is assumed that fields generally have no ungridded dimensions. However, for
! fields ending with the water tracer suffix, it is instead assumed that they have a
! single ungridded dimension of size given by shr_wtracers_get_num_tracers.

! input/output variables
character(len=*), intent(in) :: field_names(:)
type(med_field_info_type), allocatable, intent(out) :: field_info_array(:)
integer, intent(out) :: rc

! local variables
integer :: i, n_fields
logical :: is_tracer
integer :: n_tracers
integer :: localrc
character(len=*), parameter :: subname = '(med_field_info_array_from_names_wtracers_ungridded)'
! ----------------------------------------------

rc = ESMF_SUCCESS

n_fields = size(field_names)
allocate(field_info_array(n_fields))
! For now, hard-code n_tracers, since we haven't set up the tracer information; we'll
! fix this in an upcoming set of changes
n_tracers = 0

do i = 1, n_fields
call shr_string_withoutSuffix( &
in_str = field_names(i), &
suffix = WTRACERS_SUFFIX, &
has_suffix = is_tracer, &
rc = localrc)
if (localrc /= 0) then
call shr_log_error(subname//": ERROR in shr_string_withoutSuffix", rc=rc)
return
end if

if (is_tracer) then
! Field is a water tracer; assume a single ungridded dimension
field_info_array(i) = med_field_info_create( &
name=field_names(i), &
ungridded_lbound=[1], &
ungridded_ubound=[n_tracers], &
rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
else
! Not a water tracer; assume no ungridded dimensions
field_info_array(i) = med_field_info_create( &
name=field_names(i), &
rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
end if
end do

end subroutine med_field_info_array_from_names_wtracers_ungridded

subroutine med_field_info_array_from_state(state, field_info_array, rc)
! Create an array of field_info objects based on the Fields in an ESMF State

! input/output variables
type(ESMF_State), intent(in) :: state
type(med_field_info_type), allocatable, intent(out) :: field_info_array(:)
integer, intent(out) :: rc

! local variables
integer :: i, n_fields
character(ESMF_MAXSTR), allocatable :: field_names(:)
type(ESMF_Field) :: field
logical :: is_present
integer :: n_ungridded
integer, allocatable :: ungridded_lbound(:)
integer, allocatable :: ungridded_ubound(:)
character(len=*), parameter :: subname = '(med_field_info_array_from_state)'
! ----------------------------------------------

rc = ESMF_SUCCESS

call ESMF_StateGet(state, itemCount=n_fields, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
allocate(field_names(n_fields))
allocate(field_info_array(n_fields))
call ESMF_StateGet(state, itemNameList=field_names, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return

do i = 1, n_fields
call ESMF_StateGet(state, itemName=trim(field_names(i)), field=field, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return

call ESMF_AttributeGet(field, name="UngriddedUBound", convention="NUOPC", &
purpose="Instance", itemCount=n_ungridded, isPresent=is_present, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
if (.not. is_present) then
n_ungridded = 0
end if

if (n_ungridded == 0) then
field_info_array(i) = med_field_info_create( &
name=field_names(i), &
rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
else
allocate(ungridded_lbound(n_ungridded))
allocate(ungridded_ubound(n_ungridded))
call ESMF_AttributeGet(field, name="UngriddedLBound", convention="NUOPC", &
purpose="Instance", valueList=ungridded_lbound, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call ESMF_AttributeGet(field, name="UngriddedUBound", convention="NUOPC", &
purpose="Instance", valueList=ungridded_ubound, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
field_info_array(i) = med_field_info_create( &
name=field_names(i), &
ungridded_lbound=ungridded_lbound, &
ungridded_ubound=ungridded_ubound, &
rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
deallocate(ungridded_lbound)
deallocate(ungridded_ubound)
end if
end do

end subroutine med_field_info_array_from_state

end module med_field_info_mod
26 changes: 22 additions & 4 deletions mediator/med_fraction_mod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -135,6 +135,8 @@ module med_fraction_mod
use med_methods_mod , only : fldbun_reset => med_methods_FB_reset
use med_map_mod , only : med_map_field
use med_internalstate_mod , only : ncomps, samegrid_atmlnd
use med_field_info_mod , only : med_field_info_type
use med_field_info_mod , only : med_field_info_array_from_names_wtracers_ungridded, med_field_info_array_from_state

implicit none
private
Expand Down Expand Up @@ -189,6 +191,7 @@ subroutine med_fraction_init(gcomp, rc)
type(InternalState) :: is_local
type(ESMF_Field) :: field_src
type(ESMF_Field) :: field_dst
type(med_field_info_type), allocatable :: field_info_array(:)
real(R8), pointer :: frac(:)
real(R8), pointer :: ofrac(:)
real(R8), pointer :: aofrac(:)
Expand Down Expand Up @@ -255,13 +258,18 @@ subroutine med_fraction_init(gcomp, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return

! create FBFrac
call med_field_info_array_from_names_wtracers_ungridded( &
field_names = fraclist(:,n1), &
field_info_array = field_info_array, &
rc = rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
if (fieldCount == 0) then
call fldbun_init(is_local%wrap%FBfrac(n1), is_local%wrap%flds_scalar_name, &
STgeom=is_local%wrap%NStateExp(n1), fieldNameList=fraclist(:,n1), &
field_info_array=field_info_array, STgeom=is_local%wrap%NStateExp(n1), &
name='FBfrac'//trim(compname(n1)), rc=rc)
else
call fldbun_init(is_local%wrap%FBfrac(n1), is_local%wrap%flds_scalar_name, &
STgeom=is_local%wrap%NStateImp(n1), fieldNameList=fraclist(:,n1), &
field_info_array=field_info_array, STgeom=is_local%wrap%NStateImp(n1), &
name='FBfrac'//trim(compname(n1)), rc=rc)
end if
if (ChkErr(rc,__LINE__,u_FILE_u)) return
Expand Down Expand Up @@ -673,9 +681,14 @@ subroutine med_fraction_init(gcomp, rc)
if (is_local%wrap%comp_present(compice) .and. is_local%wrap%comp_present(compocn)) then
if (.not. med_map_RH_is_created(is_local%wrap%RH(compice,compocn,:),mapfcopy, rc=rc)) then
if (.not. ESMF_FieldBundleIsCreated(is_local%wrap%FBImp(compice,compocn))) then
call med_field_info_array_from_state( &
state = is_local%wrap%NStateImp(compice), &
field_info_array = field_info_array, &
rc = rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call fldbun_init(is_local%wrap%FBImp(compice,compocn), is_local%wrap%flds_scalar_name, &
field_info_array=field_info_array, &
STgeom=is_local%wrap%NStateImp(compocn), &
STflds=is_local%wrap%NStateImp(compice), &
name='FBImp'//trim(compname(compice))//'_'//trim(compname(compocn)), rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
end if
Expand All @@ -687,9 +700,14 @@ subroutine med_fraction_init(gcomp, rc)
end if
if (.not. med_map_RH_is_created(is_local%wrap%RH(compocn,compice,:),mapfcopy, rc=rc)) then
if (.not. ESMF_FieldBundleIsCreated(is_local%wrap%FBImp(compocn,compice))) then
call med_field_info_array_from_state( &
state = is_local%wrap%NStateImp(compocn), &
field_info_array = field_info_array, &
rc = rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call fldbun_init(is_local%wrap%FBImp(compocn,compice), is_local%wrap%flds_scalar_name, &
field_info_array = field_info_array, &
STgeom=is_local%wrap%NStateImp(compice), &
STflds=is_local%wrap%NStateImp(compocn), &
name='FBImp'//trim(compname(compocn))//'_'//trim(compname(compice)), rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
end if
Expand Down
Loading
Loading