Test_FieldArithmetic.pf Source File


This file depends on

sourcefile~~test_fieldarithmetic.pf~~EfferentGraph sourcefile~test_fieldarithmetic.pf Test_FieldArithmetic.pf sourcefile~field_utils_setup.f90 field_utils_setup.F90 sourcefile~test_fieldarithmetic.pf->sourcefile~field_utils_setup.f90 sourcefile~fieldbinaryoperations.f90 FieldBinaryOperations.F90 sourcefile~test_fieldarithmetic.pf->sourcefile~fieldbinaryoperations.f90 sourcefile~fieldpointerutilities.f90 FieldPointerUtilities.F90 sourcefile~test_fieldarithmetic.pf->sourcefile~fieldpointerutilities.f90 sourcefile~fieldunaryfunctions.f90 FieldUnaryFunctions.F90 sourcefile~test_fieldarithmetic.pf->sourcefile~fieldunaryfunctions.f90 sourcefile~fieldutilities.f90 FieldUtilities.F90 sourcefile~test_fieldarithmetic.pf->sourcefile~fieldutilities.f90 sourcefile~mapl_exceptionhandling.f90 MAPL_ExceptionHandling.F90 sourcefile~test_fieldarithmetic.pf->sourcefile~mapl_exceptionhandling.f90 sourcefile~field_utils_setup.f90->sourcefile~mapl_exceptionhandling.f90 sourcefile~fieldbinaryoperations.f90->sourcefile~fieldpointerutilities.f90 sourcefile~fieldbinaryoperations.f90->sourcefile~mapl_exceptionhandling.f90 sourcefile~fieldpointerutilities.f90->sourcefile~mapl_exceptionhandling.f90 sourcefile~fieldunaryfunctions.f90->sourcefile~fieldpointerutilities.f90 sourcefile~fieldunaryfunctions.f90->sourcefile~mapl_exceptionhandling.f90 sourcefile~fieldutilities.f90->sourcefile~fieldpointerutilities.f90 sourcefile~mapl_errorhandling.f90 MAPL_ErrorHandling.F90 sourcefile~fieldutilities.f90->sourcefile~mapl_errorhandling.f90 sourcefile~mapl_exceptionhandling.f90->sourcefile~mapl_errorhandling.f90 sourcefile~mapl_throw.f90 MAPL_Throw.F90 sourcefile~mapl_exceptionhandling.f90->sourcefile~mapl_throw.f90 sourcefile~mapl_errorhandling.f90->sourcefile~mapl_throw.f90

Source Code

#include "MAPL_Generic.h"

module Test_FieldArithmetic

   use field_utils_setup
   use MAPL_FieldUnaryFunctions
   use MAPL_FieldBinaryOperations
   use MAPL_FieldUtilities
   use MAPL_FieldPointerUtilities
   use ESMF
   use pfunit
   use MAPL_ExceptionHandling

   implicit none

   real(kind=ESMF_KIND_R4), parameter :: ADD_R4 = 100.0
   real(kind=ESMF_KIND_R8), parameter :: ADD_R8 = 100.0

