FieldCondensedArray.F90 Source File


This file depends on

sourcefile~~fieldcondensedarray.f90~~EfferentGraph sourcefile~fieldcondensedarray.f90 FieldCondensedArray.F90 sourcefile~fieldcondensedarray_private.f90 FieldCondensedArray_private.F90 sourcefile~fieldcondensedarray.f90->sourcefile~fieldcondensedarray_private.f90 sourcefile~fieldget.f90 FieldGet.F90 sourcefile~fieldcondensedarray.f90->sourcefile~fieldget.f90 sourcefile~fieldpointerutilities.f90 FieldPointerUtilities.F90 sourcefile~fieldcondensedarray.f90->sourcefile~fieldpointerutilities.f90 sourcefile~mapl_exceptionhandling.f90 MAPL_ExceptionHandling.F90 sourcefile~fieldcondensedarray.f90->sourcefile~mapl_exceptionhandling.f90 sourcefile~verticalstaggerloc.f90 VerticalStaggerLoc.F90 sourcefile~fieldcondensedarray.f90->sourcefile~verticalstaggerloc.f90 sourcefile~fieldcondensedarray_private.f90->sourcefile~mapl_exceptionhandling.f90 sourcefile~fieldget.f90->sourcefile~verticalstaggerloc.f90 sourcefile~errorhandling.f90 ErrorHandling.F90 sourcefile~fieldget.f90->sourcefile~errorhandling.f90 sourcefile~fieldinfo.f90 FieldInfo.F90 sourcefile~fieldget.f90->sourcefile~fieldinfo.f90 sourcefile~keywordenforcer.f90 KeywordEnforcer.F90 sourcefile~fieldget.f90->sourcefile~keywordenforcer.f90 sourcefile~ungriddeddims.f90 UngriddedDims.F90 sourcefile~fieldget.f90->sourcefile~ungriddeddims.f90 sourcefile~fieldpointerutilities.f90->sourcefile~mapl_exceptionhandling.f90 sourcefile~mapl_exceptionhandling.f90->sourcefile~errorhandling.f90 sourcefile~mapl_throw.f90 MAPL_Throw.F90 sourcefile~mapl_exceptionhandling.f90->sourcefile~mapl_throw.f90 sourcefile~errorhandling.f90->sourcefile~mapl_throw.f90 sourcefile~fieldinfo.f90->sourcefile~verticalstaggerloc.f90 sourcefile~fieldinfo.f90->sourcefile~errorhandling.f90 sourcefile~fieldinfo.f90->sourcefile~keywordenforcer.f90 sourcefile~fieldinfo.f90->sourcefile~ungriddeddims.f90 sourcefile~infoutilities.f90 InfoUtilities.F90 sourcefile~fieldinfo.f90->sourcefile~infoutilities.f90 sourcefile~mapl_esmf_infokeys.f90 MAPL_ESMF_InfoKeys.F90 sourcefile~fieldinfo.f90->sourcefile~mapl_esmf_infokeys.f90 sourcefile~ungriddeddims.f90->sourcefile~errorhandling.f90 sourcefile~ungriddeddims.f90->sourcefile~infoutilities.f90 sourcefile~lu_bound.f90 LU_Bound.F90 sourcefile~ungriddeddims.f90->sourcefile~lu_bound.f90 sourcefile~ungriddeddims.f90->sourcefile~mapl_esmf_infokeys.f90 sourcefile~ungriddeddim.f90 UngriddedDim.F90 sourcefile~ungriddeddims.f90->sourcefile~ungriddeddim.f90 sourcefile~ungriddeddimvector.f90 UngriddedDimVector.F90 sourcefile~ungriddeddims.f90->sourcefile~ungriddeddimvector.f90 sourcefile~infoutilities.f90->sourcefile~errorhandling.f90 sourcefile~infoutilities.f90->sourcefile~keywordenforcer.f90 sourcefile~infoutilities.f90->sourcefile~mapl_esmf_infokeys.f90 sourcefile~mapl_esmf_infokeys.f90->sourcefile~errorhandling.f90 sourcefile~ungriddeddim.f90->sourcefile~errorhandling.f90 sourcefile~ungriddeddim.f90->sourcefile~infoutilities.f90 sourcefile~ungriddeddim.f90->sourcefile~lu_bound.f90 sourcefile~ungriddeddimvector.f90->sourcefile~ungriddeddim.f90

Files dependent on this one

