! $Id: ESMF_FieldRedistArbUTest.F90,v 1.16 2011/07/02 05:54:03 oehmke Exp $ ! ! Earth System Modeling Framework ! Copyright 2002-2011, University Corporation for Atmospheric Research, ! Massachusetts Institute of Technology, Geophysical Fluid Dynamics ! Laboratory, University of Michigan, National Centers for Environmental ! Prediction, Los Alamos National Laboratory, Argonne National Laboratory, ! NASA Goddard Space Flight Center. ! Licensed under the University of Illinois-NCSA License. ! !============================================================================== ! program ESMF_FieldRedistArbUTest !------------------------------------------------------------------------------ #include "ESMF.h" #define ESMF_METHOD "ESMF_TESTS" !============================================================================== !BOPI ! !PROGRAM: ESMF_FieldRedistArbUTest - Unit tests for FieldRedist methods on a arbitrarily distributed grid ! ! !DESCRIPTION: ! ! The code in this file drives F90 Field Redist on Arbitrary Grid unit tests . ! The companion folder Fieldsrc contains the definitions for the ! Field methods. !EOPI !----------------------------------------------------------------------------- ! !USES: use ESMF_TestMod ! test methods use ESMF use ESMF_GridUtilMod use ESMF_LogErrMod implicit none !------------------------------------------------------------------------------ ! The following line turns the CVS identifier string into a printable variable. character(*), parameter :: version = & '$Id' type(ESMF_Grid) :: srcgrid, dstgrid, srcgrid2D, dstgrid2D type(ESMF_VM) :: vm type(ESMF_ArraySpec) :: arrayspec1D, arrayspec2D integer :: ind1d, xdim, ydim, zdim, total, x, y integer :: i, remain integer :: myPet, petCount, halfPets integer :: localCount, localCount1 integer, allocatable :: localIndices(:,:), localIndices1(:,:) integer :: localrc type(ESMF_Field) :: srcfield, dstfield, srcfield2D, dstfield2D type(ESMF_RouteHandle) :: rhandle, rhandle2D real, dimension(:,:), pointer :: fptr real, dimension(:), pointer :: fptr1D integer :: lbnd(2), ubnd(2), cnt(2) integer :: lbnd1(1), ubnd1(1), cnt1(1) logical :: correct ! cumulative result: count failures; no failures equals "all pass" integer :: result = 0 ! individual test result code integer :: rc = 1 ! individual test failure message character(ESMF_MAXSTR) :: failMsg character(512) :: name call ESMF_TestStart(ESMF_SRCLINE, rc=rc) ! Calculate localIndices and localCount for a 100x200 2D arbitrary grid with ! an optional undistributed 3rd dimenison of size 4 ! get global VM call ESMF_VMGetGlobal(vm, rc=rc) if (rc /= ESMF_SUCCESS) call ESMF_Finalize(endflag=ESMF_END_ABORT) call ESMF_VMGet(vm, petCount=petCount, localPet=myPet, rc=rc) if (rc /= ESMF_SUCCESS) call ESMF_Finalize(endflag=ESMF_END_ABORT) ! grid dimension: xdim and ydim are arbitrarily distributed xdim = 100 ydim = 200 zdim = 4 ! -------------------------------------------------------- ! Set up distribution for the source grid ! -------------------------------------------------------- ! calculate the localcount and the local indices based on the total number of PETS total = xdim*ydim halfPets = petCount/2 ! let's make the first half pet twice of the cells of the second half localCount = total/(petCount+halfPets) remain = total-localCount*(petCount+halfPets) if (myPet < halfPets) localCount = localCount*2 if (myPet == petCount-1) localCount = localCount+remain ! car deal the cells with the first half of the Pets gets two each time ! the remaining cells are given to the last Pet allocate(localIndices(localCount,2)) if (myPet < halfPets) then ind1d = myPet*2 do i=1,localCount,2 y = mod(ind1d,ydim)+1 x = ind1d/ydim+1 localIndices(i,1)=y localIndices(i,2)=x if (y 0) correct = .false. call ESMF_Test(((rc.eq.ESMF_SUCCESS) .and. correct), name, failMsg, result, ESMF_SRCLINE) #endif call ESMF_FieldDestroy(srcfield, rc=localrc) call ESMF_FieldDestroy(dstfield, rc=localrc) call ESMF_FieldDestroy(srcfield2D, rc=localrc) call ESMF_FieldDestroy(dstfield2D, rc=localrc) call ESMF_GridDestroy(srcgrid, rc=localrc) call ESMF_GridDestroy(dstgrid, rc=localrc) call ESMF_GridDestroy(srcgrid2D, rc=localrc) call ESMF_GridDestroy(dstgrid2D, rc=localrc) call ESMF_FieldRedistRelease(rhandle, rc=localrc) call ESMF_FieldRedistRelease(rhandle2D, rc=localrc) deallocate(localIndices, localIndices1) !----------------------------------------------------------------------------- call ESMF_TestEnd(result, ESMF_SRCLINE) !----------------------------------------------------------------------------- end program ESMF_FieldRedistArbUTest