XBeach
|
00001 module wave_timestep_module 00002 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00003 ! Copyright (C) 2007 UNESCO-IHE, WL|Delft Hydraulics and Delft University ! 00004 ! Dano Roelvink, Ap van Dongeren, Ad Reniers, Jamie Lescinski, ! 00005 ! Jaap van Thiel de Vries, Robert McCall ! 00006 ! ! 00007 ! d.roelvink@unesco-ihe.org ! 00008 ! UNESCO-IHE Institute for Water Education ! 00009 ! P.O. Box 3015 ! 00010 ! 2601 DA Delft ! 00011 ! The Netherlands ! 00012 ! ! 00013 ! This library is free software; you can redistribute it and/or ! 00014 ! modify it under the terms of the GNU Lesser General Public ! 00015 ! License as published by the Free Software Foundation; either ! 00016 ! version 2.1 of the License, or (at your option) any later version. ! 00017 ! ! 00018 ! This library is distributed in the hope that it will be useful, ! 00019 ! but WITHOUT ANY WARRANTY; without even the implied warranty of ! 00020 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU ! 00021 ! Lesser General Public License for more details. ! 00022 ! ! 00023 ! You should have received a copy of the GNU Lesser General Public ! 00024 ! License along with this library; if not, write to the Free Software ! 00025 ! Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 ! 00026 ! USA ! 00027 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00028 implicit none 00029 save 00030 private 00031 public:: wave 00032 00033 contains 00034 subroutine wave(s,par) 00035 00036 use params, only: parameters 00037 use spaceparams 00038 !use wave_stationary_module 00039 !use wave_directions_module 00040 use wave_stationary_directions_module 00041 use wave_instationary_module 00042 use paramsconst 00043 use wave_functions_module, only: wave_dispersion,update_means_wave_flow 00044 00045 implicit none 00046 00047 type(spacepars) :: s 00048 type(parameters) :: par 00049 00050 real*8,dimension(:,:),allocatable,save :: gamma 00051 00052 if (.not.allocated(gamma)) then 00053 allocate(gamma(s%nx+1,s%ny+1)) 00054 endif 00055 00056 ! set basic water depth for all wave calculations: note that depending on wci setting, other water depths may be used in 00057 ! dispersion routine etc. 00058 if (par%delta>0.d0) then 00059 s%hhw = max(s%hh+par%delta*s%H,par%eps) 00060 else 00061 s%hhw = max(s%hh,par%eps) ! hh can be less than eps after morphevolution? 00062 endif 00063 00064 if (par%oldhmin==1) then 00065 s%hstokes = max(s%hh,par%hmin) 00066 else 00067 gamma = s%H/s%hh 00068 where (gamma>1.d0) 00069 s%hstokes = par%deltahmin*(gamma-1.d0)*s%H+s%hh 00070 elsewhere 00071 s%hstokes = s%hh 00072 endwhere 00073 endif 00074 00075 select case (par%wavemodel) 00076 case(WAVEMODEL_STATIONARY) 00077 ! call update_means_wave_flow not required: stationary always uses instantaneous water depth s%hhw 00078 if ((abs(mod(par%t,par%wavint))<0.001d0*par%dt) .or. s%newstatbc==1) then 00079 call wave_dispersion(s,par,0) ! use instantaneous water depth (and velocity) 00080 call wave_stationary_directions(s,par,0) 00081 s%newstatbc = 0 00082 endif 00083 case(WAVEMODEL_SURFBEAT) 00084 if (par%single_dir==1) then 00085 ! always need to update depths for wave directions model, as well as in case of wci (dealt with in subroutine) 00086 call update_means_wave_flow(s,par) 00087 ! 00088 ! if necessary, update wave directions 00089 if ((abs(mod(par%t,par%wavint))<0.001d0*par%dt) .or. s%newstatbc==1 .or. par%t==par%dt) then 00090 call wave_dispersion(s,par,1) ! use s%hhws water depth (and velocity) 00091 call wave_stationary_directions(s,par,1) 00092 s%newstatbc = 0 00093 endif 00094 ! 00095 s%newstatbc = 0 ! not sure if this is needed every timestep, but no overhead to keep in ... 00096 ! 00097 ! need to call dispersion again, because dispersion is different in timestep and direction computation 00098 if (par%wci==1) then 00099 call wave_dispersion(s,par,2) ! use s%hhwcins water depth (and velocity) 00100 else 00101 call wave_dispersion(s,par,0) ! use s%hhw water depth 00102 endif 00103 call wave_instationary(s,par) 00104 else 00105 s%newstatbc = 0 00106 if (par%wci==1) then 00107 ! in this case we only need to update flow depth in case of wci 00108 call update_means_wave_flow(s,par) 00109 call wave_dispersion(s,par,2) ! use s%hhwcins water depth (and velocity) 00110 else 00111 call wave_dispersion(s,par,0) ! use s%hhw water depth 00112 endif 00113 call wave_instationary(s,par) 00114 endif 00115 end select 00116 00117 00118 end subroutine wave 00119 00120 end module wave_timestep_module