sourcefile~~fieldcondensedarray.f90~~AfferentGraph sourcefile~fieldcondensedarray.f90 FieldCondensedArray.F90 sourcefile~fixedlevelsverticalgrid.f90 FixedLevelsVerticalGrid.F90 sourcefile~fixedlevelsverticalgrid.f90->sourcefile~fieldcondensedarray.f90 sourcefile~verticalregridaction.f90 VerticalRegridAction.F90 sourcefile~verticalregridaction.f90->sourcefile~fieldcondensedarray.f90 sourcefile~can_connect_to.f90 can_connect_to.F90 sourcefile~can_connect_to.f90->sourcefile~fixedlevelsverticalgrid.f90 sourcefile~modelverticalgrid.f90 ModelVerticalGrid.F90 sourcefile~can_connect_to.f90->sourcefile~modelverticalgrid.f90 sourcefile~couplermetacomponent.f90 CouplerMetaComponent.F90 sourcefile~couplermetacomponent.f90->sourcefile~verticalregridaction.f90 sourcefile~fieldspec.f90 FieldSpec.F90 sourcefile~fieldspec.f90->sourcefile~fixedlevelsverticalgrid.f90 sourcefile~fieldspec.f90->sourcefile~verticalregridaction.f90 sourcefile~genericcoupler.f90 GenericCoupler.F90 sourcefile~genericcoupler.f90->sourcefile~verticalregridaction.f90 sourcefile~genericcoupler.f90->sourcefile~couplermetacomponent.f90 sourcefile~parse_geometry_spec.f90 parse_geometry_spec.F90 sourcefile~parse_geometry_spec.f90->sourcefile~fixedlevelsverticalgrid.f90 sourcefile~parse_geometry_spec.f90->sourcefile~modelverticalgrid.f90 sourcefile~test_fixedlevelsverticalgrid.pf Test_FixedLevelsVerticalGrid.pf sourcefile~test_fixedlevelsverticalgrid.pf->sourcefile~fixedlevelsverticalgrid.f90 sourcefile~bracketspec.f90 BracketSpec.F90 sourcefile~bracketspec.f90->sourcefile~fieldspec.f90 sourcefile~make_itemspec.f90 make_itemSpec.F90 sourcefile~make_itemspec.f90->sourcefile~fieldspec.f90 sourcefile~make_itemspec.f90->sourcefile~bracketspec.f90 sourcefile~modelverticalgrid.f90->sourcefile~fieldspec.f90 sourcefile~stateitemextension.f90 StateItemExtension.F90 sourcefile~modelverticalgrid.f90->sourcefile~stateitemextension.f90 sourcefile~stateitemextension.f90->sourcefile~genericcoupler.f90 sourcefile~test_addfieldspec.pf Test_AddFieldSpec.pf sourcefile~test_addfieldspec.pf->sourcefile~fieldspec.f90 sourcefile~test_bracketspec.pf Test_BracketSpec.pf sourcefile~test_bracketspec.pf->sourcefile~fieldspec.f90 sourcefile~test_bracketspec.pf->sourcefile~bracketspec.f90 sourcefile~test_fieldspec.pf Test_FieldSpec.pf sourcefile~test_fieldspec.pf->sourcefile~fieldspec.f90 sourcefile~can_connect_to.f90~2 can_connect_to.F90 sourcefile~can_connect_to.f90~2->sourcefile~modelverticalgrid.f90 sourcefile~can_connect_to.f90~3 can_connect_to.F90 sourcefile~can_connect_to.f90~3->sourcefile~modelverticalgrid.f90 sourcefile~extensionfamily.f90 ExtensionFamily.F90 sourcefile~extensionfamily.f90->sourcefile~stateitemextension.f90 sourcefile~initialize_advertise.f90 initialize_advertise.F90 sourcefile~initialize_advertise.f90->sourcefile~make_itemspec.f90 sourcefile~mapl_generic.f90 MAPL_Generic.F90 sourcefile~mapl_generic.f90->sourcefile~stateitemextension.f90 sourcefile~matchconnection.f90 MatchConnection.F90 sourcefile~matchconnection.f90->sourcefile~stateitemextension.f90 sourcefile~protoextdatagc.f90 ProtoExtDataGC.F90 sourcefile~protoextdatagc.f90->sourcefile~stateitemextension.f90 sourcefile~servicespec.f90 ServiceSpec.F90 sourcefile~servicespec.f90->sourcefile~stateitemextension.f90 sourcefile~simpleconnection.f90 SimpleConnection.F90 sourcefile~simpleconnection.f90->sourcefile~stateitemextension.f90 sourcefile~stateitemextensionptrvector.f90 StateItemExtensionPtrVector.F90 sourcefile~stateitemextensionptrvector.f90->sourcefile~stateitemextension.f90 sourcefile~stateitemextensionvector.f90 StateItemExtensionVector.F90 sourcefile~stateitemextensionvector.f90->sourcefile~stateitemextension.f90 sourcefile~stateregistry.f90 StateRegistry.F90 sourcefile~stateregistry.f90->sourcefile~stateitemextension.f90 sourcefile~test_extensionfamily.pf Test_ExtensionFamily.pf sourcefile~test_extensionfamily.pf->sourcefile~stateitemextension.f90 sourcefile~test_modelverticalgrid.pf Test_ModelVerticalGrid.pf sourcefile~test_modelverticalgrid.pf->sourcefile~make_itemspec.f90 sourcefile~test_modelverticalgrid.pf->sourcefile~modelverticalgrid.f90 sourcefile~test_modelverticalgrid.pf->sourcefile~stateitemextension.f90 sourcefile~test_stateregistry.pf Test_StateRegistry.pf sourcefile~test_stateregistry.pf->sourcefile~stateitemextension.f90

