#ifndef __unix__ #define DO_DLL #endif module xbeach_bmi use iso_c_binding, only: c_loc, c_int, c_char, c_double, c_ptr use iso_c_utils, only: char_array_to_string, string_to_char_array, strlen, MAXSTRINGLEN use libxbeach_module, only: par, sglobal, s, executestep, outputext, init, final #ifdef USEMPI use libxbeach_module, only: slocal #endif use logging_module, only: LEVEL_INFO, logmsg use spaceparamsdef, only: spacepars use spaceparams, only: index_allocated, indextos, index_allocate use mnemmodule, only: chartoindex, arraytype #ifdef USEMPI use spaceparams, only: space_collect_mnem, space_distribute_mnem #endif use xmpi_module, only: using_mpi, xmaster, halt_program, xmpi_ocomm implicit none save private ! initialize -> loadmodel/initmodel ! perform single timestep -> update ! finalize -> finalize ! This is assumed..... integer(c_int), parameter :: MAXDIMS = 6 contains integer(c_int) function finalize() result(ierr) bind(C, name="finalize") #ifdef DO_DLL !DEC$ ATTRIBUTES DLLEXPORT::finalize #endif ierr = 0 call logmsg(LEVEL_INFO, 'Finalize') ierr = final() end function finalize integer(c_int) function initialize(c_configfile) result(ierr) bind(C, name="initialize") #ifdef DO_DLL !DEC$ ATTRIBUTES DLLEXPORT::initialize #endif implicit none ! Variables character(kind=c_char), intent(in) :: c_configfile(*) character(len=strlen(c_configfile)) :: configfile ! Convert c string to fortran string configfile = char_array_to_string(c_configfile) ierr = init() end function initialize !> Performs a single timestep with the current model. integer(c_int) function update(dt) result(ierr) bind(C,name="update") #ifdef DO_DLL !DEC$ ATTRIBUTES DLLEXPORT::update #endif !< Custom timestep size, use -1 to use model default. real(c_double), value, intent(in) :: dt if (dt < -1) then par%t = par%t - dt / par%morfac ierr = 0 ! added by wwvv elseif (dt >= 0) then ierr = executestep(dt / par%morfac) else ierr = executestep() end if ierr = ierr + outputext() ! wwvv: added end function update ! Void function is a subroutine subroutine get_var_type(c_var_name, c_type_name) bind(C, name="get_var_type") #ifdef DO_DLL !DEC$ ATTRIBUTES DLLEXPORT :: get_var_type #endif character(kind=c_char), intent(in) :: c_var_name(*) character(kind=c_char), intent(out) :: c_type_name(MAXSTRINGLEN) character(len=strlen(c_var_name)) :: var_name character(len=MAXSTRINGLEN) :: type_name character :: typecode integer :: index type(arraytype) :: array var_name = char_array_to_string(c_var_name) index = chartoindex(var_name) if (index .eq. -1) return call indextos(s,index,array) typecode = array%type select case(typecode) case("r") type_name = "double" case("i") type_name = "int" case("c") type_name = "character" case default continue end select c_type_name = string_to_char_array(trim(type_name)) end subroutine get_var_type subroutine get_var_rank(c_var_name, rank) bind(C, name="get_var_rank") #ifdef DO_DLL !DEC$ ATTRIBUTES DLLEXPORT :: get_var_rank #endif character(kind=c_char), intent(in) :: c_var_name(*) integer(c_int), intent(out) :: rank ! The fortran name of the attribute name character(len=strlen(c_var_name)) :: var_name type(arraytype) :: array integer :: index ! Store the name var_name = char_array_to_string(c_var_name) index = chartoindex(var_name) if (index .eq. -1) return call indextos(s,index,array) rank = array%rank end subroutine get_var_rank subroutine get_var_shape(c_var_name, shape) bind(C, name="get_var_shape") #ifdef DO_DLL !DEC$ ATTRIBUTES DLLEXPORT :: get_var_shape #endif implicit none character(kind=c_char), intent(in) :: c_var_name(*) integer(c_int), intent(inout) :: shape(MAXDIMS) character(len=strlen(c_var_name)) :: var_name type(spacepars), pointer :: spointer var_name = char_array_to_string(c_var_name) shape = (/0, 0, 0, 0, 0, 0/) if (xmaster) then if (using_mpi) then spointer => sglobal else spointer => s endif select case(var_name) include 'get_var_shape.inc' ! statements like: ! case("zs") ! shape(2) = spointer%nx+1 ! shape(1) = spointer%ny+1 case default write(*,*) 'Getting shape of unknown variable', var_name call halt_program end select endif end subroutine get_var_shape subroutine get_var(c_var_name, x) bind(C, name="get_var") #ifdef DO_DLL !DEC$ ATTRIBUTES DLLEXPORT :: get_var #endif ! Return a pointer to the variable character(kind=c_char), intent(in) :: c_var_name(*) type(c_ptr), intent(inout) :: x character(len=strlen(c_var_name)) :: var_name type(spacepars), pointer :: spointer var_name = char_array_to_string(c_var_name) #ifdef USEMPI ! collect the data of var_name into sglobal: call space_collect_mnem(sglobal,slocal,par,var_name,force = .true.,w_only=.true.) #endif if (xmaster) then if (using_mpi) then spointer => sglobal else spointer => s endif select case(var_name) include 'get_var.inc' ! statements like: ! case ('zs') ! x = c_loc(spointer%zs(lbound(s%zs, 1),lbound(s%zs, 2))) case default write(*,*) 'Getting unknown variable', var_name call halt_program end select endif end subroutine get_var ! In the mpi case: ! this subroutine is called by the compute processes (i.e. all processes mpi_comm) ! the data to be distributed is available on process 0 in mpi_comm: the process ! where xmaster .eq. .true. subroutine set_var(c_var_name, xptr) bind(C, name="set_var") #ifdef DO_DLL !DEC$ ATTRIBUTES DLLEXPORT :: set_var #endif ! Return a pointer to the variable use iso_c_binding, only: c_double, c_char, c_loc, c_f_pointer character(kind=c_char), intent(in) :: c_var_name(*) type(c_ptr), value, intent(in) :: xptr real(c_double), pointer :: x_0d_double_ptr real(c_double), pointer :: x_1d_double_ptr(:) real(c_double), pointer :: x_2d_double_ptr(:,:) real(c_double), pointer :: x_3d_double_ptr(:,:,:) real(c_double), pointer :: x_4d_double_ptr(:,:,:,:) integer(c_int), pointer :: x_0d_int_ptr integer(c_int), pointer :: x_1d_int_ptr(:) integer(c_int), pointer :: x_2d_int_ptr(:,:) !integer(c_int), pointer :: x_3d_int_ptr(:,:,:) !real(c_float), pointer :: x_1d_float_ptr(:) !real(c_float), pointer :: x_2d_float_ptr(:,:) !real(c_float), pointer :: x_3d_float_ptr(:,:,:) ! The fortran name of the attribute name character(len=strlen(c_var_name)) :: var_name type(spacepars), pointer :: spointer integer :: index #ifdef USEMPI integer, pointer :: ipointer #endif var_name = char_array_to_string(c_var_name) #ifdef USEMPI ! The mpi communicator has been determined in the calling program: if (var_name .eq. 'mpicomm') then call c_f_pointer(xptr, ipointer) xmpi_ocomm = ipointer return endif #endif if (xmaster) then ! always true in serial program if (using_mpi) then spointer => sglobal else spointer => s endif index = chartoindex(var_name) if (.not. index_allocated(spointer,index)) then call index_allocate(spointer,par,index,'a') endif select case(var_name) include 'set_var.inc' ! includes statements like: !case("zs") ! call c_f_pointer(xptr, x_2d_double_ptr, shape(spointer%zs)) ! spointer%zs(:,:) = x_2d_double_ptr case default write(*,*) 'Setting unknown variable', var_name call halt_program end select endif #ifdef USEMPI ! the variable is now available om process 0 in mpi_comm. ! space_distribute_mnem is designed to dirstribute from this ! process to the compute processes (processes in mpi_comm) call space_distribute_mnem(sglobal,slocal,var_name) #endif end subroutine set_var subroutine get_current_time(time) bind(C, name="get_current_time") #ifdef DO_DLL !DEC$ ATTRIBUTES DLLEXPORT :: get_current_time #endif real(c_double) :: time time = par%t * par%morfac end subroutine get_current_time subroutine get_start_time(time) bind(C, name="get_start_time") #ifdef DO_DLL !DEC$ ATTRIBUTES DLLEXPORT :: get_start_time #endif real(c_double) :: time time = 0 end subroutine get_start_time subroutine get_end_time(time) bind(C, name="get_end_time") #ifdef DO_DLL !DEC$ ATTRIBUTES DLLEXPORT :: get_end_time #endif real(c_double) :: time time = par%tstop * par%morfac end subroutine get_end_time end module xbeach_bmi