contains

   ! Making the fields should be done in the tests themselves so because
   ! of the npes argument.
   @Before
   subroutine set_up_data(this)
      class(MpiTestMethod), intent(inout) :: this

      integer :: status, rc

      real(kind=ESMF_KIND_R4), parameter :: ADD_R4 = 100.0
      real(kind=ESMF_KIND_R8), parameter :: ADD_R8 = 100.0
      real(kind=ESMF_KIND_R4), allocatable :: y4array(:,:)
      real(kind=ESMF_KIND_R8), allocatable :: y8array(:,:)

      allocate(y4array, source=R4_ARRAY_DEFAULT)
      allocate(y8array, source=R8_ARRAY_DEFAULT)
      y4array = y4array + ADD_R4
      y8array = y8array + ADD_R8
      XR4 = mk_field(R4_ARRAY_DEFAULT, name = 'XR4', _RC)
      YR4 = mk_field(y4array, name = 'YR4', _RC)
      XR8 = mk_field(R8_ARRAY_DEFAULT, name = 'XR8', _RC)
      YR8 = mk_field(y8array, name = 'YR8', _RC)
      call ESMF_AttributeSet(xr4,name="missing_value",value=undef,_RC)
      call ESMF_AttributeSet(xr8,name="missing_value",value=undef,_RC)
      call ESMF_AttributeSet(yr4,name="missing_value",value=undef,_RC)
      call ESMF_AttributeSet(yr8,name="missing_value",value=undef,_RC)

   end subroutine set_up_data

   @after
   subroutine teardown(this)
      class(MpiTestMethod), intent(inout) :: this
   end subroutine teardown

   @Test(npes=[4])
   subroutine test_FieldAddR4(this)
      class(MpiTestMethod), intent(inout) :: this
      type(ESMF_Field) :: x
      type(ESMF_Field) :: y
      real(kind=ESMF_KIND_R4), pointer :: x_ptr(:,:), y_ptr(:,:)
      real(kind=ESMF_KIND_R4), allocatable :: result_array(:,:)
      integer :: status, rc
      real(kind=ESMF_KIND_R4), allocatable :: y4array(:,:)

      allocate(y4array, source=R4_ARRAY_DEFAULT)
      x = mk_r4field(R4_ARRAY_DEFAULT, 'XR4', _RC)
      y = mk_r4field(y4array, 'YR4', _RC)
      call ESMF_FieldGet(x , farrayPtr = x_ptr, _RC)
      call ESMF_FieldGet(y , farrayPtr = y_ptr, _RC)

      x_ptr = 2.0
      y_ptr = 3.0
      result_array = x_ptr
      result_array = 5.0
      call FieldAdd(y, x, y, _RC)
      @assertEqual(y_ptr, result_array)

   end subroutine test_FieldAddR4

   @Test(npes=[4])
   subroutine test_FieldAddR4_missing(this)
      class(MpiTestMethod), intent(inout) :: this
      type(ESMF_Field) :: x
      type(ESMF_Field) :: y
      real(kind=ESMF_KIND_R4), pointer :: x_ptr(:,:), y_ptr(:,:)
      real(kind=ESMF_KIND_R4), allocatable :: result_array(:,:)
      integer :: status, rc

      x = XR4
      y = YR4
      call ESMF_FieldGet(x , farrayPtr = x_ptr, _RC)
      call ESMF_FieldGet(y , farrayPtr = y_ptr, _RC)

      x_ptr = reshape(source=[2.0,2.0,2.0,undef],shape=[2,2])
      y_ptr = reshape(source=[undef,3.0,3.0,undef],shape=[2,2])
      result_array = x_ptr
      result_array = reshape(source=[undef,5.0,5.0,undef],shape=[2,2])
      call FieldAdd(y, x, y, _RC)
      @assertEqual(y_ptr, result_array)
   end subroutine test_FieldAddR4_missing

   @Test(npes=[4])
   subroutine test_FieldAddR8(this)
      class(MpiTestMethod), intent(inout) :: this
      type(ESMF_Field) :: x
      type(ESMF_Field) :: y
      real(kind=ESMF_KIND_R8), pointer :: x_ptr(:,:), y_ptr(:,:)
      real(kind=ESMF_KIND_R8), allocatable :: result_array(:,:)
      integer :: status, rc
      real(kind=ESMF_KIND_R8), allocatable :: y8array(:,:)

      allocate(y8array, source=R8_ARRAY_DEFAULT)
      x = mk_r8field(R8_ARRAY_DEFAULT, 'XR8', _RC)
      y = mk_r8field(y8array, 'YR8', _RC)
      call ESMF_FieldGet(x , farrayPtr = x_ptr, _RC)
      call ESMF_FieldGet(y , farrayPtr = y_ptr, _RC)

      x_ptr = 2.d0
      y_ptr = 3.d0
      result_array = x_ptr
      result_array = 5.d0
      call FieldAdd(y, x, y, _RC)
      @assertEqual(y_ptr, result_array)
   end subroutine test_FieldAddR8

   @Test(npes=[4])
   subroutine test_FieldPowR4(this)
      class(MpiTestMethod), intent(inout) :: this
      type(ESMF_Field) :: x
      real(kind=ESMF_KIND_R4), pointer :: x_ptr(:,:)
      real(kind=ESMF_KIND_R4), allocatable :: result_array(:,:)
      integer :: status, rc
      real :: expo

      x = XR4
      call ESMF_FieldGet(x , farrayPtr = x_ptr, _RC)

      x_ptr = 2.0
      expo = 4.0
      result_array = x_ptr
      result_array = 2.0**expo
      call FieldPow(x, x, expo, _RC)
      @assertEqual(x_ptr, result_array)
   end subroutine test_FieldPowR4

   @Test(npes=[4])
   subroutine test_FieldPowR8(this)
      class(MpiTestMethod), intent(inout) :: this
      type(ESMF_Field) :: x
      real(kind=ESMF_KIND_R8), pointer :: x_ptr(:,:)
      real(kind=ESMF_KIND_R8), allocatable :: result_array(:,:)
      integer :: status, rc
      real :: expo

      x = XR8
      call ESMF_FieldGet(x , farrayPtr = x_ptr, _RC)

      x_ptr = 2.d0
      expo = 4.0
      result_array = x_ptr
      result_array = 2.d0**expo
      call FieldPow(x, x, expo, _RC)
      @assertEqual(x_ptr, result_array)
   end subroutine test_FieldPowR8

   @Test(npes=[4])
   subroutine test_FieldSinR4(this)
      class(MpiTestMethod), intent(inout) :: this
      type(ESMF_Field) :: x
      real(kind=ESMF_KIND_R4), pointer :: x_ptr(:,:)
      real(kind=ESMF_KIND_R4), allocatable :: result_array(:,:)
      integer :: status, rc

      x = XR4
      call ESMF_FieldGet(x , farrayPtr = x_ptr, _RC)

      x_ptr = 2.0
      result_array = x_ptr
      result_array = sin(2.0)
      call FieldSin(x, x, _RC)
      @assertEqual(x_ptr, result_array)
   end subroutine test_FieldSinR4

   @Test(npes=[4])
   subroutine test_FieldNegR4(this)
      class(MpiTestMethod), intent(inout) :: this
      type(ESMF_Field) :: x
      real(kind=ESMF_KIND_R4), pointer :: x_ptr(:,:)
      real(kind=ESMF_KIND_R4), allocatable :: result_array(:,:)
      integer :: status, rc

      x = XR4
      call ESMF_FieldGet(x , farrayPtr = x_ptr, _RC)

      x_ptr = 2.0
      result_array = x_ptr
      result_array = -2.0
      call FieldNegate(x, _RC)
      @assertEqual(x_ptr, result_array)
   end subroutine test_FieldNegR4

end module Test_FieldArithmetic