Source Code

#include "MAPL_Generic.h"
module mapl3g_FieldCondensedArray
   use mapl3g_FieldCondensedArray_private, only: ARRAY_RANK, get_fptr_shape_private
   use mapl_FieldPointerUtilities, only: FieldGetLocalElementCount, assign_fptr
   use mapl3g_VerticalStaggerLoc
   use mapl_ExceptionHandling
   use mapl3g_FieldGet
   use ESMF, only: ESMF_Field, ESMF_FieldGet
   use ESMF, only: ESMF_KIND_R4, ESMF_KIND_R8, ESMF_KIND_I8
   use, intrinsic :: iso_c_binding, only: c_ptr, c_f_pointer
   implicit none(type, external)
   private

   public :: assign_fptr_condensed_array

   interface assign_fptr_condensed_array
      module procedure :: assign_fptr_condensed_array_r4
      module procedure :: assign_fptr_condensed_array_r8
   end interface assign_fptr_condensed_array

contains

   subroutine assign_fptr_condensed_array_r4(x, fptr, rc)
      type(ESMF_Field), intent(inout) :: x
      real(kind=ESMF_KIND_R4), pointer, intent(out) :: fptr(:,:,:)
      integer, optional, intent(out) :: rc
      integer(ESMF_KIND_I8) :: fp_shape(ARRAY_RANK)
      integer :: status

      fp_shape = get_fptr_shape(x, _RC)
      call assign_fptr(x, fp_shape, fptr, _RC)
      _RETURN(_SUCCESS)

   end subroutine assign_fptr_condensed_array_r4

   subroutine assign_fptr_condensed_array_r8(x, fptr, rc)
      type(ESMF_Field), intent(inout) :: x
      real(kind=ESMF_KIND_R8), pointer, intent(out) :: fptr(:,:,:)
      integer, optional, intent(out) :: rc
      integer(ESMF_KIND_I8) :: fp_shape(ARRAY_RANK)
      integer :: status

      fp_shape = get_fptr_shape(x, _RC)
      call assign_fptr(x, fp_shape, fptr, _RC)
      _RETURN(_SUCCESS)

   end subroutine assign_fptr_condensed_array_r8

   function get_fptr_shape(f, rc) result(fptr_shape)
      integer :: fptr_shape(ARRAY_RANK)
      type(ESMF_Field), intent(inout) :: f
      integer, optional, intent(out) :: rc
      integer :: status
      integer :: rank
      integer, allocatable :: gridToFieldMap(:)
      integer, allocatable :: localElementCount(:)
      logical :: has_vertical
      integer :: geomDimCount
      type(VerticalStaggerLoc) :: vert_staggerloc

      call ESMF_FieldGet(f, geomDimCount=geomDimCount, rank=rank, _RC)
      _ASSERT(.not. rank < 0, 'rank cannot be negative.')
      _ASSERT(.not. geomDimCount < 0, 'geomDimCount cannot be negative.')
      allocate(localElementCount(rank))
      allocate(gridToFieldMap(geomDimCount))
      call ESMF_FieldGet(f, gridToFieldMap=gridToFieldMap, _RC) 
      !  Due to an ESMF bug, getting the localElementCount must use the module function.
      !  See FieldGetLocalElementCount (specific function) comments.
      localElementCount = FieldGetLocalElementCount(f, _RC)
      call MAPL_FieldGet(f, vert_staggerloc=vert_staggerloc, _RC)
      has_vertical = (vert_staggerloc /= VERTICAL_STAGGER_NONE)
      fptr_shape = get_fptr_shape_private(gridToFieldMap, localElementCount, has_vertical, _RC)

      _RETURN(_SUCCESS)
   end function get_fptr_shape

end module mapl3g_FieldCondensedArray