populate_field Subroutine

private subroutine populate_field(Field, Grid, Memory, Grid_info, TestFunction, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_Field), intent(inout) :: Field
type(ESMF_Grid), intent(in) :: Grid
type(memory_config), intent(in) :: Memory
type(grid_specification_record), intent(in) :: Grid_info
type(test_function_record), intent(in) :: TestFunction
integer, intent(inout) :: rc

Source Code

  subroutine populate_field(Field, Grid, Memory, Grid_info, TestFunction, rc)
  !-----------------------------------------------------------------------------
  ! routine populates the source field with the test function
  !-----------------------------------------------------------------------------
  ! arguments
  type(ESMF_Field), intent(inout) :: Field
  type(ESMF_Grid), intent(in   ) :: Grid
  type(memory_config), intent(in   ) :: Memory
  type(grid_specification_record), intent(in   ) :: Grid_info
  type(test_function_record), intent(in   ) :: TestFunction
  integer, intent(inout) :: rc

  ! local ESMF types

  ! local integer variables
  integer :: i,j,k
  ! integer :: l,m
  integer, allocatable :: lbnd(:), ubnd(:)
  integer :: localDECount, lDE, GridRank

  integer :: localrc ! local error status

  ! local real variables
  real(ESMF_KIND_R8) :: a, b, kx, ly, lenk, lenl
  ! real(ESMF_KIND_R8) :: exact
  real(ESMF_KIND_R8), pointer :: coordX2D(:,:), coordY2D(:,:)
  real(ESMF_KIND_R8), pointer :: exp2D(:,:)
  real(ESMF_KIND_R8), pointer :: coordX3D(:,:,:), coordY3D(:,:,:)
  real(ESMF_KIND_R8), pointer :: coordZ3D(:,:,:)
  real(ESMF_KIND_R8), pointer :: exp3D(:,:,:)

  ! initialize return flag
  localrc = ESMF_RC_NOT_IMPL
  rc = ESMF_RC_NOT_IMPL

  !-----------------------------------------------------------------------------
  ! Get the number of local DEs from the Grid
  !-----------------------------------------------------------------------------
  call ESMF_GridGet(grid=Grid, localDECount=localDECount, rc=localrc)
  if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting local DE count from grid",  &
                            rcToReturn=rc)) return

  !-----------------------------------------------------------------------------
  ! Get the grid rank
  !-----------------------------------------------------------------------------
  call ESMF_GridGet(grid=Grid, dimCount=GridRank, rc=localrc)
  if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting Grid rank from grid",       &
                            rcToReturn=rc)) return

  !-----------------------------------------------------------------------------
  ! allocate arrays to hold the local array bounds
  !-----------------------------------------------------------------------------
  allocate( lbnd(GridRank) )
  allocate( ubnd(GridRank) )


  !-----------------------------------------------------------------------------
  ! populate Field
  !-----------------------------------------------------------------------------

  do lDE=0,localDECount-1
    select case(GridRank)
      case(1)
      !-------------------------------------------------------------------------
      ! grid rank = 1
      !-------------------------------------------------------------------------
      localrc = ESMF_FAILURE
      call ESMF_LogSetError(ESMF_FAILURE, msg="Grid rank=1 not supported ",      &
                    rcToReturn=localrc)
      return

      case(2)
      !-------------------------------------------------------------------------
      ! grid rank = 2
      !-------------------------------------------------------------------------
      call ESMF_GridGetCoord(Grid, localDE=lDE,                                &
                   staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=1,              &
                   computationalLBound=lbnd, computationalUBound=ubnd,         &
                   farrayPtr=coordX2D, rc=localrc)
      if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting grid=1 coordinates",    &
                            rcToReturn=rc)) return

      call ESMF_GridGetCoord(Grid, localDE=lDE,                                &
                   staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=2,              &
                   computationalLBound=lbnd, computationalUBound=ubnd,         &
                   farrayPtr=coordY2D, rc=localrc)
      if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting grid=2 coordinates",    &
                            rcToReturn=rc)) return

      !-------------------------------------------------------------------------
      ! Get pointers to field array
      !-------------------------------------------------------------------------
      call ESMF_FieldGet(Field, lDE, exp2D, rc=localrc)
      if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error field get", rcToReturn=rc)) return

      !-------------------------------------------------------------------------
      !  set coordinates
      !-------------------------------------------------------------------------
      do j=lbnd(2),ubnd(2)
         do i=lbnd(1),ubnd(1)
           !-------------------------------------------------------------------
           select case( trim(TestFunction%string) )
             case("CONSTANT")
               exp2D(i,j) = TestFunction%param(1)

             case("COORDINATEX")
               exp2D(i,j) =  TestFunction%param(1)*coordX2D(i,j)

             case("COORDINATEY")
               exp2D(i,j) =  TestFunction%param(1)*coordY2D(i,j)

             case("SPHERICAL_HARMONIC")
              ! (a+b) + acos(k*2pi*x/Lx) + b*sin(l*2pi*y/Ly)
              a = TestFunction%param(1)
              kx = TestFunction%param(2)
              b = TestFunction%param(3)
              lY = TestFunction%param(4)
              lenk = Grid_info%grange(1,2) - Grid_info%grange(1,1)
              lenl = Grid_info%grange(2,2) - Grid_info%grange(2,1)
               exp2D(i,j) =  abs(a+b) +a*cos(pi2*kx*coordX2D(i,j)/lenk) +        &
                             b*sin(pi2*ly*coordY2D(i,j)/lenl)

             !case("PEAK_VALLEY")
             ! (1.-x*y)*sin(k*pi*x/Lx)*cos(l*pi*y)+2.

             case default
               ! error
                call ESMF_LogSetError( ESMF_FAILURE, msg="undefined case in select statement",    &
                       rcToReturn=localrc)
                return
           end select
           !-------------------------------------------------------------------
         enddo    ! i loop
      enddo    ! j loop

      case(3)
      !-------------------------------------------------------------------------
      ! grid rank = 3
      !-------------------------------------------------------------------------
      call ESMF_GridGetCoord(Grid, localDE=lDE,                                &
                   staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=1,              &
                   computationalLBound=lbnd, computationalUBound=ubnd,         &
                   farrayPtr=coordX3D, rc=localrc)
      if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting grid=1 coordinates",    &
                            rcToReturn=rc)) return

      call ESMF_GridGetCoord(Grid, localDE=lDE,                                &
                   staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=2,              &
                   computationalLBound=lbnd, computationalUBound=ubnd,         &
                   farrayPtr=coordY3D, rc=localrc)
      if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting grid=2 coordinates",    &
                            rcToReturn=rc)) return

      call ESMF_GridGetCoord(Grid, localDE=lDE,                                &
                   staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=3,              &
                   computationalLBound=lbnd, computationalUBound=ubnd,         &
                   farrayPtr=coordZ3D, rc=localrc)
      if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting grid=3 coordinates",    &
                            rcToReturn=rc)) return

      !-------------------------------------------------------------------------
      ! Get pointers to field array
      !-------------------------------------------------------------------------
     call ESMF_FieldGet(Field, lDE, exp3D, rc=localrc)
      if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error field get", rcToReturn=rc)) return

      !-------------------------------------------------------------------------
      !  set coordinates
      !-------------------------------------------------------------------------
      do k=lbnd(3),ubnd(3)
         do j=lbnd(2),ubnd(2)
            do i=lbnd(1),ubnd(1)
           !-------------------------------------------------------------------
           select case( trim(TestFunction%string) )
             case("CONSTANT")
               exp3D(i,j,k) = TestFunction%param(1)

             case("COORDINATEX")
               exp3D(i,j,k) =  TestFunction%param(1)*coordX3D(i,j,k)

             case("COORDINATEY")
               exp3D(i,j,k) =  TestFunction%param(1)*coordY3D(i,j,k)

             case("COORDINATEZ")
               exp3D(i,j,k) =  TestFunction%param(1)*coordZ3D(i,j,k)

             case default
               ! error
                call ESMF_LogSetError( ESMF_FAILURE, msg="undefined case in select statement",    &
                       rcToReturn=localrc)
                return
           end select
           !-------------------------------------------------------------------
            ! coordX3D(i,j,k)
            ! coordY3D(i,j,k)
            ! coordZ3D(i,j,k)
            enddo    ! i loop
         enddo    ! j loop
      enddo    ! k loop

      case(4)
      !-------------------------------------------------------------------------
      ! grid rank = 4
      !-------------------------------------------------------------------------
          localrc = ESMF_FAILURE
           call ESMF_LogSetError(ESMF_FAILURE, msg="Grid rank=4 not supported ", &
                    rcToReturn=localrc)
           return

      case(5)
      !-------------------------------------------------------------------------
      ! grid rank = 5
      !-------------------------------------------------------------------------
          localrc = ESMF_FAILURE
           call ESMF_LogSetError(ESMF_FAILURE, msg="Grid rank=5 not supported ", &
                    rcToReturn=localrc)
           return

      case(6)
      !-------------------------------------------------------------------------
      ! grid rank = 6
      !-------------------------------------------------------------------------
          localrc = ESMF_FAILURE
           call ESMF_LogSetError(ESMF_FAILURE, msg="Grid rank=6 not supported ", &
                    rcToReturn=localrc)
           return

      case(7)
      !-------------------------------------------------------------------------
      ! grid rank = 7
      !-------------------------------------------------------------------------
           localrc = ESMF_FAILURE
           call ESMF_LogSetError(ESMF_FAILURE, msg="Grid rank=7 not supported ", &
                    rcToReturn=localrc)
           return
      case default
      !-------------------------------------------------------------------------
      ! error
      !-------------------------------------------------------------------------
           localrc = ESMF_FAILURE
           call ESMF_LogSetError(ESMF_FAILURE, msg="Grid rank not between 1 & 7",&
                    rcToReturn=localrc)
           return

      end select

   enddo    ! lDE
  !-----------------------------------------------------------------------------
  !-----------------------------------------------------------------------------
  ! clean up
  !-----------------------------------------------------------------------------

  deallocate( lbnd, ubnd)

  !-----------------------------------------------------------------------------
  rc = ESMF_SUCCESS
  !-----------------------------------------------------------------------------

  !-----------------------------------------------------------------------------
  end subroutine populate_field