XBeach
C:/repositories/XBeach/trunk/src/xbeachlibrary/xbeach_bmi.f90
Go to the documentation of this file.
00001 module xbeach_bmi
00002    use iso_c_binding,    only: c_loc, c_int, c_char, c_double, c_ptr
00003    use iso_c_utils,      only: char_array_to_string, string_to_char_array, strlen, MAXSTRINGLEN
00004    use libxbeach_module, only: par, sglobal, s, executestep, outputext, init, final
00005 
00006    use libxbeach_module
00007    use timestep_module,  only: compute_dt
00008    use logging_module,   only: LEVEL_INFO, logmsg
00009    use spaceparamsdef,   only: spacepars
00010    use spaceparams,      only: index_allocated, indextos, index_allocate
00011    use mnemmodule,       only: chartoindex, arraytype
00012 
00013    use spaceparams
00014 
00015    use xmpi_module
00016 
00017    implicit none
00018    save
00019    private
00020 
00021    ! initialize  -> loadmodel/initmodel
00022    ! perform single timestep -> update
00023    ! finalize -> finalize
00024 
00025 
00026 
00027    ! This is assumed.....
00028    integer(c_int), parameter :: MAXDIMS = 6
00029 
00030 
00031 contains
00032 
00033 
00034 
00035    integer(c_int) function finalize() result(ierr) bind(C, name="finalize")
00036       !DEC$ ATTRIBUTES DLLEXPORT::finalize
00037       ierr = 0
00038       call logmsg(LEVEL_INFO, 'Finalize')
00039       ierr = final()
00040    end function finalize
00041 
00042 
00043    integer(c_int) function initialize(c_configfile) result(ierr) bind(C, name="initialize")
00044       !DEC$ ATTRIBUTES DLLEXPORT::initialize
00045 
00046       implicit none
00047 
00048       ! Variables
00049       character(kind=c_char), intent(in) :: c_configfile(*)
00050       character(len=strlen(c_configfile)) :: configfile
00051 
00052       ! Convert c string to fortran string
00053       configfile = char_array_to_string(c_configfile)
00054 
00055       ierr = init()
00056 
00057    end function initialize
00058 
00059 
00061    integer(c_int) function update(dt) result(ierr) bind(C,name="update")
00062       !DEC$ ATTRIBUTES DLLEXPORT::update
00063 
00064       !< Custom timestep size, use -1 to use model default.
00065       real(c_double), value, intent(in) :: dt
00066 
00067       if (dt >= 0) then
00068          ierr = executestep(dt)
00069       else
00070          ierr = executestep()
00071       end if
00072 
00073       ! enable output to be able to use aggregated variable values (min,
00074       ! max, mean, var) through BMI. user can still control output
00075       ! through params.txt file
00076       ierr = outputext()
00077 
00078    end function update
00079 
00080 
00081    ! Void function is a subroutine
00082    subroutine get_var_type(c_var_name, c_type_name)  bind(C, name="get_var_type")
00083       !DEC$ ATTRIBUTES DLLEXPORT :: get_var_type
00084 
00085       character(kind=c_char), intent(in) :: c_var_name(*)
00086       character(kind=c_char), intent(out) :: c_type_name(MAXSTRINGLEN)
00087 
00088       character(len=strlen(c_var_name)) :: var_name
00089       character(len=MAXSTRINGLEN) :: type_name
00090 
00091       character :: typecode
00092       integer :: index
00093       type(arraytype) :: array
00094 
00095 
00096       var_name = char_array_to_string(c_var_name)
00097 
00098 
00099       index =  chartoindex(var_name)
00100       if (index .eq. -1) return
00101       call indextos(s,index,array)
00102       typecode = array%type
00103 
00104       select case(typecode)
00105        case("r")
00106          type_name = "double"
00107        case("i")
00108          type_name = "int"
00109        case("c")
00110          type_name = "character"
00111        case default
00112          continue
00113       end select
00114       c_type_name = string_to_char_array(trim(type_name))
00115 
00116    end subroutine get_var_type
00117 
00118    subroutine get_var_rank(c_var_name, rank) bind(C, name="get_var_rank")
00119       !DEC$ ATTRIBUTES DLLEXPORT :: get_var_rank
00120 
00121       character(kind=c_char), intent(in) :: c_var_name(*)
00122       integer(c_int), intent(out) :: rank
00123 
00124       ! The fortran name of the attribute name
00125       character(len=strlen(c_var_name)) :: var_name
00126       type(arraytype) :: array
00127       integer :: index
00128 
00129       ! Store the name
00130       var_name = char_array_to_string(c_var_name)
00131 
00132       index =  chartoindex(var_name)
00133       if (index .eq. -1) return
00134       call indextos(s,index,array)
00135       rank = array%rank
00136    end subroutine get_var_rank
00137 
00138    include 'get_var_shape.inc'
00139    include 'get_var.inc'
00140    include 'set_var.inc'
00141    subroutine set_current_time(xptr) bind(C, name="set_current_time")
00142       !DEC$ ATTRIBUTES DLLEXPORT :: set_current_time
00143 
00144       use iso_c_binding, only: c_ptr, c_double, c_f_pointer
00145 
00146       type(c_ptr), value, intent(in) :: xptr
00147       real(c_double), pointer :: x_0d_double_ptr
00148 
00149       call c_f_pointer(xptr, x_0d_double_ptr)
00150       if (par%morfacopt == 1) then
00151          par%t = x_0d_double_ptr / max(par%morfac, 1.d0)
00152       else
00153          par%t = x_0d_double_ptr
00154       endif
00155 
00156    end subroutine set_current_time
00157 
00158    subroutine get_current_time(time) bind(C, name="get_current_time")
00159       !DEC$ ATTRIBUTES DLLEXPORT :: get_current_time
00160 
00161       real(c_double) :: time
00162 
00163       if (par%morfacopt == 1) then
00164          time = par%t * max(par%morfac, 1.d0)
00165       else
00166          time = par%t
00167       endif
00168    end subroutine get_current_time
00169 
00170    subroutine get_start_time(time) bind(C, name="get_start_time")
00171       !DEC$ ATTRIBUTES DLLEXPORT :: get_start_time
00172 
00173       real(c_double) :: time
00174 
00175       time = 0
00176    end subroutine get_start_time
00177 
00178    subroutine get_time_step(timestep) bind(C, name="get_time_step")
00179       !DEC$ ATTRIBUTES DLLEXPORT :: get_time_step
00180 
00181       real(c_double) :: timestep
00182       integer :: ilim = 0
00183       integer :: jlim = 0
00184           integer :: it = 0
00185       real*8 :: dtref = 0.d0
00186 
00187       call compute_dt(s,par, tpar, it, ilim, jlim, dtref)
00188       timestep = par%dt
00189    end subroutine get_time_step
00190 
00191    subroutine get_end_time(time) bind(C, name="get_end_time")
00192       !DEC$ ATTRIBUTES DLLEXPORT :: get_end_time
00193 
00194       real(c_double) :: time
00195       if (par%morfacopt == 1) then
00196          time = par%tstop * max(par%morfac, 1.d0)
00197       else
00198          time = par%tstop
00199       endif
00200    end subroutine get_end_time
00201 
00202 
00203 end module xbeach_bmi
 All Classes Files Functions Variables Defines