test_slicing Subroutine

subroutine test_slicing(rc)

Arguments

Type IntentOptional Attributes Name
integer, intent(out) :: rc

Calls

proc~~test_slicing~~CallsGraph proc~test_slicing test_slicing esmf_fieldcreate esmf_fieldcreate proc~test_slicing->esmf_fieldcreate esmf_fielddestroy esmf_fielddestroy proc~test_slicing->esmf_fielddestroy esmf_fieldget esmf_fieldget proc~test_slicing->esmf_fieldget interface~esmf_gridcreatenoperidim ESMF_GridCreateNoPeriDim proc~test_slicing->interface~esmf_gridcreatenoperidim proc~esmf_griddestroy ESMF_GridDestroy proc~test_slicing->proc~esmf_griddestroy proc~esmf_logfounderror ESMF_LogFoundError proc~test_slicing->proc~esmf_logfounderror proc~esmf_logwrite ESMF_LogWrite proc~test_slicing->proc~esmf_logwrite proc~esmf_gridcreatenoperidima ESMF_GridCreateNoPeriDimA interface~esmf_gridcreatenoperidim->proc~esmf_gridcreatenoperidima proc~esmf_gridcreatenoperidimi ESMF_GridCreateNoPeriDimI interface~esmf_gridcreatenoperidim->proc~esmf_gridcreatenoperidimi proc~esmf_gridcreatenoperidimr ESMF_GridCreateNoPeriDimR interface~esmf_gridcreatenoperidim->proc~esmf_gridcreatenoperidimr proc~esmf_griddestroy->proc~esmf_logfounderror c_esmc_griddestroy c_esmc_griddestroy proc~esmf_griddestroy->c_esmc_griddestroy proc~esmf_gridgetinit ESMF_GridGetInit proc~esmf_griddestroy->proc~esmf_gridgetinit proc~esmf_imerr ESMF_IMErr proc~esmf_griddestroy->proc~esmf_imerr proc~esmf_logfounderror->proc~esmf_logwrite esmf_breakpoint esmf_breakpoint proc~esmf_logfounderror->esmf_breakpoint proc~esmf_logrc2msg ESMF_LogRc2Msg proc~esmf_logfounderror->proc~esmf_logrc2msg c_esmc_vmwtime c_esmc_vmwtime proc~esmf_logwrite->c_esmc_vmwtime proc~esmf_logclose ESMF_LogClose proc~esmf_logwrite->proc~esmf_logclose proc~esmf_logflush ESMF_LogFlush proc~esmf_logwrite->proc~esmf_logflush proc~esmf_logopenfile ESMF_LogOpenFile proc~esmf_logwrite->proc~esmf_logopenfile proc~esmf_utiliounitflush ESMF_UtilIOUnitFlush proc~esmf_logwrite->proc~esmf_utiliounitflush proc~esmf_utilstring2array ESMF_UtilString2Array proc~esmf_logwrite->proc~esmf_utilstring2array

Called by

proc~~test_slicing~~CalledByGraph proc~test_slicing test_slicing program~esmf_fieldcreategetutest ESMF_FieldCreateGetUTest program~esmf_fieldcreategetutest->proc~test_slicing

Source Code

    subroutine test_slicing(rc)
        integer, intent(out)        :: rc
        integer                     :: localrc
        type(ESMF_Field)            :: field1, field2
        type(ESMF_Grid)             :: grid
        character(160)              :: msgStrg
        real(ESMF_KIND_R4), pointer :: fptr1a(:,:,:,:,:), fptr2a(:,:,:)
        real(ESMF_KIND_R4), pointer :: fptr1b(:,:,:), fptr2b(:,:)

        rc = ESMF_SUCCESS
        localrc = ESMF_SUCCESS

        grid = ESMF_GridCreateNoPeriDim(minIndex=(/1,1/), maxIndex=(/16,20/), &
          regDecomp=(/4,1/), name="testgrid", rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        field1 = ESMF_FieldCreate(grid, ESMF_TYPEKIND_R4, name="Field#1", &
          ungriddedLBound=[1,1,1], ungriddedUBound=[10,20,30], rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        field2 = ESMF_FieldCreate(field1, name="Field#2", &
          trailingUngridSlice=[3,7], rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_FieldGet(field1, farrayPtr=fptr1a, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        write(msgStrg,*) "shape(fptr1a)=", shape(fptr1a)
        call ESMF_LogWrite(msgStrg, ESMF_LOGMSG_INFO, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        write(msgStrg,*) "lbound(fptr1a)=", lbound(fptr1a)
        call ESMF_LogWrite(msgStrg, ESMF_LOGMSG_INFO, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        write(msgStrg,*) "ubound(fptr1a)=", ubound(fptr1a)
        call ESMF_LogWrite(msgStrg, ESMF_LOGMSG_INFO, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_FieldGet(field2, farrayPtr=fptr2a, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        write(msgStrg,*) "shape(fptr2a)=", shape(fptr2a)
        call ESMF_LogWrite(msgStrg, ESMF_LOGMSG_INFO, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        write(msgStrg,*) "lbound(fptr2a)=", lbound(fptr2a)
        call ESMF_LogWrite(msgStrg, ESMF_LOGMSG_INFO, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        write(msgStrg,*) "ubound(fptr2a)=", ubound(fptr2a)
        call ESMF_LogWrite(msgStrg, ESMF_LOGMSG_INFO, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_FieldDestroy(field1, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_FieldDestroy(field2, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        field1 = ESMF_FieldCreate(grid, ESMF_TYPEKIND_R4, &
          gridToFieldMap=[0,0], name="Field#1 with replicated dims", &
          ungriddedLBound=[1,1,1], ungriddedUBound=[10,20,30], rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        field2 = ESMF_FieldCreate(field1, name="Field#2", &
          trailingUngridSlice=[4], rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_FieldGet(field1, farrayPtr=fptr1b, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        write(msgStrg,*) "shape(fptr1b)=", shape(fptr1b)
        call ESMF_LogWrite(msgStrg, ESMF_LOGMSG_INFO, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        write(msgStrg,*) "lbound(fptr1b)=", lbound(fptr1b)
        call ESMF_LogWrite(msgStrg, ESMF_LOGMSG_INFO, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        write(msgStrg,*) "ubound(fptr1b)=", ubound(fptr1b)
        call ESMF_LogWrite(msgStrg, ESMF_LOGMSG_INFO, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_FieldGet(field2, farrayPtr=fptr2b, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        write(msgStrg,*) "shape(fptr2b)=", shape(fptr2b)
        call ESMF_LogWrite(msgStrg, ESMF_LOGMSG_INFO, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        write(msgStrg,*) "lbound(fptr2b)=", lbound(fptr2b)
        call ESMF_LogWrite(msgStrg, ESMF_LOGMSG_INFO, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        write(msgStrg,*) "ubound(fptr2b)=", ubound(fptr2b)
        call ESMF_LogWrite(msgStrg, ESMF_LOGMSG_INFO, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_GridDestroy(grid, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

    end subroutine test_slicing