module buffer !----------------------------------------------------------------------- ! ! Purpose: ! LOW level handler for f90 arrays. ! ! Author: J. Edwards ! ! This file is used with genf90.pl to generate buffer.F90 ! !----------------------------------------------------------------------- use shr_kind_mod, only: r8 => shr_kind_r8, r4=> shr_kind_r4, i4=> shr_kind_i4 use cam_logfile, only: iulog use abortutils, only: endrun implicit none private ! The maximum number of dims in a fortran array #define MAXDIMS 7 type buffer_field_default_type private real(r8), pointer :: data(:,:,:,:,:,:,:) => null() end type buffer_field_default_type ! TYPE int,double,real type buffer_field_{TYPE} private {VTYPE}, pointer :: data(:,:,:,:,:,:,:) => null() end type buffer_field_{TYPE} integer(i4), parameter,public :: dtype_i4=1 real(r8), parameter,public :: dtype_r8=1_r8 real(r4), parameter,public :: dtype_r4=1_r4 interface buffer_field_allocate ! TYPE int,double,real module procedure buffer_field_allocate_{TYPE} end interface interface buffer_set_field ! TYPE int,double,real module procedure buffer_set_field_const_{TYPE} ! DIMS 1,2,3,4,5,6,7 ! TYPE int,double,real module procedure buffer_set_field_{DIMS}d_{TYPE} end interface interface buffer_get_field_ptr ! DIMS 1,2,3,4,5,6,7 ! TYPE int,double,real module procedure buffer_get_field_ptr_{DIMS}d_{TYPE} end interface public :: buffer_field_deallocate, buffer_field_allocate, buffer_set_field, buffer_get_field_ptr, buffer_field_default_type CONTAINS subroutine buffer_field_deallocate(bfg) type(buffer_field_default_type),intent(inout) :: bfg if(.not.associated(bfg%data)) then call endrun('Attempt to deallocate unassociated array ptr') end if deallocate(bfg%data) nullify(bfg%data) end subroutine buffer_field_deallocate ! TYPE int,double,real subroutine buffer_field_allocate_{TYPE} (bfg, dimsizes, dtype) type(buffer_field_default_type),intent(inout) :: bfg integer, intent(in) :: dimsizes(:) integer :: alldimsizes( MAXDIMS ) {VTYPE}, intent(in) :: dtype integer :: ierr type(buffer_field_{TYPE}) :: b1 alldimsizes(:) = 1 alldimsizes(1:size(dimsizes)) = dimsizes if(associated(bfg%data)) then call endrun('Attempt to allocate array to associated ptr') end if allocate(b1%data(alldimsizes(1),alldimsizes(2),alldimsizes(3),alldimsizes(4),& alldimsizes(5),alldimsizes(6),alldimsizes(7)),stat=ierr) if(ierr/=0) then call endrun("allocate failed") end if bfg = transfer(b1,bfg) end subroutine buffer_field_allocate_{TYPE} ! TYPE int,double,real subroutine buffer_set_field_const_{TYPE}(bfg, const, start,kount) type(buffer_field_default_type) :: bfg {VTYPE}, intent(in) :: const integer, intent(in), optional :: start(:),kount(:) type(buffer_field_{TYPE}) :: ptr integer :: i, ns, strt(7), fin(7), cnt(7) ptr = transfer(bfg, ptr) if(present(start).and.present(kount)) then strt(:) = 1 cnt = shape(ptr%data) ns = size(start) strt(1:ns) = start fin = strt+cnt-1 do i=1,ns fin(i) = strt(i)+kount(i)-1 if(strt(i)<1 .or. fin(i)>cnt(i)) then call endrun('Start plus kount exceeds dimension bounds') endif enddo ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),& strt(5):fin(5),strt(6):fin(6),strt(7):fin(7))=const else ptr%data = const endif end subroutine buffer_set_field_const_{TYPE} !========================================================================================= ! ! Given a physics_buffer chunk and an index return a pointer to a field chunk ! ! ! TYPE int,double,real ! DIMS 1,2,3,4,5,6,7 subroutine buffer_get_field_ptr_{DIMS}d_{TYPE}(bfg, field, start,kount) type(buffer_field_default_type), intent(in) :: bfg {VTYPE}, pointer :: field{DIMSTR} integer, intent(in), optional :: start(:), kount(:) type(buffer_field_{TYPE}), target :: ptr integer :: ns, strt(7), fin(7), cnt(7) ptr = transfer(bfg, ptr) strt(:) = 1 cnt = shape(ptr%data) if(present(start)) then ns = size(start) strt(1:ns) = start end if if(present(kount)) then cnt(1:ns) = kount end if fin = strt+cnt-1 #if ({DIMS}==1) field => ptr%data(strt(1):fin(1),strt(2),strt(3),strt(4),strt(5),strt(6),strt(7)) #elif ({DIMS}==2) field => ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3),strt(4),strt(5),strt(6),strt(7)) #elif ({DIMS}==3) field => ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4),strt(5),strt(6),strt(7)) #elif ({DIMS}==4) field => ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),strt(5),strt(6),strt(7)) #elif ({DIMS}==5) field => ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),strt(5):fin(5),strt(6),strt(7)) #elif ({DIMS}==6) field => ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),& strt(5):fin(5),strt(6):fin(6),strt(7)) #else field => ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),& strt(5):fin(5),strt(6):fin(6),strt(7):fin(7)) #endif end subroutine buffer_get_field_ptr_{DIMS}d_{TYPE} ! TYPE int,double,real ! DIMS 1,2,3,4,5,6,7 subroutine buffer_set_field_{DIMS}d_{TYPE}(bfg,field,start,kount) type(buffer_field_default_type) :: bfg {VTYPE},intent(in) :: field{DIMSTR} integer,intent(in),optional :: start(:),kount(:) type(buffer_field_{TYPE}) :: ptr integer :: i, nc, strt(7), fin(7), cnt(7) ptr = transfer(bfg,ptr) if(present(start).and.present(kount)) then strt(:) = 1 cnt = shape(ptr%data) nc=size(start) strt(1:nc) = start fin = strt+cnt-1 do i=1,nc fin(i) = strt(i)+kount(i)-1 if(strt(i)<1 .or. fin(i)>cnt(i)) then call endrun('Start plus kount exceeds dimension bounds') endif enddo #if ({DIMS}==1) ptr%data(strt(1):fin(1),strt(2),strt(3),strt(4),strt(5),strt(6),strt(7))=field #elif ({DIMS}==2) ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3),strt(4),strt(5),strt(6),strt(7))=field #elif ({DIMS}==3) ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4),strt(5),strt(6),strt(7))=field #elif ({DIMS}==4) ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),strt(5),strt(6),strt(7))=field #elif ({DIMS}==5) ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),strt(5):fin(5),strt(6),strt(7))=field #elif ({DIMS}==6) ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),& strt(5):fin(5),strt(6):fin(6),strt(7))=field #else ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),& strt(5):fin(5),strt(6):fin(6),strt(7):fin(7))=field #endif else #if ({DIMS}==1) ptr%data(:,1,1,1,1,1,1) = field #elif({DIMS}==2) ptr%data(:,:,1,1,1,1,1) = field #elif({DIMS}==3) ptr%data(:,:,:,1,1,1,1) = field #elif({DIMS}==4) ptr%data(:,:,:,:,1,1,1) = field #elif({DIMS}==5) ptr%data(:,:,:,:,:,1,1) = field #elif({DIMS}==6) ptr%data(:,:,:,:,:,:,1) = field #else ptr%data = field #endif end if end subroutine buffer_set_field_{DIMS}d_{TYPE} end module buffer