XBeach
|
00001 module rainfall_module 00002 implicit none 00003 save 00004 00005 logical,save,private :: constantRainfall 00006 00007 contains 00008 00009 subroutine rainfall_init(s,par) 00010 use params 00011 use spaceparams 00012 use filefunctions 00013 use logging_module 00014 use interp 00015 00016 type(spacepars) :: s 00017 type(parameters) :: par 00018 00019 integer :: i,j,it 00020 integer :: ier,fid,dummy 00021 00022 if(.not. xmaster) return 00023 00024 ! allocate global rainfall rate variable 00025 allocate(s%rainfallrate(s%nx+1,s%ny+1)) 00026 00027 if (par%rainfall==1) then 00028 ! read from file? 00029 if (par%rainfallratefile==' ') then 00030 constantRainfall = .true. 00031 allocate(s%rainfallinput(0,0,0)) 00032 allocate(s%trainfallinput(0)) 00033 s%rainfallrate = par%rainfallrate/1000.d0/3600.d0 ! convert from [mm/hr] to [m/s] 00034 else 00035 constantRainfall = .false. 00036 allocate(s%rainfallinput(s%nx+1,s%ny+1,par%nrainfallrate)) 00037 allocate(s%trainfallinput(par%nrainfallrate)) 00038 ! start file read 00039 fid = create_new_fid() 00040 ! call check_file_exist(par%rainfallratefile) ! already done at all_input subroutine level 00041 open (fid,file=par%rainfallratefile) 00042 do it=1,par%nrainfallrate 00043 read(fid,*,iostat=ier)s%trainfallinput(it) 00044 if (ier .ne. 0) then 00045 call report_file_read_error(par%rainfallratefile) 00046 endif 00047 do j=1,s%ny+1 00048 read(fid,*,iostat=ier)(s%rainfallinput(i,j,it),i=1,s%nx+1) 00049 if (ier .ne. 0) then 00050 call report_file_read_error(par%rainfallratefile) 00051 endif 00052 enddo 00053 enddo 00054 close(fid) 00055 s%rainfallinput = s%rainfallinput/1000.d0/3600.d0 ! convert source file from [mm/hr] to [m/s] 00056 ! Interpolate initial rainfall 00057 do j=1,s%ny+1 00058 do i=1,s%nx+1 00059 call LINEAR_INTERP(s%trainfallinput,s%rainfallinput(i,j,:),par%nrainfallrate, & 00060 0.d0,s%rainfallrate(i,j),dummy) 00061 enddo 00062 enddo 00063 endif 00064 else 00065 ! give MPI bcast a memory address 00066 allocate(s%rainfallinput(0,0,0)) 00067 allocate(s%trainfallinput(0)) 00068 constantRainfall = .true. 00069 s%rainfallrate = 0.d0 00070 endif 00071 end subroutine rainfall_init 00072 00073 subroutine rainfall_update(s, par) 00074 00075 use params 00076 use spaceparams 00077 use interp 00078 00079 implicit none 00080 00081 type(spacepars) :: s 00082 type(parameters) :: par 00083 00084 integer :: i,j,dummy 00085 00086 if (.not. constantRainfall) then 00087 ! interpolate from time series 00088 do j=1,s%ny+1 00089 do i=1,s%nx+1 00090 call LINEAR_INTERP(s%trainfallinput,s%rainfallinput(i,j,:),par%nrainfallrate, & 00091 par%t,s%rainfallrate(i,j),dummy) 00092 enddo 00093 enddo 00094 endif 00095 00096 end subroutine rainfall_update 00097 00098 end module rainfall_module